Class handling finite element assembly. More...

#include <asmhandler.hpp>

Public Types

typedef
GridType::LeafGridView::IndexSet 
LeafIndexSet
 A set of indices. More...
 
typedef
GridType::LeafGridView::template
Codim< 0 >::Iterator 
LeafIterator
 An iterator over grid cells. More...
 
typedef Dune::FieldVector
< double, dim
NodeValue
 A vectorial node value. More...
 

Public Member Functions

 ASMHandler (const GridType &gv_)
 The default constructor. More...
 
 ~ASMHandler ()
 Destructor. More...
 
size_t getEqns () const
 Get the number of equations in the system. More...
 
int getEquationForDof (int node, int dof)
 Get the equation number for a given dof on a given node. More...
 
MatrixgetOperator ()
 Obtain a reference to the linear operator. More...
 
VectorgetLoadVector ()
 Obtain a reference to the load vector. More...
 
void initForAssembly ()
 This function needs to be called before starting the element assembly. More...
 
template<int esize>
void addElement (const Dune::FieldMatrix< double, esize, esize > *K, const Dune::FieldVector< double, esize > *S, const LeafIterator &cell, Vector *b=NULL)
 Add an element matrix/vector to the system. More...
 
template<int comp>
void extractValues (Dune::FieldVector< double, comp > &v, const Vector &u, const LeafIterator &it)
 Extract values corresponding to cell. More...
 
void expandSolution (Vector &result, const Vector &u)
 Expand a system vector to a solution vector. More...
 
void addMPC (MPC *mpc)
 Add a MPC. More...
 
MPCgetMPC (int node, int dof)
 Look for and return a MPC for a specified node+dof. More...
 
void updateFixedNode (int node, const std::pair< Direction, NodeValue > &entry)
 Update/add a fixed node. More...
 
bool isFixed (int node)
 Check if a node is marked as fixed (in any direction) More...
 
void printOperator () const
 Print the current operator. More...
 
void printLoadVector () const
 Print the current load vector. More...
 
AdjacencyPatterngetAdjacencyPattern ()
 Access current adjacency pattern. More...
 

Static Public Attributes

static const int dim = GridType::dimension
 The dimension of the grid. More...
 

Protected Types

typedef std::pair< Direction,
NodeValue
fixEntry
 Fixed nodes. More...
 
typedef std::map< int, fixEntryfixMap
 A mapping from dof to fix value info. More...
 
typedef fixMap::iterator fixIt
 Iterator over a fixmap. More...
 

Protected Member Functions

void resolveMPCChains ()
 Resolve chained MPCs. More...
 
void resolveMPCChain (MPC *mpc)
 Internal function. Handles a single MPC. More...
 
void preprocess ()
 Internal function. Generate meqn for registered MPC/fixed nodes. More...
 
void nodeAdjacency (const LeafIterator &it, int vertexsize, int row)
 Internal function. Generate adjacency pattern for a given node. More...
 
void determineAdjacencyPattern ()
 Internal function. Calculate adjacency pattern. More...
 
template<int esize>
void addDOF (int row, int erow, const Dune::FieldMatrix< double, esize, esize > *K, const Dune::FieldVector< double, esize > *S, const LeafIndexSet &set, const LeafIterator &cell, Vector *b, double scale=1.f)
 Internal function. Assemble entries for a single DOF. More...
 

Protected Attributes

MPCMap mpcs
 The set of MPC. More...
 
std::vector< int > meqn
 Vector of (interleaved) dof<->eqn mapping. More...
 
fixMap fixedNodes
 The map holding information about our fixed nodes. More...
 
AdjacencyPattern adjacencyPattern
 Holds the adjacency pattern of the sparse matrix. More...
 
Matrix A
 The linear operator. More...
 
Vector b
 The load vector. More...
 
const GridType & gv
 A reference to the grid in use. More...
 
size_t maxeqn
 The number of equations in the system. More...
 

Detailed Description

template<class GridType>
class Opm::Elasticity::ASMHandler< GridType >

Class handling finite element assembly.

Member Typedef Documentation

template<class GridType>
typedef std::pair<Direction,NodeValue> Opm::Elasticity::ASMHandler< GridType >::fixEntry
protected

Fixed nodes.

template<class GridType>
typedef fixMap::iterator Opm::Elasticity::ASMHandler< GridType >::fixIt
protected

Iterator over a fixmap.

template<class GridType>
typedef std::map<int,fixEntry> Opm::Elasticity::ASMHandler< GridType >::fixMap
protected

A mapping from dof to fix value info.

template<class GridType>
typedef GridType::LeafGridView::IndexSet Opm::Elasticity::ASMHandler< GridType >::LeafIndexSet

