3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
Public Types | Public Member Functions | Protected Attributes | List of all members
Dumux::LocalStiffness< TypeTag, m > Class Template Referenceabstract

Base class for local assemblers. More...

#include <dumux/porousmediumflow/2p/sequential/diffusion/mimetic/localstiffness.hh>

Inheritance diagram for Dumux::LocalStiffness< TypeTag, m >:
Inheritance graph

Description

template<class TypeTag, int m>
class Dumux::LocalStiffness< TypeTag, m >

Base class for local assemblers.

This class serves as a base class for local assemblers. It provides space and access to the local stiffness matrix. The actual assembling is done in a derived class via the virtual assemble method.

Template Parameters
TypeTagThe problem TypeTag
mnumber of degrees of freedom per node (system size)

Public Types

using MBlockType = Dune::FieldMatrix< Scalar, m, m >
 
using VBlockType = Dune::FieldVector< Scalar, m >
 
using BCBlockType = std::array< BoundaryConditions::Flags, m >
 
using BoundaryTypes = GetPropType< TypeTag, Properties::BoundaryTypes >
 

Public Member Functions

virtual ~LocalStiffness ()
 
virtual void assemble (const Entity &e, int k=1)=0
 Assembles local stiffness matrix including boundary conditions for given element and order. More...
 
virtual void assemble (const Entity &e, const Dune::BlockVector< VBlockType > &localSolution, int k=1)=0
 assemble local stiffness matrix including boundary conditions for given element and order More...
 
virtual void assembleBoundaryCondition (const Entity &e, int k=1)=0
 assemble only boundary conditions for given element and order More...
 
void print (std::ostream &s, int width, int precision)
 Prints contents of local stiffness matrix. More...
 
const MBlockTypemat (int i, int j) const
 Accesses the local stiffness matrix. More...
 
const VBlockTyperhs (int i) const
 Accesses right hand side. More...
 
const BoundaryTypesbc (int i) const
 Accesses boundary condition for each dof. More...
 
void setcurrentsize (int s)
 Sets the current size of the local stiffness matrix. More...
 
int currentsize ()
 Gets the current size of the local stiffness matrix. More...
 

Protected Attributes

Dune::Matrix< MBlockTypeA
 
std::vector< VBlockTypeb
 
std::vector< BoundaryTypesbctype
 

Member Typedef Documentation

◆ BCBlockType

template<class TypeTag , int m>
using Dumux::LocalStiffness< TypeTag, m >::BCBlockType = std::array<BoundaryConditions::Flags, m>

◆ BoundaryTypes

template<class TypeTag , int m>
using Dumux::LocalStiffness< TypeTag, m >::BoundaryTypes = GetPropType<TypeTag, Properties::BoundaryTypes>

◆ MBlockType

template<class TypeTag , int m>
using Dumux::LocalStiffness< TypeTag, m >::MBlockType = Dune::FieldMatrix<Scalar, m, m>

◆ VBlockType

template<class TypeTag , int m>
using Dumux::LocalStiffness< TypeTag, m >::VBlockType = Dune::FieldVector<Scalar, m>

Constructor & Destructor Documentation

◆ ~LocalStiffness()

template<class TypeTag , int m>
virtual Dumux::LocalStiffness< TypeTag, m >::~LocalStiffness ( )
inlinevirtual

Member Function Documentation

◆ assemble() [1/2]

template<class TypeTag , int m>
virtual void Dumux::LocalStiffness< TypeTag, m >::assemble ( const Entity &  e,
const Dune::BlockVector< VBlockType > &  localSolution,
int  k = 1 
)
pure virtual

assemble local stiffness matrix including boundary conditions for given element and order

Unlike the method with only two arguments, this one additionally takes the local solution in order to allow assembly of nonlinear operators.

On exit the following things have been done:

  • The stiffness matrix for the given entity and polynomial degree has been assembled and is accessible with the mat() method.
  • The boundary conditions have been evaluated and are accessible with the bc() method. The boundary conditions are either Neumann, process or Dirichlet. Neumann indicates that the corresponding node (assuming a nodal basis) is at the Neumann boundary, process indicates that the node is at a process boundary (arising from the parallel decomposition of the mesh). Process boundaries are treated as homogeneous Dirichlet conditions, i.e. the corresponding value in the right hand side is set to 0. Finally, Dirichlet indicates that the node is at the Dirichlet boundary.
  • The right hand side has been assembled. It contains either the value of the essential boundary condition or the assembled source term and Neumann boundary condition. It is accessible via the rhs() method.
