- Home
- About DUNE
- Download
- Documentation
- Community
- Development
#include <localstiffness.hh>
Inheritance diagram for Dune::LocalStiffness< GV, RT, m >:

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.
GV A grid view type RT The field type used in the elements of the stiffness matrix m number of degrees of freedom per node (system size)
Public Member Functions | |
| virtual void | assemble (const Entity &e, int k=1)=0 |
| assemble local stiffness matrix including boundary conditions for given element and order | |
| virtual void | assemble (const Entity &e, const BlockVector< VBlockType > &localSolution, int k=1)=0 |
| assemble local stiffness matrix including boundary conditions for given element and order | |
| virtual void | assembleBoundaryCondition (const Entity &e, int k=1)=0 |
| assemble only boundary conditions for given element and order | |
| void | print (std::ostream &s, int width, int precision) |
| print contents of local stiffness matrix | |
| const MBlockType & | mat (int i, int j) const |
| access local stiffness matrix | |
| const VBlockType & | rhs (int i) const |
| access right hand side | |
| const BCBlockType & | bc (int i) const |
| access boundary condition for each dof | |
| void | setcurrentsize (int s) |
| set the current size of the local stiffness matrix | |
| int | currentsize () |
| get the current size of the local stiffness matrix | |
| virtual void Dune::LocalStiffness< GV, RT, m >::assemble | ( | const Entity & | e, | |
| int | k = 1 | |||
| ) | [pure virtual] |
assemble local stiffness matrix including boundary conditions for given element and order
On exit the following things have been done:
| [in] | e | a codim 0 entity reference |
| [in] | k | order of Lagrange basis (default is 1) |
Implemented in Dune::GroundwaterEquationLocalStiffness< GV, RT >, Dune::LinearLocalStiffness< GV, RT, m >, and Dune::LinearLocalStiffness< GV, RT, 1 >.
| virtual void Dune::LocalStiffness< GV, RT, m >::assemble | ( | const Entity & | e, | |
| const 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:
| [in] | e | a codim 0 entity reference |
| [in] | localSolution | The current solution on the entity, which is needed by nonlinear assemblers |
| [in] | k | order of Lagrange basis (default is 1) |
Implemented in Dune::LinearLocalStiffness< GV, RT, m >, and Dune::LinearLocalStiffness< GV, RT, 1 >.
| virtual void Dune::LocalStiffness< GV, RT, 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:
| [in] | e | a codim 0 entity reference |
| [in] | k | order of Lagrange basis (default is 1) |
Implemented in Dune::GroundwaterEquationLocalStiffness< GV, RT >.
| const MBlockType& Dune::LocalStiffness< GV, RT, m >::mat | ( | int | i, | |
| int | j | |||
| ) | const [inline] |
access local stiffness matrix
Access elements of the local stiffness matrix. Elements are undefined without prior call to the assemble method.
| const VBlockType& Dune::LocalStiffness< GV, RT, m >::rhs | ( | int | i | ) | const [inline] |
access right hand side
Access elements of the right hand side vector. Elements are undefined without prior call to the assemble method.
| const BCBlockType& Dune::LocalStiffness< GV, RT, m >::bc | ( | int | i | ) | const [inline] |
access boundary condition for each dof
Access boundary condition type for each degree of freedom. Elements are undefined without prior call to the assemble method.