A set of indices.

template<class GridType>
typedef GridType::LeafGridView::template Codim<0>::Iterator Opm::Elasticity::ASMHandler< GridType >::LeafIterator

An iterator over grid cells.

template<class GridType>
typedef Dune::FieldVector<double,dim> Opm::Elasticity::ASMHandler< GridType >::NodeValue

A vectorial node value.

Constructor & Destructor Documentation

template<class GridType>
Opm::Elasticity::ASMHandler< GridType >::ASMHandler ( const GridType &  gv_)
inline

The default constructor.

Parameters
[in]gv_The grid the operator is assembled over
template<class GridType>
Opm::Elasticity::ASMHandler< GridType >::~ASMHandler ( )
inline

Member Function Documentation

template<class GridType >
template<int esize>
void Opm::Elasticity::ASMHandler< GridType >::addDOF ( int  row,
int  erow,
const Dune::FieldMatrix< double, esize, esize > *  K,
const Dune::FieldVector< double, esize > *  S,
const LeafIndexSet set,
const LeafIterator cell,
Vector b,
double  scale = 1.f 
)
protected

Internal function. Assemble entries for a single DOF.

Parameters
[in]rowThe row in the global matrix
[in]erowThe row in the element matrix
[in]KPointer to the element matrix. Can be NULL
[in]SPointer to the element load vector. Can be NULL
[in]setThe index set
[in]cellAn iterator pointing to the cell we're assembling for
[in]bVector to add contributions to
[in]scaleScale for elements. Used with MPC couplings

References Opm::Elasticity::MPC::DOF::coeff, Opm::Elasticity::MPC::DOF::dof, Opm::Elasticity::MPC::getMaster(), Opm::Elasticity::MPC::getNoMaster(), Opm::Elasticity::j, and Opm::Elasticity::MPC::DOF::node.

template<class GridType >
template<int esize>
void Opm::Elasticity::ASMHandler< GridType >::addElement ( const Dune::FieldMatrix< double, esize, esize > *  K,
const Dune::FieldVector< double, esize > *  S,
const LeafIterator cell,
Vector b = NULL 
)

Add an element matrix/vector to the system.

Parameters
[in]KPointer to the element matrix. Can be NULL
[in]SPointer to the element load vector. Can be NULL
[in]cellAn iterator pointing to the cell we're assembling for
[in]bVector to add contributions to. If not given, we use the internal vector

References Opm::Elasticity::MPC::DOF::coeff, Opm::Elasticity::MPC::DOF::dof, Opm::Elasticity::MPC::getMaster(), Opm::Elasticity::MPC::getNoMaster(), Opm::Elasticity::MPC::DOF::node, and Opm::Elasticity::XYZ.

template<class GridType >
void Opm::Elasticity::ASMHandler< GridType >::addMPC ( MPC mpc)

Add a MPC.

Parameters
[in]mpcPointer to the MPC to add.
Note
This class handles destruction

References Opm::Elasticity::MPC::DOF::dof, Opm::Elasticity::MPC::getSlave(), and Opm::Elasticity::MPC::DOF::node.

template<class GridType >
void Opm::Elasticity::ASMHandler< GridType >::determineAdjacencyPattern ( )
protected
template<class GridType >
void Opm::Elasticity::ASMHandler< GridType >::expandSolution ( Vector result,
const Vector u 
)
template<class GridType >
template<int comp>
void Opm::Elasticity::ASMHandler< GridType >::extractValues ( Dune::FieldVector< double, comp > &  v,
const Vector u,
const LeafIterator it 
)

Extract values corresponding to cell.

Parameters
[in]uThe global load vector
[in]itAn iterator to the cell we want to extract values for
[out]vVector holding the values requested

References Opm::Elasticity::MPC::DOF::coeff, Opm::Elasticity::MPC::DOF::dof, Opm::Elasticity::MPC::getMaster(), Opm::Elasticity::MPC::getNoMaster(), and Opm::Elasticity::MPC::DOF::node.

template<class GridType>
AdjacencyPattern& Opm::Elasticity::ASMHandler< GridType >::getAdjacencyPattern ( )
inline

Access current adjacency pattern.

Can be used to add extra entries, such as other blocks

References Opm::Elasticity::ASMHandler< GridType >::adjacencyPattern.

template<class GridType>
size_t Opm::Elasticity::ASMHandler< GridType >::getEqns ( ) const
inline

Get the number of equations in the system.

Returns
The number of equations in the system

References Opm::Elasticity::ASMHandler< GridType >::maxeqn.

template<class GridType>
int Opm::Elasticity::ASMHandler< GridType >::getEquationForDof ( int  node,
int  dof 
)
inline

Get the equation number for a given dof on a given node.

