version 3.8
helper.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_DISCRETIZATION_CC_MPFA_HELPER_HH
13#define DUMUX_DISCRETIZATION_CC_MPFA_HELPER_HH
14
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>
20
21#include <dumux/common/math.hh>
22
23namespace Dumux {
24
29template<class GridGeometry, int dim, int dimWorld>
31
36template<class GridGeometry>
37class MpfaDimensionHelper<GridGeometry, /*dim*/2, /*dimWorld*/2>
38{
39 using GridView = typename GridGeometry::GridView;
40 using CoordScalar = typename GridView::ctype;
41 using Element = typename GridView::template Codim<0>::Entity;
42
43 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
44
45 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
46 using ScvfCornerVector = typename SubControlVolumeFace::Traits::CornerStorage;
47
48 // Container to store the positions of intersections required for scvf
49 // corner computation. In 2d, these are the center plus the two corners
50 using ScvfPositionsOnIntersection = std::array<GlobalPosition, 3>;
51
52public:
57 template< class ScvBasis >
58 static ScvBasis calculateInnerNormals(const ScvBasis& scvBasis)
59 {
60 static const Dune::FieldMatrix<CoordScalar, 2, 2> R = {{0.0, 1.0}, {-1.0, 0.0}};
61
62 ScvBasis innerNormals;
63 R.mv(scvBasis[1], innerNormals[0]);
64 R.mv(scvBasis[0], innerNormals[1]);
65
66 // adjust sign depending on basis being a RHS
67 if (!isRightHandSystem(scvBasis))
68 innerNormals[0] *= -1.0;
69 else
70 innerNormals[1] *= -1.0;
71
72 return innerNormals;
73 }
74
80 template< class ScvBasis >
81 static CoordScalar calculateDetX(const ScvBasis& scvBasis)
82 {
83 using std::abs;
84 return abs(crossProduct<CoordScalar>(scvBasis[0], scvBasis[1]));
85 }
86
92 static std::size_t getGlobalNumScvf(const GridView& gridView)
93 {
94 assert(gridView.size(Dune::GeometryTypes::triangle)
95 + gridView.size(Dune::GeometryTypes::quadrilateral) == gridView.size(0));
96
97 return gridView.size(Dune::GeometryTypes::triangle)
98 * getNumLocalScvfs(Dune::GeometryTypes::triangle)
99 + gridView.size(Dune::GeometryTypes::quadrilateral)
100 * getNumLocalScvfs(Dune::GeometryTypes::quadrilateral);
101 }
102
107 template< class ScvBasis >
108 static bool isRightHandSystem(const ScvBasis& scvBasis)
109 { return !std::signbit(crossProduct<CoordScalar>(scvBasis[0], scvBasis[1])); }
110
121 template<class ElementGeometry, class ReferenceElement>
122 static ScvfPositionsOnIntersection computeScvfCornersOnIntersection(const ElementGeometry& eg,
123 const ReferenceElement& refElement,
124 unsigned int indexInElement,
125 unsigned int numCorners)
126 {
127 ScvfPositionsOnIntersection p;
128
129 // compute facet center and corners
130 p[0] = 0.0;
131 for (unsigned int c = 0; c < numCorners; ++c)
132 {
133 // codim = 1, dim = 2
134 p[c+1] = eg.global(refElement.position(refElement.subEntity(indexInElement, 1, c, 2), 2));
135 p[0] += p[c+1];
136 }
137 p[0] /= numCorners;
138
139 return p;
140 }
141
150 static ScvfCornerVector getScvfCorners(const ScvfPositionsOnIntersection& p,
151 unsigned int numIntersectionCorners,
152 unsigned int cornerIdx)
153 {
154 // make sure the given input is admissible
155 assert(cornerIdx < 2 && "provided index exceeds the number of corners of facets in 2d");
156
157 // create & return the scvf corner vector
158 return cornerIdx == 0 ? ScvfCornerVector({p[0], p[1]})
159 : ScvfCornerVector({p[0], p[2]});
160 }
161
166 static CoordScalar computeScvfArea(const ScvfCornerVector& scvfCorners)
167 { return (scvfCorners[1]-scvfCorners[0]).two_norm(); }
168
173 static std::size_t getNumLocalScvfs(const Dune::GeometryType& gt)
174 {
175 if (gt == Dune::GeometryTypes::triangle)
176 return 6;
177 else if (gt == Dune::GeometryTypes::quadrilateral)
178 return 8;
179 else
180 DUNE_THROW(Dune::NotImplemented, "Mpfa for 2d geometry type " << gt);
181 }
182};
183
189template<class GridGeometry>
190class MpfaDimensionHelper<GridGeometry, /*dim*/2, /*dimWorld*/3>
191: public MpfaDimensionHelper<GridGeometry, 2, 2>
192{
193 using GridView = typename GridGeometry::GridView;
194 using CoordScalar = typename GridView::ctype;
195public:
196
201 template< class ScvBasis >
202 static ScvBasis calculateInnerNormals(const ScvBasis& scvBasis)
203 {
204 // compute vector normal to the basis plane
205 const auto normal = [&scvBasis] ()
206 {
207 auto n = crossProduct<CoordScalar>(scvBasis[0], scvBasis[1]);
208 n /= n.two_norm();
209 return n;
210 } ();
211
212 // compute inner normals using the normal vector
213 ScvBasis innerNormals;
214 innerNormals[0] = crossProduct<CoordScalar>(scvBasis[1], normal);
215 innerNormals[1] = crossProduct<CoordScalar>(normal, scvBasis[0]);
216
217 return innerNormals;
218 }
219
228 template< class ScvBasis >
229 static CoordScalar calculateDetX(const ScvBasis& scvBasis)
230 {
231 using std::abs;
232 return abs(crossProduct<CoordScalar>(scvBasis[0], scvBasis[1]).two_norm());
233 }
234
242 template< class ScvBasis >
243 static constexpr bool isRightHandSystem(const ScvBasis& scvBasis) { return true; }
244};
250template<class GridGeometry>
251class MpfaDimensionHelper<GridGeometry, /*dim*/3, /*dimWorld*/3>
252{
253 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
254 using ScvfCornerVector = typename SubControlVolumeFace::Traits::CornerStorage;
255
256 // Be picky about the dimensions
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;
261
262 // container to store the positions of intersections required for
263 // scvf corner computation. Maximum number of points needed is 9
264 // for the supported geometry types (quadrilateral facet)
265 using ScvfPositionsOnIntersection = std::array<GlobalPosition, 9>;
266
267public:
268
274 template< class ScvBasis >
275 static ScvBasis calculateInnerNormals(const ScvBasis& scvBasis)
276 {
277 ScvBasis innerNormals;
278
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]);
282
283 if (!isRightHandSystem(scvBasis))
284 std::for_each(innerNormals.begin(), innerNormals.end(), [] (auto& n) { n *= -1.0; });
285
286 return innerNormals;
287 }
288
295 template< class ScvBasis >
296 static CoordScalar calculateDetX(const ScvBasis& scvBasis)
297 {
298 using std::abs;
299 return abs(tripleProduct<CoordScalar>(scvBasis[0], scvBasis[1], scvBasis[2]));
300 }
301
308 static std::size_t getGlobalNumScvf(const GridView& gridView)
309 {
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));
314
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);
323 }
324
329 template< class ScvBasis >
330 static bool isRightHandSystem(const ScvBasis& scvBasis)
331 { return !std::signbit(tripleProduct<CoordScalar>(scvBasis[0], scvBasis[1], scvBasis[2])); }
332
343 template<class ElementGeometry, class ReferenceElement>
344 static ScvfPositionsOnIntersection computeScvfCornersOnIntersection(const ElementGeometry& eg,
345 const ReferenceElement& refElement,
346 unsigned int indexInElement,
347 unsigned int numCorners)
348 {
349 // The size of ScvfPositionsOnIntersection doesn't allow for faces with more than four corners!
350 ScvfPositionsOnIntersection p;
351 if (numCorners > 4)
352 DUNE_THROW(Dune::InvalidStateException, "Mpfa implementation cannot handle faces with more than 4 corners");
353
354 // compute facet center and corners
355 p[0] = 0.0;
356 for (unsigned int c = 0; c < numCorners; ++c)
357 {
358 // codim = 1, dim = 3
359 p[c+1] = eg.global(refElement.position(refElement.subEntity(indexInElement, 1, c, 3), 3));
360 p[0] += p[c+1];
361 }
362 p[0] /= numCorners;
363
364 // proceed according to number of corners
365 switch (numCorners)
366 {
367 case 3: // triangle
368 {
369 // add edge midpoints, there are 3 for triangles
370 p[numCorners+1] = p[2] + p[1];
371 p[numCorners+1] /= 2;
372 p[numCorners+2] = p[3] + p[1];
373 p[numCorners+2] /= 2;
374 p[numCorners+3] = p[3] + p[2];
375 p[numCorners+3] /= 2;
376 return p;
377 }
378 case 4: // quadrilateral
379 {
380 // add edge midpoints, there are 4 for quadrilaterals
381 p[numCorners+1] = p[3] + p[1];
382 p[numCorners+1] /= 2;
383 p[numCorners+2] = p[4] + p[2];
384 p[numCorners+2] /= 2;
385 p[numCorners+3] = p[2] + p[1];
386 p[numCorners+3] /= 2;
387 p[numCorners+4] = p[4] + p[3];
388 p[numCorners+4] /= 2;
389 return p;
390 }
391 default:
392 DUNE_THROW(Dune::NotImplemented,
393 "Mpfa scvf corners for dim = 3, dimWorld = 3, corners = " << numCorners);
394 }
395 }
396
405 static ScvfCornerVector getScvfCorners(const ScvfPositionsOnIntersection& p,
406 unsigned int numIntersectionCorners,
407 unsigned int cornerIdx)
408 {
409 // proceed according to number of corners
410 // we assume the ordering according to the above method computeScvfCornersOnIntersection()
411 switch (numIntersectionCorners)
412 {
413 case 3: // triangle
414 {
415 // Only build the maps the first time we encounter a triangle
416 static const std::uint8_t vo = 1; // vertex offset in point vector p
417 static const std::uint8_t eo = 4; // edge offset in point vector p
418 static const std::uint8_t map[3][4] =
419 {
420 {0, eo+1, eo+0, vo+0},
421 {0, eo+0, eo+2, vo+1},
422 {0, eo+2, eo+1, vo+2}
423 };
424
425 return ScvfCornerVector( {p[map[cornerIdx][0]],
426 p[map[cornerIdx][1]],
427 p[map[cornerIdx][2]],
428 p[map[cornerIdx][3]]} );
429 }
430 case 4: // quadrilateral
431 {
432 // Only build the maps the first time we encounter a quadrilateral
433 static const std::uint8_t vo = 1; // vertex offset in point vector p
434 static const std::uint8_t eo = 5; // face offset in point vector p
435 static const std::uint8_t map[4][4] =
436 {
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}
441 };
442
443 return ScvfCornerVector( {p[map[cornerIdx][0]],
444 p[map[cornerIdx][1]],
445 p[map[cornerIdx][2]],
446 p[map[cornerIdx][3]]} );
447 }
448 default:
449 DUNE_THROW(Dune::NotImplemented,
450 "Mpfa scvf corners for dim = 3, dimWorld = 3, corners = " << numIntersectionCorners);
451 }
452 }
453
458 static CoordScalar computeScvfArea(const ScvfCornerVector& scvfCorners)
459 {
460 // after Wolfram alpha quadrilateral area
461 return 0.5*Dumux::crossProduct(scvfCorners[3]-scvfCorners[0], scvfCorners[2]-scvfCorners[1]).two_norm();
462 }
463
469 static std::size_t getNumLocalScvfs(const Dune::GeometryType& gt)
470 {
471 if (gt == Dune::GeometryTypes::tetrahedron)
472 return 12;
473 else if (gt == Dune::GeometryTypes::pyramid)
474 return 16;
475 else if (gt == Dune::GeometryTypes::prism)
476 return 18;
477 else if (gt == Dune::GeometryTypes::hexahedron)
478 return 24;
479 else
480 DUNE_THROW(Dune::NotImplemented, "Mpfa for 3d geometry type " << gt);
481 }
482};
483
490template<class GridGeometry>
491class CCMpfaHelper : public MpfaDimensionHelper<GridGeometry,
492 GridGeometry::GridView::dimension,
493 GridGeometry::GridView::dimensionworld>
494{
495 using PrimaryInteractionVolume = typename GridGeometry::GridIVIndexSets::PrimaryInteractionVolume;
496 using SecondaryInteractionVolume = typename GridGeometry::GridIVIndexSets::SecondaryInteractionVolume;
497
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;
502
503 using GridView = typename GridGeometry::GridView;
504 static constexpr int dim = GridView::dimension;
505
506 using Element = typename GridView::template Codim<0>::Entity;
507 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
508 using CoordScalar = typename GridView::ctype;
509
510public:
517 template< class Scalar >
518 static GlobalPosition getScvfIntegrationPoint(const ScvfCornerVector& scvfCorners, Scalar q)
519 {
520 // ordering -> first corner: facet center, last corner: vertex
521 if (q == 0.0)
522 return scvfCorners[0];
523
524 auto ip = scvfCorners.back() - scvfCorners.front();
525 ip *= q;
526 ip += scvfCorners[0];
527 return ip;
528 }
529
536 static std::vector<bool> findGhostVertices(const GridView& gridView, const VertexMapper& vertexMapper)
537 {
538 std::vector<bool> ghostVertices(gridView.size(dim), false);
539
540 // if not run in parallel, skip the rest
541 if (Dune::MPIHelper::getCommunication().size() == 1)
542 return ghostVertices;
543
544 // mpfa methods cannot yet handle ghost cells
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.");
547
548 // mpfa methods have to have overlapping cells
549 if (gridView.overlapSize(0) == 0)
550 DUNE_THROW(Dune::InvalidStateException, "Grid no overlapping cells. This is required by mpfa methods in parallel.");
551
552 for (const auto& element : elements(gridView))
553 {
554 for (const auto& is : intersections(gridView, element))
555 {
556 if (!is.neighbor() && !is.boundary())
557 {
558 const auto refElement = referenceElement(element);
559
560 for (int isVertex = 0; isVertex < is.geometry().corners(); ++isVertex)
561 {
562 const auto vIdxLocal = refElement.subEntity(is.indexInInside(), 1, isVertex, dim);
563 const auto vIdxGlobal = vertexMapper.subIndex(element, vIdxLocal, dim);
564 ghostVertices[vIdxGlobal] = true;
565 }
566 }
567 }
568 }
569
570 return ghostVertices;
571 }
572
575 static constexpr bool considerSecondaryIVs()
576 { return !std::is_same<PrimaryInteractionVolume, SecondaryInteractionVolume>::value; }
577
579 template<typename V, typename T>
580 static bool vectorContainsValue(const V& vector, const T value)
581 { return std::find(vector.begin(), vector.end(), value) != vector.end(); }
582};
583
584} // end namespace Dumux
585
586#endif
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 &gt)
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 &gt)
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:642
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
Definition: adapt.hh:17