24#ifndef DUMUX_DISCRETIZATION_CC_MPFA_HELPER_HH
25#define DUMUX_DISCRETIZATION_CC_MPFA_HELPER_HH
27#include <dune/common/exceptions.hh>
28#include <dune/common/fvector.hh>
29#include <dune/common/fmatrix.hh>
30#include <dune/common/parallel/mpihelper.hh>
31#include <dune/geometry/type.hh>
32#include <dune/geometry/referenceelements.hh>
42template<
class Gr
idGeometry,
int dim,
int dimWorld>
49template<
class Gr
idGeometry>
52 using GridView =
typename GridGeometry::GridView;
53 using CoordScalar =
typename GridView::ctype;
54 using Element =
typename GridView::template Codim<0>::Entity;
56 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
58 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
59 using ScvfCornerVector =
typename SubControlVolumeFace::Traits::CornerStorage;
63 using ScvfPositionsOnIntersection = std::array<GlobalPosition, 3>;
70 template<
class ScvBasis >
73 static const Dune::FieldMatrix<CoordScalar, 2, 2> R = {{0.0, 1.0}, {-1.0, 0.0}};
75 ScvBasis innerNormals;
76 R.mv(scvBasis[1], innerNormals[0]);
77 R.mv(scvBasis[0], innerNormals[1]);
80 if (!isRightHandSystem(scvBasis))
81 innerNormals[0] *= -1.0;
83 innerNormals[1] *= -1.0;
93 template<
class ScvBasis >
97 return abs(crossProduct<CoordScalar>(scvBasis[0], scvBasis[1]));
107 assert(gridView.size(Dune::GeometryTypes::triangle)
108 + gridView.size(Dune::GeometryTypes::quadrilateral) == gridView.size(0));
110 return gridView.size(Dune::GeometryTypes::triangle)
111 * getNumLocalScvfs(Dune::GeometryTypes::triangle)
112 + gridView.size(Dune::GeometryTypes::quadrilateral)
113 * getNumLocalScvfs(Dune::GeometryTypes::quadrilateral);
120 template<
class ScvBasis >
122 {
return !std::signbit(crossProduct<CoordScalar>(scvBasis[0], scvBasis[1])); }
134 template<
class ElementGeometry,
class ReferenceElement>
136 const ReferenceElement& refElement,
137 unsigned int indexInElement,
138 unsigned int numCorners)
140 ScvfPositionsOnIntersection p;
144 for (
unsigned int c = 0; c < numCorners; ++c)
147 p[c+1] = eg.global(refElement.position(refElement.subEntity(indexInElement, 1, c, 2), 2));
164 unsigned int numIntersectionCorners,
165 unsigned int cornerIdx)
168 assert(cornerIdx < 2 &&
"provided index exceeds the number of corners of facets in 2d");
171 return cornerIdx == 0 ? ScvfCornerVector({p[0], p[1]})
172 : ScvfCornerVector({p[0], p[2]});
180 {
return (scvfCorners[1]-scvfCorners[0]).two_norm(); }
188 if (gt == Dune::GeometryTypes::triangle)
190 else if (gt == Dune::GeometryTypes::quadrilateral)
193 DUNE_THROW(Dune::NotImplemented,
"Mpfa for 2d geometry type " << gt);
202template<
class Gr
idGeometry>
206 using GridView =
typename GridGeometry::GridView;
207 using CoordScalar =
typename GridView::ctype;
214 template<
class ScvBasis >
218 const auto normal = [&scvBasis] ()
220 auto n = crossProduct<CoordScalar>(scvBasis[0], scvBasis[1]);
226 ScvBasis innerNormals;
227 innerNormals[0] = crossProduct<CoordScalar>(scvBasis[1], normal);
228 innerNormals[1] = crossProduct<CoordScalar>(normal, scvBasis[0]);
241 template<
class ScvBasis >
245 return abs(crossProduct<CoordScalar>(scvBasis[0], scvBasis[1]).two_norm());
255 template<
class ScvBasis >
263template<
class Gr
idGeometry>
266 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
267 using ScvfCornerVector =
typename SubControlVolumeFace::Traits::CornerStorage;
270 using GridView =
typename GridGeometry::GridView;
271 using Element =
typename GridView::template Codim<0>::Entity;
272 using CoordScalar =
typename GridView::ctype;
273 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
278 using ScvfPositionsOnIntersection = std::array<GlobalPosition, 9>;
287 template<
class ScvBasis >
290 ScvBasis innerNormals;
292 innerNormals[0] = crossProduct<CoordScalar>(scvBasis[1], scvBasis[2]);
293 innerNormals[1] = crossProduct<CoordScalar>(scvBasis[2], scvBasis[0]);
294 innerNormals[2] = crossProduct<CoordScalar>(scvBasis[0], scvBasis[1]);
296 if (!isRightHandSystem(scvBasis))
297 std::for_each(innerNormals.begin(), innerNormals.end(), [] (
auto& n) { n *= -1.0; });
308 template<
class ScvBasis >
312 return abs(tripleProduct<CoordScalar>(scvBasis[0], scvBasis[1], scvBasis[2]));
323 assert(gridView.size(Dune::GeometryTypes::tetrahedron)
324 + gridView.size(Dune::GeometryTypes::pyramid)
325 + gridView.size(Dune::GeometryTypes::prism)
326 + gridView.size(Dune::GeometryTypes::hexahedron) == gridView.size(0));
328 return gridView.size(Dune::GeometryTypes::tetrahedron)
329 * getNumLocalScvfs(Dune::GeometryTypes::tetrahedron)
330 + gridView.size(Dune::GeometryTypes::pyramid)
331 * getNumLocalScvfs(Dune::GeometryTypes::pyramid)
332 + gridView.size(Dune::GeometryTypes::prism)
333 * getNumLocalScvfs(Dune::GeometryTypes::prism)
334 + gridView.size(Dune::GeometryTypes::hexahedron)
335 * getNumLocalScvfs(Dune::GeometryTypes::hexahedron);
342 template<
class ScvBasis >
344 {
return !std::signbit(tripleProduct<CoordScalar>(scvBasis[0], scvBasis[1], scvBasis[2])); }
356 template<
class ElementGeometry,
class ReferenceElement>
358 const ReferenceElement& refElement,
359 unsigned int indexInElement,
360 unsigned int numCorners)
363 ScvfPositionsOnIntersection p;
365 DUNE_THROW(Dune::InvalidStateException,
"Mpfa implementation cannot handle faces with more than 4 corners");
369 for (
unsigned int c = 0; c < numCorners; ++c)
372 p[c+1] = eg.global(refElement.position(refElement.subEntity(indexInElement, 1, c, 3), 3));
383 p[numCorners+1] = p[2] + p[1];
384 p[numCorners+1] /= 2;
385 p[numCorners+2] = p[3] + p[1];
386 p[numCorners+2] /= 2;
387 p[numCorners+3] = p[3] + p[2];
388 p[numCorners+3] /= 2;
394 p[numCorners+1] = p[3] + p[1];
395 p[numCorners+1] /= 2;
396 p[numCorners+2] = p[4] + p[2];
397 p[numCorners+2] /= 2;
398 p[numCorners+3] = p[2] + p[1];
399 p[numCorners+3] /= 2;
400 p[numCorners+4] = p[4] + p[3];
401 p[numCorners+4] /= 2;
405 DUNE_THROW(Dune::NotImplemented,
406 "Mpfa scvf corners for dim = 3, dimWorld = 3, corners = " << numCorners);
419 unsigned int numIntersectionCorners,
420 unsigned int cornerIdx)
424 switch (numIntersectionCorners)
429 static const std::uint8_t vo = 1;
430 static const std::uint8_t eo = 4;
431 static const std::uint8_t map[3][4] =
433 {0, eo+1, eo+0, vo+0},
434 {0, eo+0, eo+2, vo+1},
435 {0, eo+2, eo+1, vo+2}
438 return ScvfCornerVector( {p[map[cornerIdx][0]],
439 p[map[cornerIdx][1]],
440 p[map[cornerIdx][2]],
441 p[map[cornerIdx][3]]} );
446 static const std::uint8_t vo = 1;
447 static const std::uint8_t eo = 5;
448 static const std::uint8_t map[4][4] =
450 {0, eo+0, eo+2, vo+0},
451 {0, eo+2, eo+1, vo+1},
452 {0, eo+3, eo+0, vo+2},
453 {0, eo+1, eo+3, vo+3}
456 return ScvfCornerVector( {p[map[cornerIdx][0]],
457 p[map[cornerIdx][1]],
458 p[map[cornerIdx][2]],
459 p[map[cornerIdx][3]]} );
462 DUNE_THROW(Dune::NotImplemented,
463 "Mpfa scvf corners for dim = 3, dimWorld = 3, corners = " << numIntersectionCorners);
474 return 0.5*
Dumux::crossProduct(scvfCorners[3]-scvfCorners[0], scvfCorners[2]-scvfCorners[1]).two_norm();
484 if (gt == Dune::GeometryTypes::tetrahedron)
486 else if (gt == Dune::GeometryTypes::pyramid)
488 else if (gt == Dune::GeometryTypes::prism)
490 else if (gt == Dune::GeometryTypes::hexahedron)
493 DUNE_THROW(Dune::NotImplemented,
"Mpfa for 3d geometry type " << gt);
503template<
class Gr
idGeometry>
505 GridGeometry::GridView::dimension,
506 GridGeometry::GridView::dimensionworld>
508 using PrimaryInteractionVolume =
typename GridGeometry::GridIVIndexSets::PrimaryInteractionVolume;
509 using SecondaryInteractionVolume =
typename GridGeometry::GridIVIndexSets::SecondaryInteractionVolume;
511 using VertexMapper =
typename GridGeometry::VertexMapper;
512 using FVElementGeometry =
typename GridGeometry::LocalView;
513 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
514 using ScvfCornerVector =
typename SubControlVolumeFace::Traits::CornerStorage;
516 using GridView =
typename GridGeometry::GridView;
517 static constexpr int dim = GridView::dimension;
519 using Element =
typename GridView::template Codim<0>::Entity;
520 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
521 using CoordScalar =
typename GridView::ctype;
522 using ReferenceElements =
typename Dune::ReferenceElements<CoordScalar, dim>;
531 template<
class Scalar >
536 return scvfCorners[0];
538 auto ip = scvfCorners.back() - scvfCorners.front();
540 ip += scvfCorners[0];
546 static std::vector<bool>
findGhostVertices(
const GridView& gridView,
const VertexMapper& vertexMapper)
548 std::vector<bool> ghostVertices(gridView.size(dim),
false);
551 if (Dune::MPIHelper::getCollectiveCommunication().size() == 1)
552 return ghostVertices;
555 if (gridView.ghostSize(0) > 0)
556 DUNE_THROW(Dune::InvalidStateException,
"Mpfa methods in parallel do not work with ghost cells. Use overlap cells instead.");
559 if (gridView.overlapSize(0) == 0)
560 DUNE_THROW(Dune::InvalidStateException,
"Grid no overlaping cells. This is required by mpfa methods in parallel.");
562 for (
const auto& element : elements(gridView))
566 if (!is.neighbor() && !is.boundary())
568 const auto refElement = ReferenceElements::general(element.type());
570 for (
int isVertex = 0; isVertex < is.geometry().corners(); ++isVertex)
572 const auto vIdxLocal = refElement.subEntity(is.indexInInside(), 1, isVertex, dim);
573 const auto vIdxGlobal = vertexMapper.subIndex(element, vIdxLocal, dim);
574 ghostVertices[vIdxGlobal] =
true;
580 return ghostVertices;
586 {
return !std::is_same<PrimaryInteractionVolume, SecondaryInteractionVolume>::value; }
589 template<
typename V,
typename T>
591 {
return std::find(vector.begin(), vector.end(), value) != vector.end(); }
Define some often used mathematical functions.
Dune::FieldVector< Scalar, 3 > crossProduct(const Dune::FieldVector< Scalar, 3 > &vec1, const Dune::FieldVector< Scalar, 3 > &vec2)
Cross product of two vectors in three-dimensional Euclidean space.
Definition: math.hh:631
Dune::IteratorRange< typename MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper >::Intersections::const_iterator > intersections(const MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper > &glue)
Range generator to iterate with range-based for loops over all intersections as follows: for (const a...
Definition: glue.hh:62
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Dimension-specific helper class to get data required for mpfa scheme.
Definition: helper.hh:43
static std::size_t getNumLocalScvfs(const Dune::GeometryType >)
Calculates the number of scvfs in a given element geometry type.
Definition: helper.hh:186
static ScvfPositionsOnIntersection computeScvfCornersOnIntersection(const ElementGeometry &eg, const ReferenceElement &refElement, unsigned int indexInElement, unsigned int numCorners)
Returns a vector containing the positions on a given intersection that are relevant for scvf corner c...
Definition: helper.hh:135
static ScvBasis calculateInnerNormals(const ScvBasis &scvBasis)
Calculates the inner normal vectors to a given scv basis.
Definition: helper.hh:71
static CoordScalar computeScvfArea(const ScvfCornerVector &scvfCorners)
Calculates the area of an scvf.
Definition: helper.hh:179
static ScvfCornerVector getScvfCorners(const ScvfPositionsOnIntersection &p, unsigned int numIntersectionCorners, unsigned int cornerIdx)
Returns the corners of the sub control volume face constructed in a corner (vertex) of an intersectio...
Definition: helper.hh:163
static CoordScalar calculateDetX(const ScvBasis &scvBasis)
Calculates the determinant of an scv basis. This is equal to the cross product for dim = dimWorld = 2...
Definition: helper.hh:94
static std::size_t getGlobalNumScvf(const GridView &gridView)
Returns the global number of scvfs in the grid. Assumes the grid to be made up of only basic geometry...
Definition: helper.hh:105
static bool isRightHandSystem(const ScvBasis &scvBasis)
Checks whether or not a given scv basis forms a right hand system.
Definition: helper.hh:121
static constexpr bool isRightHandSystem(const ScvBasis &scvBasis)
Checks whether or not a given scv basis forms a right hand system. Note that for dim = 2 < dimWorld =...
Definition: helper.hh:256
static ScvBasis calculateInnerNormals(const ScvBasis &scvBasis)
Calculates the inner normal vectors to a given scv basis.
Definition: helper.hh:215
static CoordScalar calculateDetX(const ScvBasis &scvBasis)
Calculates the determinant of an scv basis. For dim = 2 < dimWorld = 3 this is actually not the deter...
Definition: helper.hh:242
static CoordScalar calculateDetX(const ScvBasis &scvBasis)
Calculates the determinant of an scv basis. This is equal to the cross product for dim = dimWorld = 2...
Definition: helper.hh:309
static ScvfPositionsOnIntersection computeScvfCornersOnIntersection(const ElementGeometry &eg, const ReferenceElement &refElement, unsigned int indexInElement, unsigned int numCorners)
Returns a vector containing the positions on a given intersection that are relevant for scvf corner c...
Definition: helper.hh:357
static ScvBasis calculateInnerNormals(const ScvBasis &scvBasis)
Calculates the inner normal vectors to a given scv basis.
Definition: helper.hh:288
static bool isRightHandSystem(const ScvBasis &scvBasis)
Checks whether or not a given scv basis forms a right hand system.
Definition: helper.hh:343
static std::size_t getGlobalNumScvf(const GridView &gridView)
Returns the global number of scvfs in the grid. Assumes the grid to be made up of only basic geometry...
Definition: helper.hh:321
static ScvfCornerVector getScvfCorners(const ScvfPositionsOnIntersection &p, unsigned int numIntersectionCorners, unsigned int cornerIdx)
Returns the corners of the sub control volume face constructed in a corner (vertex) of an intersectio...
Definition: helper.hh:418
static std::size_t getNumLocalScvfs(const Dune::GeometryType >)
Calculates the number of scvfs in a given element geometry type.
Definition: helper.hh:482
static CoordScalar computeScvfArea(const ScvfCornerVector &scvfCorners)
Calculates the area of an scvf.
Definition: helper.hh:471
Helper class to get the required information on an interaction volume.
Definition: helper.hh:507
static constexpr bool considerSecondaryIVs()
Definition: helper.hh:585
static bool vectorContainsValue(const V &vector, const T value)
returns whether or not a value exists in a vector
Definition: helper.hh:590
static std::vector< bool > findGhostVertices(const GridView &gridView, const VertexMapper &vertexMapper)
Definition: helper.hh:546
static GlobalPosition getScvfIntegrationPoint(const ScvfCornerVector &scvfCorners, Scalar q)
Calculates the integration point on an scvf.
Definition: helper.hh:532