24#ifndef DUMUX_CROPERATOR2PADAPTIVE_HH
25#define DUMUX_CROPERATOR2PADAPTIVE_HH
33#include<dune/common/timer.hh>
34#include<dune/common/fvector.hh>
35#include<dune/common/fmatrix.hh>
36#include<dune/common/exceptions.hh>
37#include<dune/geometry/type.hh>
38#include<dune/grid/common/grid.hh>
39#include<dune/grid/common/mcmgmapper.hh>
41#include<dune/istl/bvector.hh>
42#include<dune/istl/operators.hh>
43#include<dune/istl/bcrsmatrix.hh>
67template<
class TypeTag>
72 enum {dim=GridView::dimension};
73 using Element =
typename GridView::template Codim<0>::Entity;
74 using IS =
typename GridView::IndexSet;
75 using BlockType = Dune::FieldMatrix<Scalar, 1, 1>;
76 using MatrixType = Dune::BCRSMatrix<BlockType>;
77 using MBlockType =
typename MatrixType::block_type;
78 using rowiterator =
typename MatrixType::RowIterator;
79 using coliterator =
typename MatrixType::ColIterator;
80 using BCBlockType = std::array<BoundaryConditions::Flags, 1>;
81 using SatType = Dune::BlockVector< Dune::FieldVector<double, 1> >;
84 using Indices =
typename GET_PROP_TYPE(TypeTag, ModelTraits)::Indices;
88 pressureEqIdx = Indices::pressureEqIdx,
103 A_.setBuildMode(MatrixType::random);
114 A_.setSize(size(), size());
157 template <
class LocalStiffness,
class Vector>
169template<
class TypeTag>
173 for (
int i = 0; i < size(); i++)
177 std::vector<bool> visited(size(),
false);
179 int numElem = gridView_.size(0);
181 for (
int elemIdx = 0; elemIdx < numElem; elemIdx++)
183 int numFaces = intersectionMapper_.size(elemIdx);
184 for (
int fIdx = 0; fIdx < numFaces; fIdx++)
186 int faceIdxGlobal = intersectionMapper_.subIndex(elemIdx, fIdx);
187 if (!visited[faceIdxGlobal])
189 A_.incrementrowsize(faceIdxGlobal);
190 visited[faceIdxGlobal] =
true;
192 for (
int k = 0; k < numFaces-1; k++) {
193 A_.incrementrowsize(faceIdxGlobal);
203 for (
int i = 0; i < size(); i++)
206 for (
int elemIdx = 0; elemIdx < numElem; elemIdx++)
208 int numFaces = intersectionMapper_.size(elemIdx);
209 for (
int fIdx = 0; fIdx < numFaces; fIdx++)
211 int faceIdxGlobalI = intersectionMapper_.subIndex(elemIdx, fIdx);
212 if (!visited[faceIdxGlobalI])
214 A_.addindex(faceIdxGlobalI,faceIdxGlobalI);
215 visited[faceIdxGlobalI] =
true;
217 for (
int k = 0; k < numFaces; k++)
219 int faceIdxGlobalJ = intersectionMapper_.subIndex(elemIdx, k);
220 A_.addindex(faceIdxGlobalI, faceIdxGlobalJ);
230template<
class TypeTag>
231template <
class LocalStiffness,
class Vector>
236 if (u.N()!=A_.M() || f.N()!=A_.N())
237 DUNE_THROW(Dune::MathError,
"CROperatorAssemblerTwoPAdaptive::assemble(): size mismatch");
243 std::vector<BCBlockType> essential(intersectionMapper_.size());
244 for (
typename std::vector<BCBlockType>::size_type i=0; i<essential.size(); i++)
248 std::vector<int> local2Global(2*dim);
251 for (
const auto& element : elements(gridView_))
257 int eIdxGlobal = intersectionMapper_.index(element);
259 unsigned int numFaces = intersectionMapper_.size(eIdxGlobal);
260 local2Global.resize(numFaces);
262 for (
unsigned int i = 0; i < numFaces; i++)
264 int idx = intersectionMapper_.subIndex(element, i);
265 local2Global[i] = idx;
269 for (
unsigned int i=0; i<numFaces; i++)
272 for (
unsigned int j=0; j<numFaces; j++)
275 A_[local2Global[i]][local2Global[j]] += loc.
mat(i,j);
282 f[local2Global[i]][0] = loc.
rhs(i)[0];
285 f[local2Global[i]][0] += loc.
rhs(i)[0];
290 for (
const auto& element : elements(gridView_))
292 int eIdxGlobal = intersectionMapper_.index(element);
294 unsigned int numFaces = intersectionMapper_.size(eIdxGlobal);
295 local2Global.resize(numFaces);
297 for (
unsigned int i = 0; i < numFaces; i++)
299 int idx = intersectionMapper_.subIndex(element, i);
300 local2Global[i] = idx;
303 loc.completeRHS(element, local2Global, f);
307 rowiterator endi=A_.end();
308 for (rowiterator i=A_.begin(); i!=endi; ++i)
311 if ((
int) i.index() >= (
int) size())
313 coliterator endj=(*i).end();
314 for (coliterator j=(*i).begin(); j!=endj; ++j)
317 if (j.index()==i.index())
327 coliterator endj=(*i).end();
328 for (coliterator j=(*i).begin(); j!=endj; ++j)
329 if (j.index()==i.index())
337 u[i.index()][0] = f[i.index()][0];
#define GET_PROP_TYPE(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:283
Definition of boundary condition types, extend if necessary.
defines intersection mappers.
Base class for assembling local stiffness matrices.
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
@ dirichlet
Dirichlet boundary.
Definition: boundaryconditions.hh:41
@ neumann
Neumann boundary.
Definition: boundaryconditions.hh:39
bool isDirichlet(unsigned eqIdx) const
Returns true if an equation is used to specify a Dirichlet condition.
Definition: common/boundarytypes.hh:236
defines an intersection mapper for mapping of global DOFs assigned to faces which also works for adap...
Definition: intersectionmapper.hh:185
unsigned int size() const
Definition: intersectionmapper.hh:292
void update()
Definition: intersectionmapper.hh:307
Extends CROperatorBase by a generic methods to assemble global stiffness matrix from local stiffness ...
Definition: croperatoradaptive.hh:69
void updateMatrix()
Definition: croperatoradaptive.hh:170
const IntersectionMapper & intersectionMapper()
Definition: croperatoradaptive.hh:130
MatrixType RepresentationType
Definition: croperatoradaptive.hh:98
CROperatorAssemblerTwoPAdaptive(const GridView &gridview)
Definition: croperatoradaptive.hh:100
void assemble(LocalStiffness &loc, Vector &u, Vector &f)
Assembles global stiffness matrix.
Definition: croperatoradaptive.hh:232
void adapt()
Definition: croperatoradaptive.hh:111
void initialize()
Definition: croperatoradaptive.hh:106
const IS & is_
Definition: croperatoradaptive.hh:164
const RepresentationType & operator*() const
Returns const reference to operator matrix.
Definition: croperatoradaptive.hh:119
const IntersectionMapper & intersectionMapper() const
Definition: croperatoradaptive.hh:135
const GridView gridView_
Definition: croperatoradaptive.hh:163
IntersectionMapper intersectionMapper_
Definition: croperatoradaptive.hh:165
RepresentationType A_
Definition: croperatoradaptive.hh:166
Base class for local assemblers.
Definition: localstiffness.hh:60
const VBlockType & rhs(int i) const
Accesses right hand side.
Definition: localstiffness.hh:210
const MBlockType & mat(int i, int j) const
Accesses the local stiffness matrix.
Definition: localstiffness.hh:197
const BoundaryTypes & bc(int i) const
Accesses boundary condition for each dof.
Definition: localstiffness.hh:223
virtual void assemble(const Entity &e, int k=1)=0
Assembles local stiffness matrix including boundary conditions for given element and order.
Base file for properties related to sequential IMPET algorithms.