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>
41template<
class Gr
idGeometry,
int dim,
int dimWorld>
48template<
class Gr
idGeometry>
51 using GridView =
typename GridGeometry::GridView;
52 using CoordScalar =
typename GridView::ctype;
53 using Element =
typename GridView::template Codim<0>::Entity;
55 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
57 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
58 using ScvfCornerVector =
typename SubControlVolumeFace::Traits::CornerStorage;
62 using ScvfPositionsOnIntersection = std::array<GlobalPosition, 3>;
69 template<
class ScvBasis >
72 static const Dune::FieldMatrix<CoordScalar, 2, 2> R = {{0.0, 1.0}, {-1.0, 0.0}};
74 ScvBasis innerNormals;
75 R.mv(scvBasis[1], innerNormals[0]);
76 R.mv(scvBasis[0], innerNormals[1]);
79 if (!isRightHandSystem(scvBasis))
80 innerNormals[0] *= -1.0;
82 innerNormals[1] *= -1.0;
92 template<
class ScvBasis >
96 return abs(crossProduct<CoordScalar>(scvBasis[0], scvBasis[1]));
106 assert(gridView.size(Dune::GeometryTypes::triangle)
107 + gridView.size(Dune::GeometryTypes::quadrilateral) == gridView.size(0));
109 return gridView.size(Dune::GeometryTypes::triangle)
110 * getNumLocalScvfs(Dune::GeometryTypes::triangle)
111 + gridView.size(Dune::GeometryTypes::quadrilateral)
112 * getNumLocalScvfs(Dune::GeometryTypes::quadrilateral);
119 template<
class ScvBasis >
121 {
return !std::signbit(crossProduct<CoordScalar>(scvBasis[0], scvBasis[1])); }
133 template<
class ElementGeometry,
class ReferenceElement>
135 const ReferenceElement& refElement,
136 unsigned int indexInElement,
137 unsigned int numCorners)
139 ScvfPositionsOnIntersection p;
143 for (
unsigned int c = 0; c < numCorners; ++c)
146 p[c+1] = eg.global(refElement.position(refElement.subEntity(indexInElement, 1, c, 2), 2));
163 unsigned int numIntersectionCorners,
164 unsigned int cornerIdx)
167 assert(cornerIdx < 2 &&
"provided index exceeds the number of corners of facets in 2d");
170 return cornerIdx == 0 ? ScvfCornerVector({p[0], p[1]})
171 : ScvfCornerVector({p[0], p[2]});
179 {
return (scvfCorners[1]-scvfCorners[0]).two_norm(); }
187 if (gt == Dune::GeometryTypes::triangle)
189 else if (gt == Dune::GeometryTypes::quadrilateral)
192 DUNE_THROW(Dune::NotImplemented,
"Mpfa for 2d geometry type " << gt);
201template<
class Gr
idGeometry>
205 using GridView =
typename GridGeometry::GridView;
206 using CoordScalar =
typename GridView::ctype;
213 template<
class ScvBasis >
217 const auto normal = [&scvBasis] ()
219 auto n = crossProduct<CoordScalar>(scvBasis[0], scvBasis[1]);
225 ScvBasis innerNormals;
226 innerNormals[0] = crossProduct<CoordScalar>(scvBasis[1],
normal);
227 innerNormals[1] = crossProduct<CoordScalar>(
normal, scvBasis[0]);
240 template<
class ScvBasis >
244 return abs(crossProduct<CoordScalar>(scvBasis[0], scvBasis[1]).two_norm());
254 template<
class ScvBasis >
262template<
class Gr
idGeometry>
265 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
266 using ScvfCornerVector =
typename SubControlVolumeFace::Traits::CornerStorage;
269 using GridView =
typename GridGeometry::GridView;
270 using Element =
typename GridView::template Codim<0>::Entity;
271 using CoordScalar =
typename GridView::ctype;
272 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
277 using ScvfPositionsOnIntersection = std::array<GlobalPosition, 9>;
286 template<
class ScvBasis >
289 ScvBasis innerNormals;
291 innerNormals[0] = crossProduct<CoordScalar>(scvBasis[1], scvBasis[2]);
292 innerNormals[1] = crossProduct<CoordScalar>(scvBasis[2], scvBasis[0]);
293 innerNormals[2] = crossProduct<CoordScalar>(scvBasis[0], scvBasis[1]);
295 if (!isRightHandSystem(scvBasis))
296 std::for_each(innerNormals.begin(), innerNormals.end(), [] (
auto& n) { n *= -1.0; });
307 template<
class ScvBasis >
311 return abs(tripleProduct<CoordScalar>(scvBasis[0], scvBasis[1], scvBasis[2]));
322 assert(gridView.size(Dune::GeometryTypes::tetrahedron)
323 + gridView.size(Dune::GeometryTypes::pyramid)
324 + gridView.size(Dune::GeometryTypes::prism)
325 + gridView.size(Dune::GeometryTypes::hexahedron) == gridView.size(0));
327 return gridView.size(Dune::GeometryTypes::tetrahedron)
328 * getNumLocalScvfs(Dune::GeometryTypes::tetrahedron)
329 + gridView.size(Dune::GeometryTypes::pyramid)
330 * getNumLocalScvfs(Dune::GeometryTypes::pyramid)
331 + gridView.size(Dune::GeometryTypes::prism)
332 * getNumLocalScvfs(Dune::GeometryTypes::prism)
333 + gridView.size(Dune::GeometryTypes::hexahedron)
334 * getNumLocalScvfs(Dune::GeometryTypes::hexahedron);
341 template<
class ScvBasis >
343 {
return !std::signbit(tripleProduct<CoordScalar>(scvBasis[0], scvBasis[1], scvBasis[2])); }
355 template<
class ElementGeometry,
class ReferenceElement>
357 const ReferenceElement& refElement,
358 unsigned int indexInElement,
359 unsigned int numCorners)
362 ScvfPositionsOnIntersection p;
364 DUNE_THROW(Dune::InvalidStateException,
"Mpfa implementation cannot handle faces with more than 4 corners");
368 for (
unsigned int c = 0; c < numCorners; ++c)
371 p[c+1] = eg.global(refElement.position(refElement.subEntity(indexInElement, 1, c, 3), 3));
382 p[numCorners+1] = p[2] + p[1];
383 p[numCorners+1] /= 2;
384 p[numCorners+2] = p[3] + p[1];
385 p[numCorners+2] /= 2;
386 p[numCorners+3] = p[3] + p[2];
387 p[numCorners+3] /= 2;
393 p[numCorners+1] = p[3] + p[1];
394 p[numCorners+1] /= 2;
395 p[numCorners+2] = p[4] + p[2];
396 p[numCorners+2] /= 2;
397 p[numCorners+3] = p[2] + p[1];
398 p[numCorners+3] /= 2;
399 p[numCorners+4] = p[4] + p[3];
400 p[numCorners+4] /= 2;
404 DUNE_THROW(Dune::NotImplemented,
405 "Mpfa scvf corners for dim = 3, dimWorld = 3, corners = " << numCorners);
418 unsigned int numIntersectionCorners,
419 unsigned int cornerIdx)
423 switch (numIntersectionCorners)
428 static const std::uint8_t vo = 1;
429 static const std::uint8_t eo = 4;
430 static const std::uint8_t map[3][4] =
432 {0, eo+1, eo+0, vo+0},
433 {0, eo+0, eo+2, vo+1},
434 {0, eo+2, eo+1, vo+2}
437 return ScvfCornerVector( {p[map[cornerIdx][0]],
438 p[map[cornerIdx][1]],
439 p[map[cornerIdx][2]],
440 p[map[cornerIdx][3]]} );
445 static const std::uint8_t vo = 1;
446 static const std::uint8_t eo = 5;
447 static const std::uint8_t map[4][4] =
449 {0, eo+0, eo+2, vo+0},
450 {0, eo+2, eo+1, vo+1},
451 {0, eo+3, eo+0, vo+2},
452 {0, eo+1, eo+3, vo+3}
455 return ScvfCornerVector( {p[map[cornerIdx][0]],
456 p[map[cornerIdx][1]],
457 p[map[cornerIdx][2]],
458 p[map[cornerIdx][3]]} );
461 DUNE_THROW(Dune::NotImplemented,
462 "Mpfa scvf corners for dim = 3, dimWorld = 3, corners = " << numIntersectionCorners);
473 return 0.5*
Dumux::crossProduct(scvfCorners[3]-scvfCorners[0], scvfCorners[2]-scvfCorners[1]).two_norm();
483 if (gt == Dune::GeometryTypes::tetrahedron)
485 else if (gt == Dune::GeometryTypes::pyramid)
487 else if (gt == Dune::GeometryTypes::prism)
489 else if (gt == Dune::GeometryTypes::hexahedron)
492 DUNE_THROW(Dune::NotImplemented,
"Mpfa for 3d geometry type " << gt);
502template<
class Gr
idGeometry>
504 GridGeometry::GridView::dimension,
505 GridGeometry::GridView::dimensionworld>
507 using PrimaryInteractionVolume =
typename GridGeometry::GridIVIndexSets::PrimaryInteractionVolume;
508 using SecondaryInteractionVolume =
typename GridGeometry::GridIVIndexSets::SecondaryInteractionVolume;
510 using VertexMapper =
typename GridGeometry::VertexMapper;
511 using FVElementGeometry =
typename GridGeometry::LocalView;
512 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
513 using ScvfCornerVector =
typename SubControlVolumeFace::Traits::CornerStorage;
515 using GridView =
typename GridGeometry::GridView;
516 static constexpr int dim = GridView::dimension;
518 using Element =
typename GridView::template Codim<0>::Entity;
519 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
520 using CoordScalar =
typename GridView::ctype;
529 template<
class Scalar >
534 return scvfCorners[0];
536 auto ip = scvfCorners.back() - scvfCorners.front();
538 ip += scvfCorners[0];
548 static std::vector<bool>
findGhostVertices(
const GridView& gridView,
const VertexMapper& vertexMapper)
550 std::vector<bool> ghostVertices(gridView.size(dim),
false);
553 if (Dune::MPIHelper::getCollectiveCommunication().size() == 1)
554 return ghostVertices;
557 if (gridView.ghostSize(0) > 0)
558 DUNE_THROW(Dune::InvalidStateException,
"Mpfa methods in parallel do not work with ghost cells. Use overlap cells instead.");
561 if (gridView.overlapSize(0) == 0)
562 DUNE_THROW(Dune::InvalidStateException,
"Grid no overlaping cells. This is required by mpfa methods in parallel.");
564 for (
const auto& element : elements(gridView))
566 for (
const auto& is : intersections(gridView, element))
568 if (!is.neighbor() && !is.boundary())
570 const auto refElement = referenceElement(element);
572 for (
int isVertex = 0; isVertex < is.geometry().corners(); ++isVertex)
574 const auto vIdxLocal = refElement.subEntity(is.indexInInside(), 1, isVertex, dim);
575 const auto vIdxGlobal = vertexMapper.subIndex(element, vIdxLocal, dim);
576 ghostVertices[vIdxGlobal] =
true;
582 return ghostVertices;
588 {
return !std::is_same<PrimaryInteractionVolume, SecondaryInteractionVolume>::value; }
591 template<
typename V,
typename T>
593 {
return std::find(vector.begin(), vector.end(), value) != vector.end(); }
Define some often used mathematical functions.
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: geometry/normal.hh:36
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:640
Dimension-specific helper class to get data required for mpfa scheme.
Definition: helper.hh:42
static std::size_t getNumLocalScvfs(const Dune::GeometryType >)
Calculates the number of scvfs in a given element geometry type.
Definition: helper.hh:185
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:134
static ScvBasis calculateInnerNormals(const ScvBasis &scvBasis)
Calculates the inner normal vectors to a given scv basis.
Definition: helper.hh:70
static CoordScalar computeScvfArea(const ScvfCornerVector &scvfCorners)
Calculates the area of an scvf.
Definition: helper.hh:178
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:162
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:93
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:104
static bool isRightHandSystem(const ScvBasis &scvBasis)
Checks whether or not a given scv basis forms a right hand system.
Definition: helper.hh:120
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:255
static ScvBasis calculateInnerNormals(const ScvBasis &scvBasis)
Calculates the inner normal vectors to a given scv basis.
Definition: helper.hh:214
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:241
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:308
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:356
static ScvBasis calculateInnerNormals(const ScvBasis &scvBasis)
Calculates the inner normal vectors to a given scv basis.
Definition: helper.hh:287
static bool isRightHandSystem(const ScvBasis &scvBasis)
Checks whether or not a given scv basis forms a right hand system.
Definition: helper.hh:342
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:320
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:417
static std::size_t getNumLocalScvfs(const Dune::GeometryType >)
Calculates the number of scvfs in a given element geometry type.
Definition: helper.hh:481
static CoordScalar computeScvfArea(const ScvfCornerVector &scvfCorners)
Calculates the area of an scvf.
Definition: helper.hh:470
Helper class to get the required information on an interaction volume.
Definition: helper.hh:506
static constexpr bool considerSecondaryIVs()
Definition: helper.hh:587
static bool vectorContainsValue(const V &vector, const T value)
returns whether or not a value exists in a vector
Definition: helper.hh:592
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:548
static GlobalPosition getScvfIntegrationPoint(const ScvfCornerVector &scvfCorners, Scalar q)
Calculates the integration point on an scvf.
Definition: helper.hh:530