12#ifndef DUMUX_DISCRETIZATION_CC_MPFA_HELPER_HH
13#define DUMUX_DISCRETIZATION_CC_MPFA_HELPER_HH
15#include <dune/common/exceptions.hh>
16#include <dune/common/fvector.hh>
17#include <dune/common/fmatrix.hh>
18#include <dune/common/parallel/mpihelper.hh>
19#include <dune/geometry/type.hh>
29template<
class Gr
idGeometry,
int dim,
int dimWorld>
36template<
class Gr
idGeometry>
39 using GridView =
typename GridGeometry::GridView;
40 using CoordScalar =
typename GridView::ctype;
41 using Element =
typename GridView::template Codim<0>::Entity;
43 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
45 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
46 using ScvfCornerVector =
typename SubControlVolumeFace::Traits::CornerStorage;
50 using ScvfPositionsOnIntersection = std::array<GlobalPosition, 3>;
57 template<
class ScvBasis >
60 static const Dune::FieldMatrix<CoordScalar, 2, 2> R = {{0.0, 1.0}, {-1.0, 0.0}};
62 ScvBasis innerNormals;
63 R.mv(scvBasis[1], innerNormals[0]);
64 R.mv(scvBasis[0], innerNormals[1]);
67 if (!isRightHandSystem(scvBasis))
68 innerNormals[0] *= -1.0;
70 innerNormals[1] *= -1.0;
80 template<
class ScvBasis >
84 return abs(crossProduct<CoordScalar>(scvBasis[0], scvBasis[1]));
94 assert(gridView.size(Dune::GeometryTypes::triangle)
95 + gridView.size(Dune::GeometryTypes::quadrilateral) == gridView.size(0));
97 return gridView.size(Dune::GeometryTypes::triangle)
98 * getNumLocalScvfs(Dune::GeometryTypes::triangle)
99 + gridView.size(Dune::GeometryTypes::quadrilateral)
100 * getNumLocalScvfs(Dune::GeometryTypes::quadrilateral);
107 template<
class ScvBasis >
109 {
return !std::signbit(crossProduct<CoordScalar>(scvBasis[0], scvBasis[1])); }
121 template<
class ElementGeometry,
class ReferenceElement>
123 const ReferenceElement& refElement,
124 unsigned int indexInElement,
127 ScvfPositionsOnIntersection p;
134 p[c+1] = eg.global(refElement.position(refElement.subEntity(indexInElement, 1, c, 2), 2));
151 unsigned int numIntersectionCorners,
152 unsigned int cornerIdx)
155 assert(cornerIdx < 2 &&
"provided index exceeds the number of corners of facets in 2d");
158 return cornerIdx == 0 ? ScvfCornerVector({p[0], p[1]})
159 : ScvfCornerVector({p[0], p[2]});
167 {
return (scvfCorners[1]-scvfCorners[0]).two_norm(); }
175 if (gt == Dune::GeometryTypes::triangle)
177 else if (gt == Dune::GeometryTypes::quadrilateral)
180 DUNE_THROW(Dune::NotImplemented,
"Mpfa for 2d geometry type " << gt);
189template<
class Gr
idGeometry>
193 using GridView =
typename GridGeometry::GridView;
194 using CoordScalar =
typename GridView::ctype;
201 template<
class ScvBasis >
205 const auto normal = [&scvBasis] ()
207 auto n = crossProduct<CoordScalar>(scvBasis[0], scvBasis[1]);
213 ScvBasis innerNormals;
214 innerNormals[0] = crossProduct<CoordScalar>(scvBasis[1],
normal);
215 innerNormals[1] = crossProduct<CoordScalar>(
normal, scvBasis[0]);
228 template<
class ScvBasis >
232 return abs(crossProduct<CoordScalar>(scvBasis[0], scvBasis[1]).two_norm());
242 template<
class ScvBasis >
250template<
class Gr
idGeometry>
253 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
254 using ScvfCornerVector =
typename SubControlVolumeFace::Traits::CornerStorage;
257 using GridView =
typename GridGeometry::GridView;
258 using Element =
typename GridView::template Codim<0>::Entity;
259 using CoordScalar =
typename GridView::ctype;
260 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
265 using ScvfPositionsOnIntersection = std::array<GlobalPosition, 9>;
274 template<
class ScvBasis >
277 ScvBasis innerNormals;
279 innerNormals[0] = crossProduct<CoordScalar>(scvBasis[1], scvBasis[2]);
280 innerNormals[1] = crossProduct<CoordScalar>(scvBasis[2], scvBasis[0]);
281 innerNormals[2] = crossProduct<CoordScalar>(scvBasis[0], scvBasis[1]);
283 if (!isRightHandSystem(scvBasis))
284 std::for_each(innerNormals.begin(), innerNormals.end(), [] (
auto& n) { n *= -1.0; });
295 template<
class ScvBasis >
299 return abs(tripleProduct<CoordScalar>(scvBasis[0], scvBasis[1], scvBasis[2]));
310 assert(gridView.size(Dune::GeometryTypes::tetrahedron)
311 + gridView.size(Dune::GeometryTypes::pyramid)
312 + gridView.size(Dune::GeometryTypes::prism)
313 + gridView.size(Dune::GeometryTypes::hexahedron) == gridView.size(0));
315 return gridView.size(Dune::GeometryTypes::tetrahedron)
316 * getNumLocalScvfs(Dune::GeometryTypes::tetrahedron)
317 + gridView.size(Dune::GeometryTypes::pyramid)
318 * getNumLocalScvfs(Dune::GeometryTypes::pyramid)
319 + gridView.size(Dune::GeometryTypes::prism)
320 * getNumLocalScvfs(Dune::GeometryTypes::prism)
321 + gridView.size(Dune::GeometryTypes::hexahedron)
322 * getNumLocalScvfs(Dune::GeometryTypes::hexahedron);
329 template<
class ScvBasis >
331 {
return !std::signbit(tripleProduct<CoordScalar>(scvBasis[0], scvBasis[1], scvBasis[2])); }
343 template<
class ElementGeometry,
class ReferenceElement>
345 const ReferenceElement& refElement,
346 unsigned int indexInElement,
350 ScvfPositionsOnIntersection p;
352 DUNE_THROW(Dune::InvalidStateException,
"Mpfa implementation cannot handle faces with more than 4 corners");
359 p[c+1] = eg.global(refElement.position(refElement.subEntity(indexInElement, 1, c, 3), 3));
392 DUNE_THROW(Dune::NotImplemented,
393 "Mpfa scvf corners for dim = 3, dimWorld = 3, corners = " <<
numCorners);
406 unsigned int numIntersectionCorners,
407 unsigned int cornerIdx)
411 switch (numIntersectionCorners)
416 static const std::uint8_t vo = 1;
417 static const std::uint8_t eo = 4;
418 static const std::uint8_t map[3][4] =
420 {0, eo+1, eo+0, vo+0},
421 {0, eo+0, eo+2, vo+1},
422 {0, eo+2, eo+1, vo+2}
425 return ScvfCornerVector( {p[map[cornerIdx][0]],
426 p[map[cornerIdx][1]],
427 p[map[cornerIdx][2]],
428 p[map[cornerIdx][3]]} );
433 static const std::uint8_t vo = 1;
434 static const std::uint8_t eo = 5;
435 static const std::uint8_t map[4][4] =
437 {0, eo+0, eo+2, vo+0},
438 {0, eo+2, eo+1, vo+1},
439 {0, eo+3, eo+0, vo+2},
440 {0, eo+1, eo+3, vo+3}
443 return ScvfCornerVector( {p[map[cornerIdx][0]],
444 p[map[cornerIdx][1]],
445 p[map[cornerIdx][2]],
446 p[map[cornerIdx][3]]} );
449 DUNE_THROW(Dune::NotImplemented,
450 "Mpfa scvf corners for dim = 3, dimWorld = 3, corners = " << numIntersectionCorners);
461 return 0.5*
Dumux::crossProduct(scvfCorners[3]-scvfCorners[0], scvfCorners[2]-scvfCorners[1]).two_norm();
471 if (gt == Dune::GeometryTypes::tetrahedron)
473 else if (gt == Dune::GeometryTypes::pyramid)
475 else if (gt == Dune::GeometryTypes::prism)
477 else if (gt == Dune::GeometryTypes::hexahedron)
480 DUNE_THROW(Dune::NotImplemented,
"Mpfa for 3d geometry type " << gt);
490template<
class Gr
idGeometry>
492 GridGeometry::GridView::dimension,
493 GridGeometry::GridView::dimensionworld>
495 using PrimaryInteractionVolume =
typename GridGeometry::GridIVIndexSets::PrimaryInteractionVolume;
496 using SecondaryInteractionVolume =
typename GridGeometry::GridIVIndexSets::SecondaryInteractionVolume;
498 using VertexMapper =
typename GridGeometry::VertexMapper;
499 using FVElementGeometry =
typename GridGeometry::LocalView;
500 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
501 using ScvfCornerVector =
typename SubControlVolumeFace::Traits::CornerStorage;
503 using GridView =
typename GridGeometry::GridView;
504 static constexpr int dim = GridView::dimension;
506 using Element =
typename GridView::template Codim<0>::Entity;
507 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
508 using CoordScalar =
typename GridView::ctype;
517 template<
class Scalar >
522 return scvfCorners[0];
524 auto ip = scvfCorners.back() - scvfCorners.front();
526 ip += scvfCorners[0];
536 static std::vector<bool>
findGhostVertices(
const GridView& gridView,
const VertexMapper& vertexMapper)
538 std::vector<bool> ghostVertices(gridView.size(dim),
false);
541 if (Dune::MPIHelper::getCommunication().size() == 1)
542 return ghostVertices;
545 if (gridView.ghostSize(0) > 0)
546 DUNE_THROW(Dune::InvalidStateException,
"Mpfa methods in parallel do not work with ghost cells. Use overlap cells instead.");
549 if (gridView.overlapSize(0) == 0)
550 DUNE_THROW(Dune::InvalidStateException,
"Grid no overlapping cells. This is required by mpfa methods in parallel.");
552 for (
const auto& element : elements(gridView))
554 for (
const auto& is : intersections(gridView, element))
556 if (!is.neighbor() && !is.boundary())
558 const auto refElement = referenceElement(element);
560 for (
int isVertex = 0; isVertex < is.geometry().corners(); ++isVertex)
562 const auto vIdxLocal = refElement.subEntity(is.indexInInside(), 1, isVertex, dim);
563 const auto vIdxGlobal = vertexMapper.subIndex(element, vIdxLocal, dim);
564 ghostVertices[vIdxGlobal] =
true;
570 return ghostVertices;
576 {
return !std::is_same<PrimaryInteractionVolume, SecondaryInteractionVolume>::value; }
579 template<
typename V,
typename T>
581 {
return std::find(vector.begin(), vector.end(), value) != vector.end(); }
Helper class to get the required information on an interaction volume.
Definition: helper.hh:494
static constexpr bool considerSecondaryIVs()
Definition: helper.hh:575
static bool vectorContainsValue(const V &vector, const T value)
returns whether or not a value exists in a vector
Definition: helper.hh:580
static std::vector< bool > findGhostVertices(const GridView &gridView, const VertexMapper &vertexMapper)
Returns a vector which maps true to each vertex on processor boundaries and false otherwise.
Definition: helper.hh:536
static GlobalPosition getScvfIntegrationPoint(const ScvfCornerVector &scvfCorners, Scalar q)
Calculates the integration point on an scvf.
Definition: helper.hh:518
static std::size_t getNumLocalScvfs(const Dune::GeometryType >)
Calculates the number of scvfs in a given element geometry type.
Definition: helper.hh:173
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:122
static ScvBasis calculateInnerNormals(const ScvBasis &scvBasis)
Calculates the inner normal vectors to a given scv basis.
Definition: helper.hh:58
static CoordScalar computeScvfArea(const ScvfCornerVector &scvfCorners)
Calculates the area of an scvf.
Definition: helper.hh:166
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:150
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:81
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:92
static bool isRightHandSystem(const ScvBasis &scvBasis)
Checks whether or not a given scv basis forms a right hand system.
Definition: helper.hh:108
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:243
static ScvBasis calculateInnerNormals(const ScvBasis &scvBasis)
Calculates the inner normal vectors to a given scv basis.
Definition: helper.hh:202
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:229
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:296
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:344
static ScvBasis calculateInnerNormals(const ScvBasis &scvBasis)
Calculates the inner normal vectors to a given scv basis.
Definition: helper.hh:275
static bool isRightHandSystem(const ScvBasis &scvBasis)
Checks whether or not a given scv basis forms a right hand system.
Definition: helper.hh:330
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:308
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:405
static std::size_t getNumLocalScvfs(const Dune::GeometryType >)
Calculates the number of scvfs in a given element geometry type.
Definition: helper.hh:469
static CoordScalar computeScvfArea(const ScvfCornerVector &scvfCorners)
Calculates the area of an scvf.
Definition: helper.hh:458
Dimension-specific helper class to get data required for mpfa scheme.
Definition: helper.hh:30
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:671
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.hh:26
Define some often used mathematical functions.
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:220