Parameters
[in]nodeThe node number
[in]dofThe DOF
Returns
The equation number if active, -1 otherwise

References Opm::Elasticity::ASMHandler< GridType >::meqn.

template<class GridType>
Vector& Opm::Elasticity::ASMHandler< GridType >::getLoadVector ( )
inline

Obtain a reference to the load vector.

Returns
Reference to load vector

References Opm::Elasticity::ASMHandler< GridType >::b.

template<class GridType >
MPC * Opm::Elasticity::ASMHandler< GridType >::getMPC ( int  node,
int  dof 
)

Look for and return a MPC for a specified node+dof.

Parameters
[in]nodeThe requested node
[in]dofThe requested DOF at given node
Returns
The MPC for the node/dof if found, else NULL
template<class GridType>
Matrix& Opm::Elasticity::ASMHandler< GridType >::getOperator ( )
inline

Obtain a reference to the linear operator.

Returns
Reference to linear operator

References Opm::Elasticity::ASMHandler< GridType >::A.

template<class GridType >
void Opm::Elasticity::ASMHandler< GridType >::initForAssembly ( )

This function needs to be called before starting the element assembly.

References Opm::Elasticity::MatrixOps::fromAdjacency(), Opm::Elasticity::X, Opm::Elasticity::Y, and Opm::Elasticity::Z.

template<class GridType>
bool Opm::Elasticity::ASMHandler< GridType >::isFixed ( int  node)
inline

Check if a node is marked as fixed (in any direction)

Parameters
[in]nodeThe node to query for

References Opm::Elasticity::ASMHandler< GridType >::fixedNodes.

template<class GridType >
void Opm::Elasticity::ASMHandler< GridType >::nodeAdjacency ( const LeafIterator it,
int  vertexsize,
int  row 
)
protected

Internal function. Generate adjacency pattern for a given node.

Parameters
[in]itIterator pointing to the cell in of the node
[in]vertexsizeNumber of vertices in the cell
[in]rowThe equation number/row in matrix

References Opm::Elasticity::MPC::DOF::dof, Opm::Elasticity::MPC::getMaster(), Opm::Elasticity::MPC::getNoMaster(), Opm::Elasticity::j, and Opm::Elasticity::MPC::DOF::node.

template<class GridType >
void Opm::Elasticity::ASMHandler< GridType >::preprocess ( )
protected

Internal function. Generate meqn for registered MPC/fixed nodes.

template<class GridType >
void Opm::Elasticity::ASMHandler< GridType >::printLoadVector ( ) const

Print the current load vector.

template<class GridType >
void Opm::Elasticity::ASMHandler< GridType >::printOperator ( ) const

Print the current operator.

References Opm::Elasticity::MatrixOps::print().

template<class GridType>
void Opm::Elasticity::ASMHandler< GridType >::resolveMPCChains ( )
inlineprotected
template<class GridType >
void Opm::Elasticity::ASMHandler< GridType >::updateFixedNode ( int  node,
const std::pair< Direction, NodeValue > &  entry 
)

Update/add a fixed node.

Parameters
[in]nodeThe node number
[in]entryThe fixed values

Member Data Documentation

template<class GridType>
Matrix Opm::Elasticity::ASMHandler< GridType >::A
protected

The linear operator.

Referenced by Opm::Elasticity::ASMHandler< GridType >::getOperator().

template<class GridType>
AdjacencyPattern Opm::Elasticity::ASMHandler< GridType >::adjacencyPattern
protected

Holds the adjacency pattern of the sparse matrix.

Referenced by Opm::Elasticity::ASMHandler< GridType >::getAdjacencyPattern().

template<class GridType>
Vector Opm::Elasticity::ASMHandler< GridType >::b
protected
template<class GridType>
const int Opm::Elasticity::ASMHandler< GridType >::dim = GridType::dimension
static

The dimension of the grid.

template<class GridType>
fixMap Opm::Elasticity::ASMHandler< GridType >::fixedNodes
protected

The map holding information about our fixed nodes.

Referenced by Opm::Elasticity::ASMHandler< GridType >::isFixed().

template<class GridType>
const GridType& Opm::Elasticity::ASMHandler< GridType >::gv
protected

A reference to the grid in use.

template<class GridType>
size_t Opm::Elasticity::ASMHandler< GridType >::maxeqn
protected

The number of equations in the system.

Referenced by Opm::Elasticity::ASMHandler< GridType >::getEqns().

template<class GridType>
std::vector<int> Opm::Elasticity::ASMHandler< GridType >::meqn
protected

Vector of (interleaved) dof<->eqn mapping.

Referenced by Opm::Elasticity::ASMHandler< GridType >::getEquationForDof().


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