Parameters
ea codim 0 entity reference
localSolutionThe current solution on the entity, which is needed by nonlinear assemblers
korder of Lagrange basis (default is 1)

Implemented in Dumux::LinearLocalStiffness< TypeTag, m >.

◆ assemble() [2/2]

template<class TypeTag , int m>
virtual void Dumux::LocalStiffness< TypeTag, m >::assemble ( const Entity &  e,
int  k = 1 
)
pure virtual

Assembles local stiffness matrix including boundary conditions for given element and order.

On exit the following things have been done:

  • The stiffness matrix for the given entity and polynomial degree has been assembled and is accessible with the mat() method.
  • The boundary conditions have been evaluated and are accessible with the bc() method. The boundary conditions are either neumann, process or dirichlet. Neumann indicates that the corresponding node (assuming a nodal basis) is at the Neumann boundary, process indicates that the node is at a process boundary (arising from the parallel decomposition of the mesh). Process boundaries are treated as homogeneous Dirichlet conditions, i.e. the corresponding value in the right hand side is set to 0. Finally, Dirichlet indicates that the node is at the Dirichlet boundary.
  • The right hand side has been assembled. It contains either the value of the essential boundary condition or the assembled source term and neumann boundary condition. It is accessible via the rhs() method.
Parameters
ea codim 0 entity reference
korder of Lagrange basis (default is 1)

Implemented in Dumux::LinearLocalStiffness< TypeTag, m >.

◆ assembleBoundaryCondition()

template<class TypeTag , int m>
virtual void Dumux::LocalStiffness< TypeTag, m >::assembleBoundaryCondition ( const Entity &  e,
int  k = 1 
)
pure virtual

assemble only boundary conditions for given element and order

On exit the following things have been done:

  • The boundary conditions have been evaluated and are accessible with the bc() method. The boundary conditions are either Neumann, process or Dirichlet. Neumann indicates that the corresponding node (assuming a nodal basis) is at the Neumann boundary, process indicates that the node is at a process boundary (arising from the parallel decomposition of the mesh). Process boundaries are treated as homogeneous Dirichlet conditions, i.e. the corresponding value in the right hand side is set to 0. Finally, Dirichlet indicates that the node is at the Dirichlet boundary.
  • The right hand side has been assembled as far as boundary conditions are concerned. It contains either the value of the essential boundary condition or the assembled Neumann boundary condition. It is accessible via the rhs() method.
Parameters
ea codim 0 entity reference
korder of Lagrange basis (default is 1)

◆ bc()

template<class TypeTag , int m>
const BoundaryTypes & Dumux::LocalStiffness< TypeTag, m >::bc ( int  i) const
inline

Accesses boundary condition for each dof.

Access boundary condition type for each degree of freedom. Elements are undefined without prior call to the assemble method.

Parameters
iindex i

◆ currentsize()

template<class TypeTag , int m>
int Dumux::LocalStiffness< TypeTag, m >::currentsize ( )
inline

Gets the current size of the local stiffness matrix.

◆ mat()

template<class TypeTag , int m>
const MBlockType & Dumux::LocalStiffness< TypeTag, m >::mat ( int  i,
int  j 
) const
inline

Accesses the local stiffness matrix.

Access elements of the local stiffness matrix. Elements are undefined without prior call to the assemble method.

Parameters
iindex i
jindex j

◆ print()

template<class TypeTag , int m>
void Dumux::LocalStiffness< TypeTag, m >::print ( std::ostream &  s,
int  width,
int  precision 
)
inline

Prints contents of local stiffness matrix.

Parameters
soutput stream
widththe width
precisionthe precision

◆ rhs()

template<class TypeTag , int m>
const VBlockType & Dumux::LocalStiffness< TypeTag, m >::rhs ( int  i) const
inline

Accesses right hand side.

Access elements of the right hand side vector. Elements are undefined without prior call to the assemble method.

Parameters
iindex i

◆ setcurrentsize()

template<class TypeTag , int m>
void Dumux::LocalStiffness< TypeTag, m >::setcurrentsize ( int  s)
inline

Sets the current size of the local stiffness matrix.

Parameters
ssize which is to be set

Member Data Documentation

◆ A

template<class TypeTag , int m>
Dune::Matrix<MBlockType> Dumux::LocalStiffness< TypeTag, m >::A
protected

◆ b

template<class TypeTag , int m>
std::vector<VBlockType> Dumux::LocalStiffness< TypeTag, m >::b
protected

◆ bctype

template<class TypeTag , int m>
std::vector<BoundaryTypes> Dumux::LocalStiffness< TypeTag, m >::bctype
protected

The documentation for this class was generated from the following file: