asmhandler.hpp
Go to the documentation of this file.
1 //==============================================================================
11 //==============================================================================
12 #ifndef ASMHANDLER_HPP_
13 #define ASMHANDLER_HPP_
14 
15 #include <dune/geometry/referenceelements.hh>
16 #include <dune/common/fmatrix.hh>
17 #include <dune/istl/bvector.hh>
18 #include <dune/common/fvector.hh>
20 #include <opm/elasticity/mpc.hh>
22 
23 namespace Opm {
24 namespace Elasticity {
25 
26 
28  template<class GridType>
29 class ASMHandler {
30  public:
32  static const int dim = GridType::dimension;
33 
35  typedef typename GridType::LeafGridView::IndexSet LeafIndexSet;
36 
38  typedef typename GridType::LeafGridView::template Codim<0>::Iterator LeafIterator;
39 
42  ASMHandler(const GridType& gv_) : gv(gv_), maxeqn(0)
43  {
44  }
45 
48  {
49  for (MPCMap::iterator it=mpcs.begin(); it != mpcs.end();++it)
50  delete it->second;
51  }
52 
54  typedef Dune::FieldVector<double,dim> NodeValue;
55 
58  size_t getEqns() const
59  {
60  return maxeqn;
61  }
62 
67  int getEquationForDof(int node, int dof)
68  {
69  return meqn[node*dim+dof];
70  }
71 
75  {
76  return A;
77  }
78 
82  {
83  return b;
84  }
85 
88  void initForAssembly();
89 
96  template<int esize>
97  void addElement(const Dune::FieldMatrix<double,esize,esize>* K,
98  const Dune::FieldVector<double,esize>* S,
99  const LeafIterator& cell,
100  Vector* b=NULL);
101 
106  template<int comp>
107  void extractValues(Dune::FieldVector<double,comp>& v,
108  const Vector& u, const LeafIterator& it);
109 
111  void expandSolution(Vector& result, const Vector& u);
112 
116  void addMPC(MPC* mpc);
117 
122  MPC* getMPC(int node, int dof);
123 
127  void updateFixedNode(int node,
128  const std::pair<Direction,NodeValue>& entry);
129 
132  bool isFixed(int node)
133  {
134  return (fixedNodes.find(node) != fixedNodes.end());
135  }
136 
138  void printOperator() const;
139 
141  void printLoadVector() const;
142 
146  {
147  return adjacencyPattern;
148  }
149  protected:
152  {
153  for (MPCMap::iterator it = mpcs.begin();
154  it != mpcs.end();++it)
155  resolveMPCChain(it->second);
156  }
157 
160  void resolveMPCChain(MPC* mpc);
161 
163  void preprocess();
164 
169  void nodeAdjacency(const LeafIterator& it, int vertexsize, int row);
170 
173 
183  template<int esize>
184  void addDOF(int row, int erow,
185  const Dune::FieldMatrix<double,esize,esize>* K,
186  const Dune::FieldVector<double,esize>* S,
187  const LeafIndexSet& set,
188  const LeafIterator& cell,
189  Vector* b,
190  double scale=1.f);
191 
194 
196  std::vector<int> meqn;
197 
199  typedef std::pair<Direction,NodeValue> fixEntry;
200 
202  typedef std::map<int,fixEntry> fixMap;
203 
205  typedef typename fixMap::iterator fixIt;
206 
208  fixMap fixedNodes;
209 
212 
215 
218 
220  const GridType& gv;
221 
223  size_t maxeqn;
224  private:
226  ASMHandler(const ASMHandler&) {}
228  ASMHandler& operator=(const ASMHandler&) {}
229 };
230 
231 }
232 }
233 
234 #include "asmhandler_impl.hpp"
235 
236 #endif
size_t getEqns() const
Get the number of equations in the system.
Definition: asmhandler.hpp:58
fixMap fixedNodes
The map holding information about our fixed nodes.
Definition: asmhandler.hpp:208
size_t maxeqn
The number of equations in the system.
Definition: asmhandler.hpp:223
Class handling finite element assembly.
Definition: asmhandler.hpp:29
void printLoadVector() const
Print the current load vector.
Definition: asmhandler_impl.hpp:243
~ASMHandler()
Destructor.
Definition: asmhandler.hpp:47
Definition: applier.hpp:18
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.
Definition: asmhandler_impl.hpp:52
Vector b
The load vector.
Definition: asmhandler.hpp:217
Matrix A
The linear operator.
Definition: asmhandler.hpp:214
void preprocess()
Internal function. Generate meqn for registered MPC/fixed nodes.
Definition: asmhandler_impl.hpp:317
GridType::LeafGridView::template Codim< 0 >::Iterator LeafIterator
An iterator over grid cells.
Definition: asmhandler.hpp:38
const GridType & gv
A reference to the grid in use.
Definition: asmhandler.hpp:220
std::vector< std::set< int > > AdjacencyPattern
For storing matrix adjacency/sparsity patterns.
Definition: matrixops.hpp:26
void extractValues(Dune::FieldVector< double, comp > &v, const Vector &u, const LeafIterator &it)
Extract values corresponding to cell.
Definition: asmhandler_impl.hpp:114
static const int dim
The dimension of the grid.
Definition: asmhandler.hpp:32
std::pair< Direction, NodeValue > fixEntry
Fixed nodes.
Definition: asmhandler.hpp:199
Vector & getLoadVector()
Obtain a reference to the load vector.
Definition: asmhandler.hpp:81
A class for representing a general multi-point constraint equation.
Definition: mpc.hh:58
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.
Definition: asmhandler_impl.hpp:85
ASMHandler(const GridType &gv_)
The default constructor.
Definition: asmhandler.hpp:42
Dune::FieldVector< double, dim > NodeValue
A vectorial node value.
Definition: asmhandler.hpp:54
std::map< int, MPC * > MPCMap
A mapping from dof to MPCs.
Definition: mpc.hh:159
Dune::BCRSMatrix< Dune::FieldMatrix< double, 1, 1 > > Matrix
A sparse matrix holding our operator.
Definition: matrixops.hpp:23
void expandSolution(Vector &result, const Vector &u)
Expand a system vector to a solution vector.
Definition: asmhandler_impl.hpp:147
MPCMap mpcs
The set of MPC.
Definition: asmhandler.hpp:193
Helper class with some matrix operations.
void resolveMPCChain(MPC *mpc)
Internal function. Handles a single MPC.
Definition: asmhandler_impl.hpp:252
bool isFixed(int node)
Check if a node is marked as fixed (in any direction)
Definition: asmhandler.hpp:132
void updateFixedNode(int node, const std::pair< Direction, NodeValue > &entry)
Update/add a fixed node.
Definition: asmhandler_impl.hpp:216
Logging helper utilities.
void determineAdjacencyPattern()
Internal function. Calculate adjacency pattern.
Definition: asmhandler_impl.hpp:377
AdjacencyPattern adjacencyPattern
Holds the adjacency pattern of the sparse matrix.
Definition: asmhandler.hpp:211
void printOperator() const
Print the current operator.
Definition: asmhandler_impl.hpp:237
void nodeAdjacency(const LeafIterator &it, int vertexsize, int row)
Internal function. Generate adjacency pattern for a given node.
Definition: asmhandler_impl.hpp:353
Matrix & getOperator()
Obtain a reference to the linear operator.
Definition: asmhandler.hpp:74
Representation of multi-point constraint (MPC) equations.
void initForAssembly()
This function needs to be called before starting the element assembly.
Definition: asmhandler_impl.hpp:22
std::vector< int > meqn
Vector of (interleaved) dof<->eqn mapping.
Definition: asmhandler.hpp:196
std::map< int, fixEntry > fixMap
A mapping from dof to fix value info.
Definition: asmhandler.hpp:202
void addMPC(MPC *mpc)
Add a MPC.
Definition: asmhandler_impl.hpp:189
MPC * getMPC(int node, int dof)
Look for and return a MPC for a specified node+dof.
Definition: asmhandler_impl.hpp:207
Class handling finite element assembly - template implementations.
void resolveMPCChains()
Resolve chained MPCs.
Definition: asmhandler.hpp:151
Dune::BlockVector< Dune::FieldVector< double, 1 > > Vector
A vector holding our RHS.
Definition: matrixops.hpp:29
GridType::LeafGridView::IndexSet LeafIndexSet
A set of indices.
Definition: asmhandler.hpp:35
int getEquationForDof(int node, int dof)
Get the equation number for a given dof on a given node.
Definition: asmhandler.hpp:67
AdjacencyPattern & getAdjacencyPattern()
Access current adjacency pattern.
Definition: asmhandler.hpp:145
fixMap::iterator fixIt
Iterator over a fixmap.
Definition: asmhandler.hpp:205