vcfvstencil.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 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef EWOMS_VCFV_STENCIL_HH
29#define EWOMS_VCFV_STENCIL_HH
30
31#include <dune/common/version.hh>
32#include <dune/geometry/referenceelements.hh>
33
34#include <dune/grid/common/intersectioniterator.hh>
35#include <dune/grid/common/mcmgmapper.hh>
36
38
39#if HAVE_DUNE_LOCALFUNCTIONS
40#include <dune/localfunctions/lagrange/pqkfactory.hh>
41#include <vector>
42#endif // HAVE_DUNE_LOCALFUNCTIONS
43
44#include <array>
45#include <algorithm>
46#include <cassert>
47#include <cmath>
48#include <stdexcept>
49#include <string>
50#include <utility>
51
52namespace Opm {
53
57enum class ElementType {
58 none,
59 simplex,
60 cube,
61};
62
66template <class Scalar, unsigned dim, ElementType basicGeomType>
67class VcfvScvGeometries;
68
70// local geometries for 1D elements
72template <class Scalar>
73class VcfvScvGeometries<Scalar, /*dim=*/1, ElementType::cube>
74{
75 enum { dim = 1 };
76 enum { numScv = 2 };
77
78public:
79 using ScvLocalGeometry = QuadrialteralQuadratureGeometry<Scalar, dim>;
80
81 static void init()
82 {
83 // 1D LINE SEGMENTS
84 const Scalar scvCorners[numScv][ScvLocalGeometry::numCorners][dim] = {
85 { // corners of the first sub control volume
86 { 0.0 },
87 { 0.5 }
88 },
89
90 { // corners of the second sub control volume
91 { 0.5 },
92 { 1.0 }
93 }
94 };
95 for (unsigned scvIdx = 0; scvIdx < numScv; ++scvIdx) {
96 scvGeoms_[scvIdx].setCorners(scvCorners[scvIdx], ScvLocalGeometry::numCorners);
97 }
98 }
99
100 static const ScvLocalGeometry& get(unsigned scvIdx)
101 { return scvGeoms_[scvIdx]; }
102
103private:
104 static std::array<ScvLocalGeometry, numScv> scvGeoms_;
105};
106
107template <class Scalar>
108class VcfvScvGeometries<Scalar, /*dim=*/1, ElementType::simplex>
109{
110 enum { dim = 1 };
111 enum { numScv = 2 };
112
113public:
114 using ScvLocalGeometry = QuadrialteralQuadratureGeometry<Scalar, dim>;
115
116 static const ScvLocalGeometry& get(unsigned)
117 {
118 throw std::logic_error("Not implemented: VcfvScvGeometries<Scalar, 1, ElementType::simplex>");
119 }
120};
121
123// local geometries for 2D elements
125template <class Scalar>
126class VcfvScvGeometries<Scalar, /*dim=*/2, ElementType::simplex>
127{
128 enum { dim = 2 };
129 enum { numScv = 3 };
130
131public:
132 using ScvLocalGeometry = QuadrialteralQuadratureGeometry<Scalar, dim>;
133
134 static const ScvLocalGeometry& get(unsigned scvIdx)
135 { return scvGeoms_[scvIdx]; }
136
137 static void init()
138 {
139 // 2D SIMPLEX
140 const Scalar scvCorners[numScv][ScvLocalGeometry::numCorners][dim] =
141 {
142 { // SCV 0 corners
143 { 0.0, 0.0 },
144 { 1.0/2, 0.0 },
145 { 1.0/3, 1.0/3 },
146 { 0.0, 1.0/2 },
147 },
148
149 { // SCV 1 corners
150 { 1.0/2, 0.0 },
151 { 1.0, 0.0 },
152 { 1.0/3, 1.0/3 },
153 { 1.0/2, 1.0/2 },
154 },
155
156 { // SCV 2 corners
157 { 0.0, 1.0/2 },
158 { 1.0/3, 1.0/3 },
159 { 0.0, 1.0 },
160 { 1.0/2, 1.0/2 },
161 }
162 };
163
164 for (unsigned scvIdx = 0; scvIdx < numScv; ++scvIdx) {
165 scvGeoms_[scvIdx].setCorners(scvCorners[scvIdx], ScvLocalGeometry::numCorners);
166 }
167 }
168
169private:
170 static std::array<ScvLocalGeometry, numScv> scvGeoms_;
171};
172
173template <class Scalar>
174class VcfvScvGeometries<Scalar, /*dim=*/2, ElementType::cube>
175{
176 enum { dim = 2 };
177 enum { numScv = 4 };
178
179public:
180 using ScvLocalGeometry = QuadrialteralQuadratureGeometry<Scalar, dim>;
181
182 static const ScvLocalGeometry& get(unsigned scvIdx)
183 { return scvGeoms_[scvIdx]; }
184
185 static void init()
186 {
187 // 2D CUBE
188 const Scalar scvCorners[numScv][ScvLocalGeometry::numCorners][dim] =
189 {
190 { // SCV 0 corners
191 { 0.0, 0.0 },
192 { 0.5, 0.0 },
193 { 0.0, 0.5 },
194 { 0.5, 0.5 }
195 },
196
197 { // SCV 1 corners
198 { 0.5, 0.0 },
199 { 1.0, 0.0 },
200 { 0.5, 0.5 },
201 { 1.0, 0.5 }
202 },
203
204 { // SCV 2 corners
205 { 0.0, 0.5 },
206 { 0.5, 0.5 },
207 { 0.0, 1.0 },
208 { 0.5, 1.0 }
209 },
210
211 { // SCV 3 corners
212 { 0.5, 0.5 },
213 { 1.0, 0.5 },
214 { 0.5, 1.0 },
215 { 1.0, 1.0 }
216 }
217 };
218
219 for (unsigned scvIdx = 0; scvIdx < numScv; ++scvIdx) {
220 scvGeoms_[scvIdx].setCorners(scvCorners[scvIdx], ScvLocalGeometry::numCorners);
221 }
222 }
223
224 static std::array<ScvLocalGeometry, numScv> scvGeoms_;
225};
226
228// local geometries for 3D elements
230template <class Scalar>
231class VcfvScvGeometries<Scalar, /*dim=*/3, ElementType::simplex>
232{
233 enum { dim = 3 };
234 enum { numScv = 4 };
235
236public:
237 using ScvLocalGeometry = QuadrialteralQuadratureGeometry<Scalar, dim>;
238
239 static const ScvLocalGeometry& get(unsigned scvIdx)
240 { return scvGeoms_[scvIdx]; }
241
242 static void init()
243 {
244 // 3D SIMPLEX
245 const Scalar scvCorners[numScv][ScvLocalGeometry::numCorners][dim] =
246 {
247 { // SCV 0 corners
248 { 0.0, 0.0, 0.0 },
249 { 1.0/2, 0.0, 0.0 },
250 { 0.0, 1.0/2, 0.0 },
251 { 1.0/3, 1.0/3, 0.0 },
252
253 { 0.0, 0.0, 0.5 },
254 { 1.0/3, 0.0, 1.0/3 },
255 { 0.0, 1.0/3, 1.0/3 },
256 { 1.0/4, 1.0/4, 1.0/4 },
257 },
258
259 { // SCV 1 corners
260 { 1.0/2, 0.0, 0.0 },
261 { 1.0, 0.0, 0.0 },
262 { 1.0/3, 1.0/3, 0.0 },
263 { 1.0/2, 1.0/2, 0.0 },
264
265 { 1.0/3, 0.0, 1.0/3 },
266 { 1.0/2, 0.0, 1.0/2 },
267 { 1.0/4, 1.0/4, 1.0/4 },
268 { 1.0/3, 1.0/3, 1.0/3 },
269 },
270
271 { // SCV 2 corners
272 { 0.0, 1.0/2, 0.0 },
273 { 1.0/3, 1.0/3, 0.0 },
274 { 0.0, 1.0, 0.0 },
275 { 1.0/2, 1.0/2, 0.0 },
276
277 { 0.0, 1.0/3, 1.0/3 },
278 { 1.0/4, 1.0/4, 1.0/4 },
279 { 0.0, 1.0/2, 1.0/2 },
280 { 1.0/3, 1.0/3, 1.0/3 },
281 },
282
283 { // SCV 3 corners
284 { 0.0, 0.0, 1.0/2 },
285 { 1.0/3, 0.0, 1.0/3 },
286 { 0.0, 1.0/3, 1.0/3 },
287 { 1.0/4, 1.0/4, 1.0/4 },
288
289 { 0.0, 0.0, 1.0 },
290 { 1.0/2, 0.0, 1.0/2 },
291 { 0.0, 1.0/2, 1.0/2 },
292 { 1.0/3, 1.0/3, 1.0/3 },
293 }
294 };
295
296 for (unsigned scvIdx = 0; scvIdx < numScv; ++scvIdx) {
297 scvGeoms_[scvIdx].setCorners(scvCorners[scvIdx], ScvLocalGeometry::numCorners);
298 }
299 }
300
301private:
302 static std::array<ScvLocalGeometry, numScv> scvGeoms_;
303};
304
305template <class Scalar>
306class VcfvScvGeometries<Scalar, /*dim=*/3, ElementType::cube>
307{
308 enum { dim = 3 };
309 enum { numScv = 8 };
310
311public:
312 using ScvLocalGeometry = QuadrialteralQuadratureGeometry<Scalar, dim>;
313
314 static const ScvLocalGeometry& get(unsigned scvIdx)
315 { return scvGeoms_[scvIdx]; }
316
317 static void init()
318 {
319 // 3D CUBE
320 const Scalar scvCorners[numScv][ScvLocalGeometry::numCorners][dim] =
321 {
322 { // SCV 0 corners
323 { 0.0, 0.0, 0.0 },
324 { 1.0/2, 0.0, 0.0 },
325 { 0.0, 1.0/2, 0.0 },
326 { 1.0/2, 1.0/2, 0.0 },
327
328 { 0.0, 0.0, 1.0/2 },
329 { 1.0/2, 0.0, 1.0/2 },
330 { 0.0, 1.0/2, 1.0/2 },
331 { 1.0/2, 1.0/2, 1.0/2 },
332 },
333
334 { // SCV 1 corners
335 { 1.0/2, 0.0, 0.0 },
336 { 1.0, 0.0, 0.0 },
337 { 1.0/2, 1.0/2, 0.0 },
338 { 1.0, 1.0/2, 0.0 },
339
340 { 1.0/2, 0.0, 1.0/2 },
341 { 1.0, 0.0, 1.0/2 },
342 { 1.0/2, 1.0/2, 1.0/2 },
343 { 1.0, 1.0/2, 1.0/2 },
344 },
345
346 { // SCV 2 corners
347 { 0.0, 1.0/2, 0.0 },
348 { 1.0/2, 1.0/2, 0.0 },
349 { 0.0, 1.0, 0.0 },
350 { 1.0/2, 1.0, 0.0 },
351
352 { 0.0, 1.0/2, 1.0/2 },
353 { 1.0/2, 1.0/2, 1.0/2 },
354 { 0.0, 1.0, 1.0/2 },
355 { 1.0/2, 1.0, 1.0/2 },
356 },
357
358 { // SCV 3 corners
359 { 1.0/2, 1.0/2, 0.0 },
360 { 1.0, 1.0/2, 0.0 },
361 { 1.0/2, 1.0, 0.0 },
362 { 1.0, 1.0, 0.0 },
363
364 { 1.0/2, 1.0/2, 1.0/2 },
365 { 1.0, 1.0/2, 1.0/2 },
366 { 1.0/2, 1.0, 1.0/2 },
367 { 1.0, 1.0, 1.0/2 },
368 },
369
370 { // SCV 4 corners
371 { 0.0, 0.0, 1.0/2 },
372 { 1.0/2, 0.0, 1.0/2 },
373 { 0.0, 1.0/2, 1.0/2 },
374 { 1.0/2, 1.0/2, 1.0/2 },
375
376 { 0.0, 0.0, 1.0 },
377 { 1.0/2, 0.0, 1.0 },
378 { 0.0, 1.0/2, 1.0 },
379 { 1.0/2, 1.0/2, 1.0 },
380 },
381
382 { // SCV 5 corners
383 { 1.0/2, 0.0, 1.0/2 },
384 { 1.0, 0.0, 1.0/2 },
385 { 1.0/2, 1.0/2, 1.0/2 },
386 { 1.0, 1.0/2, 1.0/2 },
387
388 { 1.0/2, 0.0, 1.0 },
389 { 1.0, 0.0, 1.0 },
390 { 1.0/2, 1.0/2, 1.0 },
391 { 1.0, 1.0/2, 1.0 },
392 },
393
394 { // SCV 6 corners
395 { 0.0, 1.0/2, 1.0/2 },
396 { 1.0/2, 1.0/2, 1.0/2 },
397 { 0.0, 1.0, 1.0/2 },
398 { 1.0/2, 1.0, 1.0/2 },
399
400 { 0.0, 1.0/2, 1.0 },
401 { 1.0/2, 1.0/2, 1.0 },
402 { 0.0, 1.0, 1.0 },
403 { 1.0/2, 1.0, 1.0 },
404 },
405
406 { // SCV 7 corners
407 { 1.0/2, 1.0/2, 1.0/2 },
408 { 1.0, 1.0/2, 1.0/2 },
409 { 1.0/2, 1.0, 1.0/2 },
410 { 1.0, 1.0, 1.0/2 },
411
412 { 1.0/2, 1.0/2, 1.0 },
413 { 1.0, 1.0/2, 1.0 },
414 { 1.0/2, 1.0, 1.0 },
415 { 1.0, 1.0, 1.0 },
416 },
417 };
418
419 for (unsigned scvIdx = 0; scvIdx < numScv; ++scvIdx) {
420 scvGeoms_[scvIdx].setCorners(scvCorners[scvIdx], ScvLocalGeometry::numCorners);
421 }
422 }
423private:
424 static std::array<ScvLocalGeometry, numScv> scvGeoms_;
425};
426
450template <class Scalar, class GridView>
452{
453 enum { dim = GridView::dimension };
454 enum { dimWorld = GridView::dimensionworld };
455 enum { maxNC = dim < 3 ? 4 : 8 };
456 enum { maxNE = dim < 3 ? 4 : 12 };
457 enum { maxNF = dim < 3 ? 1 : 6 };
458 enum { maxBF = dim < 3 ? 8 : 24 };
459 using CoordScalar = typename GridView::ctype;
460 using Element = typename GridView::Traits::template Codim<0>::Entity;
461
462public:
463 using Entity = typename GridView::Traits::template Codim<dim>::Entity;
464
465private:
466 using Geometry = typename Element::Geometry;
467 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
468 using GlobalPosition = Dune::FieldVector<CoordScalar, dimWorld>;
469 using LocalPosition = Dune::FieldVector<CoordScalar, dim>;
470 using IntersectionIterator = typename GridView::IntersectionIterator;
471
473
474#if HAVE_DUNE_LOCALFUNCTIONS
475 using LocalFiniteElementCache = Dune::PQkLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
476 using LocalFiniteElement = typename LocalFiniteElementCache::FiniteElementType;
477 using LocalBasisTraits = typename LocalFiniteElement::Traits::LocalBasisType::Traits;
478 using ShapeJacobian = typename LocalBasisTraits::JacobianType;
479#endif // HAVE_DUNE_LOCALFUNCTIONS
480
481 Scalar quadrilateralArea(const GlobalPosition& p0,
482 const GlobalPosition& p1,
483 const GlobalPosition& p2,
484 const GlobalPosition& p3)
485 {
486 return 0.5 * std::abs((p3[0] - p1[0]) * (p2[1] - p0[1]) -
487 (p3[1] - p1[1]) * (p2[0] - p0[0]));
488 }
489
490 void crossProduct(DimVector& c, const DimVector& a, const DimVector& b)
491 {
492 c[0] = a[1] * b[2] - a[2] * b[1];
493 c[1] = a[2] * b[0] - a[0] * b[2];
494 c[2] = a[0] * b[1] - a[1] * b[0];
495 }
496
497 Scalar pyramidVolume(const GlobalPosition& p0,
498 const GlobalPosition& p1,
499 const GlobalPosition& p2,
500 const GlobalPosition& p3,
501 const GlobalPosition& p4)
502 {
503 DimVector a(p2); a -= p0;
504 DimVector b(p3); b -= p1;
505
506 DimVector n;
507 crossProduct(n, a, b);
508
509 a = p4; a -= p0;
510
511 return 1.0 / 6.0 * (n * a);
512 }
513
514 Scalar prismVolume(const GlobalPosition& p0,
515 const GlobalPosition& p1,
516 const GlobalPosition& p2,
517 const GlobalPosition& p3,
518 const GlobalPosition& p4,
519 const GlobalPosition& p5)
520 {
521 DimVector a(p4);
522 for (unsigned k = 0; k < dimWorld; ++k) {
523 a[k] -= p0[k];
524 }
525 DimVector b(p1);
526 for (unsigned k = 0; k < dimWorld; ++k) {
527 b[k] -= p3[k];
528 }
529 DimVector m;
530 crossProduct(m, a, b);
531
532 for (unsigned k = 0; k < dimWorld; ++k) {
533 a[k] = p1[k] - p0[k];
534 }
535 for (unsigned k = 0; k < dimWorld; ++k) {
536 b[k] = p2[k] - p0[k];
537 }
538 DimVector n;
539 crossProduct(n, a, b);
540 n += m;
541
542 for (unsigned k = 0; k < dimWorld; ++k) {
543 a[k] = p5[k] - p0[k];
544 }
545
546 return std::abs(1.0 / 6.0 * (n * a));
547 }
548
549 Scalar hexahedronVolume(const GlobalPosition& p0,
550 const GlobalPosition& p1,
551 const GlobalPosition& p2,
552 const GlobalPosition& p3,
553 const GlobalPosition& p4,
554 const GlobalPosition& p5,
555 const GlobalPosition& p6,
556 const GlobalPosition& p7)
557 {
558 return prismVolume(p0, p1, p2, p4, p5, p6) +
559 prismVolume(p0, p2, p3, p4, p6, p7);
560 }
561
562 void normalOfQuadrilateral3D(DimVector& normal,
563 const GlobalPosition& p0,
564 const GlobalPosition& p1,
565 const GlobalPosition& p2,
566 const GlobalPosition& p3)
567 {
568 DimVector a(p2);
569 for (unsigned k = 0; k < dimWorld; ++k) {
570 a[k] -= p0[k];
571 }
572 DimVector b(p3);
573 for (unsigned k = 0; k < dimWorld; ++k) {
574 b[k] -= p1[k];
575 }
576
577 crossProduct(normal, a, b);
578 normal *= 0.5;
579 }
580
581 Scalar quadrilateralArea3D(const GlobalPosition& p0,
582 const GlobalPosition& p1,
583 const GlobalPosition& p2,
584 const GlobalPosition& p3)
585 {
586 DimVector normal;
587 normalOfQuadrilateral3D(normal, p0, p1, p2, p3);
588 return normal.two_norm();
589 }
590
591 void getFaceIndices(unsigned numElemVertices, unsigned k, unsigned& leftFace, unsigned& rightFace)
592 {
593 static const unsigned edgeToFaceTet[2][6] = {
594 {1, 0, 3, 2, 1, 3},
595 {0, 2, 0, 1, 3, 2}
596 };
597 static const unsigned edgeToFacePyramid[2][8] = {
598 {1, 2, 3, 4, 1, 3, 4, 2},
599 {0, 0, 0, 0, 3, 2, 1, 4}
600 };
601 static const unsigned edgeToFacePrism[2][9] = {
602 {1, 0, 2, 0, 1, 2, 4, 4, 4},
603 {0, 2, 1, 3, 3, 3, 0, 1, 2}
604 };
605 static const unsigned edgeToFaceHex[2][12] = {
606 {0, 2, 3, 1, 4, 1, 2, 4, 0, 5, 5, 3},
607 {2, 1, 0, 3, 0, 4, 4, 3, 5, 1, 2, 5}
608 };
609
610 switch (numElemVertices) {
611 case 4:
612 leftFace = edgeToFaceTet[0][k];
613 rightFace = edgeToFaceTet[1][k];
614 break;
615 case 5:
616 leftFace = edgeToFacePyramid[0][k];
617 rightFace = edgeToFacePyramid[1][k];
618 break;
619 case 6:
620 leftFace = edgeToFacePrism[0][k];
621 rightFace = edgeToFacePrism[1][k];
622 break;
623 case 8:
624 leftFace = edgeToFaceHex[0][k];
625 rightFace = edgeToFaceHex[1][k];
626 break;
627 default:
628 throw std::logic_error("Not implemented: VcfvStencil::getFaceIndices for "
629 + std::to_string(numElemVertices) + " vertices");
630 break;
631 }
632 }
633
634 void getEdgeIndices(unsigned numElemVertices, unsigned face, unsigned vert,
635 unsigned& leftEdge, unsigned& rightEdge)
636 {
637 static const int faceAndVertexToLeftEdgeTet[4][4] = {
638 { 0, 0, 2, -1},
639 { 0, 0, -1, 3},
640 { 1, -1, 1, 3},
641 {-1, 2, 2, 4}
642 };
643 static const int faceAndVertexToRightEdgeTet[4][4] = {
644 { 1, 2, 1, -1},
645 { 3, 4, -1, 4},
646 { 3, -1, 5, 5},
647 {-1, 4, 5, 5}
648 };
649 static const int faceAndVertexToLeftEdgePyramid[5][5] = {
650 { 0, 2, 3, 1, -1},
651 { 0, -1, 0, -1, 4},
652 {-1, 1, -1, 1, 5},
653 { 2, 2, -1, -1, 4},
654 {-1, -1, 3, 3, 7}
655 };
656 static const int faceAndVertexToRightEdgePyramid[5][5] = {
657 { 2, 1, 0, 3, -1},
658 { 4, -1, 6, -1, 6},
659 {-1, 5, -1, 7, 7},
660 { 4, 5, -1, -1, 5},
661 {-1, -1, 6, 7, 6}
662 };
663 static const int faceAndVertexToLeftEdgePrism[5][6] = {
664 { 3, 3, -1, 0, 1, -1},
665 { 4, -1, 4, 0, -1, 2},
666 {-1, 5, 5, -1, 1, 2},
667 { 3, 3, 5, -1, -1, -1},
668 {-1, -1, -1, 6, 6, 8}
669 };
670 static const int faceAndVertexToRightEdgePrism[5][6] = {
671 { 0, 1, -1, 6, 6, -1},
672 { 0, -1, 2, 7, -1, 7},
673 {-1, 1, 2, -1, 8, 8},
674 { 4, 5, 4, -1, -1, -1},
675 {-1, -1, -1, 7, 8, 7}
676 };
677 static const int faceAndVertexToLeftEdgeHex[6][8] = {
678 { 0, -1, 4, -1, 8, -1, 2, -1},
679 {-1, 5, -1, 3, -1, 1, -1, 9},
680 { 6, 1, -1, -1, 0, 10, -1, -1},
681 {-1, -1, 2, 7, -1, -1, 11, 3},
682 { 4, 6, 7, 5, -1, -1, -1, -1},
683 {-1, -1, -1, -1, 10, 9, 8, 11}
684 };
685 static const int faceAndVertexToRightEdgeHex[6][8] = {
686 { 4, -1, 2, -1, 0, -1, 8, -1},
687 {-1, 1, -1, 5, -1, 9, -1, 3},
688 { 0, 6, -1, -1, 10, 1, -1, -1},
689 {-1, -1, 7, 3, -1, -1, 2, 11},
690 { 6, 5, 4, 7, -1, -1, -1, -1},
691 {-1, -1, -1, -1, 8, 10, 11, 9}
692 };
693
694 switch (numElemVertices) {
695 case 4:
696 leftEdge = static_cast<unsigned>(faceAndVertexToLeftEdgeTet[face][vert]);
697 rightEdge = static_cast<unsigned>(faceAndVertexToRightEdgeTet[face][vert]);
698 break;
699 case 5:
700 leftEdge = static_cast<unsigned>(faceAndVertexToLeftEdgePyramid[face][vert]);
701 rightEdge = static_cast<unsigned>(faceAndVertexToRightEdgePyramid[face][vert]);
702 break;
703 case 6:
704 leftEdge = static_cast<unsigned>(faceAndVertexToLeftEdgePrism[face][vert]);
705 rightEdge = static_cast<unsigned>(faceAndVertexToRightEdgePrism[face][vert]);
706 break;
707 case 8:
708 leftEdge = static_cast<unsigned>(faceAndVertexToLeftEdgeHex[face][vert]);
709 rightEdge = static_cast<unsigned>(faceAndVertexToRightEdgeHex[face][vert]);
710 break;
711 default:
712 throw std::logic_error("Not implemented: VcfvStencil::getFaceIndices for "
713 +std::to_string(numElemVertices)+" vertices");
714 break;
715 }
716 }
717
718public:
720 using Mapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
721
723 {
724 public:
725 const GlobalPosition center() const
726 { return global(localGeometry_->center()); }
727
728 const GlobalPosition corner(unsigned cornerIdx) const
729 { return global(localGeometry_->corner(cornerIdx)); }
730
731 const GlobalPosition global(const LocalPosition& localPos) const
732 { return element_->geometry().global(localPos); }
733
735 { return *localGeometry_; }
736
738 const Element* element_;
739 };
740
742 {
743 const GlobalPosition& globalPos() const
744 { return global; }
745
746 const GlobalPosition center() const
747 { return geometry_.center(); }
748
749 Scalar volume() const
750 { return volume_; }
751
753 { return geometry_.localGeometry(); }
754
755 const ScvGeometry& geometry() const
756 { return geometry_; }
757
759 LocalPosition local;
760
762 GlobalPosition global;
763
765 Scalar volume_;
766
769
771 Dune::FieldVector<DimVector, maxNC> gradCenter;
772 };
773
775 {
776 const DimVector& normal() const
777 { return normal_; }
778
779 unsigned short interiorIndex() const
780 { return i; }
781
782 unsigned short exteriorIndex() const
783 { return j; }
784
785 Scalar area() const
786 { return area_; }
787
788 const LocalPosition& localPos() const
789 { return ipLocal_; }
790
791 const GlobalPosition& integrationPos() const
792 { return ipGlobal_; }
793
795 unsigned short i,j;
796
798 LocalPosition ipLocal_;
799
801 GlobalPosition ipGlobal_;
802
804 DimVector normal_;
805
807 Scalar area_;
808 };
809
812
813 VcfvStencil(const GridView& gridView, const Mapper& mapper)
814 : gridView_(gridView)
815 , vertexMapper_(mapper )
816 , element_(*gridView.template begin</*codim=*/0>())
817 {
818 // try to check if the mapper really maps the vertices
819 assert(static_cast<int>(gridView.size(/*codim=*/dimWorld)) == static_cast<int>(mapper.size()));
820
821 static bool localGeometriesInitialized = false;
822 if (!localGeometriesInitialized) {
823 localGeometriesInitialized = true;
824
825 VcfvScvGeometries<Scalar, /*dim=*/1, ElementType::cube>::init();
826 VcfvScvGeometries<Scalar, /*dim=*/2, ElementType::cube>::init();
827 VcfvScvGeometries<Scalar, /*dim=*/2, ElementType::simplex>::init();
828 VcfvScvGeometries<Scalar, /*dim=*/3, ElementType::cube>::init();
829 VcfvScvGeometries<Scalar, /*dim=*/3, ElementType::simplex>::init();
830 }
831 }
832
838 void updateTopology(const Element& e)
839 {
840 element_ = e;
841
842 numVertices = e.subEntities(/*codim=*/dim);
843 numEdges = e.subEntities(/*codim=*/dim-1);
844 if constexpr (dim == 3) {
845 numFaces = e.subEntities(/*codim=*/1);
846 }
847 else {
848 numFaces = 0;
849 }
850
851 numBoundarySegments_ = 0; // TODO: really required here(?)
852
853 // compute the local and global coordinates of the element
854 const Geometry& geometry = e.geometry();
855 geometryType_ = geometry.type();
856 const auto& referenceElement = Dune::ReferenceElements<CoordScalar,dim>::general(geometryType_);
857 // should be the same, but let's let the compiler know to quell (false) diagnostics
858 assert(numVertices == static_cast<unsigned>(referenceElement.size(dim)));
859 numVertices = std::min(numVertices, static_cast<unsigned>(referenceElement.size(dim)));
860 for (unsigned vertexIdx = 0; vertexIdx < numVertices; ++vertexIdx) {
861 subContVol[vertexIdx].local = referenceElement.position(static_cast<int>(vertexIdx), dim);
862 subContVol[vertexIdx].global = geometry.corner(static_cast<int>(vertexIdx));
863 }
864 }
865
866 void updatePrimaryTopology(const Element& element)
867 {
868 // since all degrees of freedom in a stencil are "primary" DOFs for the
869 // vertex-centered finite volume method, there's no difference to
870 // updateTopology()
871 updateTopology(element);
872 }
873
874 void update(const Element& e)
875 {
877
878 const Geometry& geometry = e.geometry();
879 geometryType_ = geometry.type();
880
881 const auto& referenceElement = Dune::ReferenceElements<CoordScalar,dim>::general(geometryType_);
882
883 elementVolume = geometry.volume();
884 elementLocal = referenceElement.position(0, 0);
885 elementGlobal = geometry.global(elementLocal);
886
887 // corners:
888 for (unsigned vert = 0; vert < numVertices; ++vert) {
889 subContVol[vert].local = referenceElement.position(static_cast<int>(vert), dim);
890 subContVol[vert].global = geometry.global(subContVol[vert].local);
891 }
892
893 // edges:
894 for (unsigned edge = 0; edge < numEdges; ++edge) {
895 edgeCoord[edge] = geometry.global(referenceElement.position(static_cast<int>(edge), dim - 1));
896 }
897
898 // faces:
899 for (unsigned face = 0; face < numFaces; ++face) {
900 faceCoord[face] = geometry.global(referenceElement.position(static_cast<int>(face), 1));
901 }
902
903 // fill sub control volume data use specialization for this
904 // \todo maybe it would be a good idea to migrate everything
905 // which is dependend of the grid's dimension to
906 // _VcfvFVElemGeomHelper in order to benefit from more aggressive
907 // compiler optimizations...
908 fillSubContVolData_();
909
910 // fill sub control volume face data:
911 for (unsigned k = 0; k < numEdges; ++k) { // begin loop over edges / sub control volume faces
912 unsigned short i =
913 static_cast<unsigned short>(referenceElement.subEntity(static_cast<int>(k), dim - 1, 0, dim));
914 unsigned short j =
915 static_cast<unsigned short>(referenceElement.subEntity(static_cast<int>(k), dim - 1, 1, dim));
916 if (numEdges == 4 && (i == 2 || j == 2)) {
917 std::swap(i, j);
918 }
919 subContVolFace[k].i = i;
920 subContVolFace[k].j = j;
921
922 // calculate the local integration point and
923 // the face normal. note that since dim is a
924 // constant which is known at compile time
925 // the compiler can optimize away all if
926 // cases which don't apply.
927 LocalPosition ipLocal;
928 DimVector diffVec;
929 if constexpr (dim == 1) {
930 subContVolFace[k].ipLocal_ = 0.5;
931 subContVolFace[k].normal_ = 1.0;
932 subContVolFace[k].area_ = 1.0;
933 ipLocal = subContVolFace[k].ipLocal_;
934 }
935 else if constexpr (dim == 2) {
936 ipLocal = referenceElement.position(static_cast<int>(k), dim - 1) + elementLocal;
937 ipLocal *= 0.5;
938 subContVolFace[k].ipLocal_ = ipLocal;
939 for (unsigned m = 0; m < dimWorld; ++m) {
940 diffVec[m] = elementGlobal[m] - edgeCoord[k][m];
941 }
942 subContVolFace[k].normal_[0] = diffVec[1];
943 subContVolFace[k].normal_[1] = -diffVec[0];
944
945 for (unsigned m = 0; m < dimWorld; ++m) {
946 diffVec[m] = subContVol[j].global[m] - subContVol[i].global[m];
947 }
948 // make sure the normal points to the right direction
949 if (subContVolFace[k].normal_ * diffVec < 0) {
950 subContVolFace[k].normal_ *= -1;
951 }
952
953 subContVolFace[k].area_ = subContVolFace[k].normal_.two_norm();
954 subContVolFace[k].normal_ /= subContVolFace[k].area_;
955 }
956 else if constexpr (dim == 3) {
957 unsigned leftFace;
958 unsigned rightFace;
959 getFaceIndices(numVertices, k, leftFace, rightFace);
960 ipLocal = referenceElement.position(static_cast<int>(k), dim - 1) +
961 elementLocal +
962 referenceElement.position(static_cast<int>(leftFace), 1) +
963 referenceElement.position(static_cast<int>(rightFace), 1);
964 ipLocal *= 0.25;
965 subContVolFace[k].ipLocal_ = ipLocal;
966 normalOfQuadrilateral3D(subContVolFace[k].normal_,
967 edgeCoord[k], faceCoord[rightFace],
968 elementGlobal, faceCoord[leftFace]);
969 subContVolFace[k].area_ = subContVolFace[k].normal_.two_norm();
970 subContVolFace[k].normal_ /= subContVolFace[k].area_;
971 }
972
973 // get the global integration point and the Jacobian inverse
974 subContVolFace[k].ipGlobal_ = geometry.global(ipLocal);
975 } // end loop over edges / sub control volume faces
976
977 // fill boundary face data:
978 for (const auto& intersection : intersections(gridView_, e)) {
979 if (!intersection.boundary()) {
980 continue;
981 }
982
983 const unsigned face = static_cast<unsigned>(intersection.indexInInside());
984 const unsigned numVerticesOfFace = static_cast<unsigned>(referenceElement.size(static_cast<int>(face), 1, dim));
985 for (unsigned vertInFace = 0; vertInFace < numVerticesOfFace; ++vertInFace) {
986 const unsigned short vertInElement = static_cast<unsigned short>(referenceElement.subEntity(static_cast<int>(face), 1, static_cast<int>(vertInFace), dim));
987 const unsigned bfIdx = numBoundarySegments_;
988 ++numBoundarySegments_;
989
990 if constexpr (dim == 1) {
991 boundaryFace_[bfIdx].ipLocal_ = referenceElement.position(static_cast<int>(vertInElement), dim);
992 boundaryFace_[bfIdx].area_ = 1.0;
993 }
994 else if constexpr (dim == 2) {
995 boundaryFace_[bfIdx].ipLocal_ =
996 referenceElement.position(static_cast<int>(vertInElement), dim) +
997 referenceElement.position(static_cast<int>(face), 1);
998 boundaryFace_[bfIdx].ipLocal_ *= 0.5;
999 boundaryFace_[bfIdx].area_ = 0.5 * intersection.geometry().volume();
1000 }
1001 else if constexpr (dim == 3) {
1002 unsigned leftEdge;
1003 unsigned rightEdge;
1004 getEdgeIndices(numVertices, face, vertInElement, leftEdge, rightEdge);
1005 boundaryFace_[bfIdx].ipLocal_ =
1006 referenceElement.position(static_cast<int>(vertInElement), dim) +
1007 referenceElement.position(static_cast<int>(face), 1) +
1008 referenceElement.position(static_cast<int>(leftEdge), dim - 1) +
1009 referenceElement.position(static_cast<int>(rightEdge), dim - 1);
1010 boundaryFace_[bfIdx].ipLocal_ *= 0.25;
1011 boundaryFace_[bfIdx].area_ =
1012 quadrilateralArea3D(subContVol[vertInElement].global,
1013 edgeCoord[rightEdge],
1014 faceCoord[face],
1015 edgeCoord[leftEdge]);
1016 }
1017 else {
1018 throw std::logic_error("Not implemented:VcfvStencil for dim = " + std::to_string(dim));
1019 }
1020
1021 boundaryFace_[bfIdx].ipGlobal_ = geometry.global(boundaryFace_[bfIdx].ipLocal_);
1022 boundaryFace_[bfIdx].i = vertInElement;
1023 boundaryFace_[bfIdx].j = vertInElement;
1024
1025 // ASSUME constant normal on the segment of the boundary face
1026 boundaryFace_[bfIdx].normal_ = intersection.centerUnitOuterNormal();
1027 }
1028 }
1029
1031 }
1032
1033 void updateScvGeometry(const Element& element)
1034 {
1035 auto geomType = element.geometry().type();
1036
1037 // get the local geometries of the sub control volumes
1038 if (geomType.isTriangle() || geomType.isTetrahedron()) {
1039 for (unsigned vertIdx = 0; vertIdx < numVertices; ++vertIdx) {
1040 subContVol[vertIdx].geometry_.element_ = &element;
1041 subContVol[vertIdx].geometry_.localGeometry_ =
1042 &VcfvScvGeometries<Scalar, dim, ElementType::simplex>::get(vertIdx);
1043 }
1044 }
1045 else if (geomType.isLine() || geomType.isQuadrilateral() || geomType.isHexahedron()) {
1046 for (unsigned vertIdx = 0; vertIdx < numVertices; ++vertIdx) {
1047 subContVol[vertIdx].geometry_.element_ = &element;
1048 subContVol[vertIdx].geometry_.localGeometry_ =
1049 &VcfvScvGeometries<Scalar, dim, ElementType::cube>::get(vertIdx);
1050 }
1051 }
1052 else {
1053 throw std::logic_error("Not implemented: SCV geometries for non hexahedron elements");
1054 }
1055 }
1056
1057#if HAVE_DUNE_LOCALFUNCTIONS
1058 void updateCenterGradients()
1059 {
1060 const auto& localFiniteElement = feCache_.get(element_.type());
1061 const auto& geom = element_.geometry();
1062
1063 std::vector<ShapeJacobian> localJac;
1064
1065 for (unsigned scvIdx = 0; scvIdx < numVertices; ++scvIdx) {
1066 const auto& localCenter = subContVol[scvIdx].localGeometry().center();
1067 localFiniteElement.localBasis().evaluateJacobian(localCenter, localJac);
1068 const auto& globalPos = subContVol[scvIdx].geometry().center();
1069
1070 const auto jacInvT = geom.jacobianInverseTransposed(globalPos);
1071 for (unsigned vert = 0; vert < numVertices; ++vert) {
1072 jacInvT.mv(localJac[vert][0], subContVol[scvIdx].gradCenter[vert]);
1073 }
1074 }
1075 }
1076#endif
1077
1078 unsigned numDof() const
1079 { return numVertices; }
1080
1081 unsigned numPrimaryDof() const
1082 { return numDof(); }
1083
1084 Dune::PartitionType partitionType(unsigned scvIdx) const
1085 { return element_.template subEntity</*codim=*/dim>(scvIdx)->partitionType(); }
1086
1087 const SubControlVolume& subControlVolume(unsigned scvIdx) const
1088 {
1089 assert(scvIdx < numDof());
1090 return subContVol[scvIdx];
1091 }
1092
1093 unsigned numInteriorFaces() const
1094 { return numEdges; }
1095
1096 unsigned numBoundaryFaces() const
1097 { return numBoundarySegments_; }
1098
1099 const SubControlVolumeFace& interiorFace(unsigned faceIdx) const
1100 { return subContVolFace[faceIdx]; }
1101
1102 const BoundaryFace& boundaryFace(unsigned bfIdx) const
1103 { return boundaryFace_[bfIdx]; }
1104
1109 unsigned globalSpaceIndex(unsigned dofIdx) const
1110 {
1111 assert(dofIdx < numDof());
1112
1113 return static_cast<unsigned>(vertexMapper_.subIndex(element_, static_cast<int>(dofIdx), /*codim=*/dim));
1114 }
1115
1120 Entity entity(unsigned dofIdx) const
1121 {
1122 assert(dofIdx < numDof());
1123 return element_.template subEntity<dim>(static_cast<int>(dofIdx));
1124 }
1125
1126private:
1127#if __GNUC__ || __clang__
1128#pragma GCC diagnostic push
1129#pragma GCC diagnostic ignored "-Wpragmas"
1130#pragma GCC diagnostic ignored "-Warray-bounds"
1131#endif
1132 void fillSubContVolData_()
1133 {
1134 if constexpr (dim == 1) {
1135 // 1D
1136 subContVol[0].volume_ = 0.5 * elementVolume;
1137 subContVol[1].volume_ = 0.5 * elementVolume;
1138 }
1139 else if constexpr (dim == 2) {
1140 switch (numVertices) {
1141 case 3: // 2D, triangle
1142 subContVol[0].volume_ = elementVolume / 3;
1143 subContVol[1].volume_ = elementVolume / 3;
1144 subContVol[2].volume_ = elementVolume / 3;
1145 break;
1146 case 4: // 2D, quadrilinear
1147 subContVol[0].volume_ =
1148 quadrilateralArea(subContVol[0].global,
1149 edgeCoord[2],
1150 elementGlobal,
1151 edgeCoord[0]);
1152 subContVol[1].volume_ =
1153 quadrilateralArea(subContVol[1].global,
1154 edgeCoord[1],
1155 elementGlobal,
1156 edgeCoord[2]);
1157 subContVol[2].volume_ =
1158 quadrilateralArea(subContVol[2].global,
1159 edgeCoord[0],
1160 elementGlobal,
1161 edgeCoord[3]);
1162 subContVol[3].volume_ =
1163 quadrilateralArea(subContVol[3].global,
1164 edgeCoord[3],
1165 elementGlobal,
1166 edgeCoord[1]);
1167 break;
1168 default:
1169 throw std::logic_error("Not implemented:VcfvStencil dim = " + std::to_string(dim) +
1170 ", numVertices = " + std::to_string(numVertices));
1171 }
1172 }
1173 else if constexpr (dim == 3) {
1174 switch (numVertices) {
1175 case 4: // 3D, tetrahedron
1176 for (unsigned k = 0; k < numVertices; k++) {
1177 subContVol[k].volume_ = elementVolume / 4.0;
1178 }
1179 break;
1180 case 5: // 3D, pyramid
1181 subContVol[0].volume_ =
1182 hexahedronVolume(subContVol[0].global,
1183 edgeCoord[2],
1184 faceCoord[0],
1185 edgeCoord[0],
1186 edgeCoord[4],
1187 faceCoord[3],
1188 elementGlobal,
1189 faceCoord[1]);
1190 subContVol[1].volume_ =
1191 hexahedronVolume(subContVol[1].global,
1192 edgeCoord[1],
1193 faceCoord[0],
1194 edgeCoord[2],
1195 edgeCoord[5],
1196 faceCoord[2],
1197 elementGlobal,
1198 faceCoord[3]);
1199 subContVol[2].volume_ =
1200 hexahedronVolume(subContVol[2].global,
1201 edgeCoord[0],
1202 faceCoord[0],
1203 edgeCoord[3],
1204 edgeCoord[6],
1205 faceCoord[1],
1206 elementGlobal,
1207 faceCoord[4]);
1208 subContVol[3].volume_ =
1209 hexahedronVolume(subContVol[3].global,
1210 edgeCoord[3],
1211 faceCoord[0],
1212 edgeCoord[1],
1213 edgeCoord[7],
1214 faceCoord[4],
1215 elementGlobal,
1216 faceCoord[2]);
1217 subContVol[4].volume_ = elementVolume -
1218 subContVol[0].volume_ -
1219 subContVol[1].volume_ -
1220 subContVol[2].volume_ -
1221 subContVol[3].volume_;
1222 break;
1223 case 6: // 3D, prism
1224 subContVol[0].volume_ =
1225 hexahedronVolume(subContVol[0].global,
1226 edgeCoord[3],
1227 faceCoord[3],
1228 edgeCoord[4],
1229 edgeCoord[0],
1230 faceCoord[0],
1231 elementGlobal,
1232 faceCoord[1]);
1233 subContVol[1].volume_ =
1234 hexahedronVolume(subContVol[1].global,
1235 edgeCoord[5],
1236 faceCoord[3],
1237 edgeCoord[3],
1238 edgeCoord[1],
1239 faceCoord[2],
1240 elementGlobal,
1241 faceCoord[0]);
1242 subContVol[2].volume_ =
1243 hexahedronVolume(subContVol[2].global,
1244 edgeCoord[4],
1245 faceCoord[3],
1246 edgeCoord[5],
1247 edgeCoord[2],
1248 faceCoord[1],
1249 elementGlobal,
1250 faceCoord[2]);
1251 subContVol[3].volume_ =
1252 hexahedronVolume(edgeCoord[0],
1253 faceCoord[0],
1254 elementGlobal,
1255 faceCoord[1],
1256 subContVol[3].global,
1257 edgeCoord[6],
1258 faceCoord[4],
1259 edgeCoord[7]);
1260 subContVol[4].volume_ =
1261 hexahedronVolume(edgeCoord[1],
1262 faceCoord[2],
1263 elementGlobal,
1264 faceCoord[0],
1265 subContVol[4].global,
1266 edgeCoord[8],
1267 faceCoord[4],
1268 edgeCoord[6]);
1269 subContVol[5].volume_ =
1270 hexahedronVolume(edgeCoord[2],
1271 faceCoord[1],
1272 elementGlobal,
1273 faceCoord[2],
1274 subContVol[5].global,
1275 edgeCoord[7],
1276 faceCoord[4],
1277 edgeCoord[8]);
1278 break;
1279 case 8: // 3D, hexahedron
1280 subContVol[0].volume_ =
1281 hexahedronVolume(subContVol[0].global,
1282 edgeCoord[6],
1283 faceCoord[4],
1284 edgeCoord[4],
1285 edgeCoord[0],
1286 faceCoord[2],
1287 elementGlobal,
1288 faceCoord[0]);
1289 subContVol[1].volume_ =
1290 hexahedronVolume(subContVol[1].global,
1291 edgeCoord[5],
1292 faceCoord[4],
1293 edgeCoord[6],
1294 edgeCoord[1],
1295 faceCoord[1],
1296 elementGlobal,
1297 faceCoord[2]);
1298 subContVol[2].volume_ =
1299 hexahedronVolume(subContVol[2].global,
1300 edgeCoord[4],
1301 faceCoord[4],
1302 edgeCoord[7],
1303 edgeCoord[2],
1304 faceCoord[0],
1305 elementGlobal,
1306 faceCoord[3]);
1307 subContVol[3].volume_ =
1308 hexahedronVolume(subContVol[3].global,
1309 edgeCoord[7],
1310 faceCoord[4],
1311 edgeCoord[5],
1312 edgeCoord[3],
1313 faceCoord[3],
1314 elementGlobal,
1315 faceCoord[1]);
1316 subContVol[4].volume_ =
1317 hexahedronVolume(edgeCoord[0],
1318 faceCoord[2],
1319 elementGlobal,
1320 faceCoord[0],
1321 subContVol[4].global,
1322 edgeCoord[10],
1323 faceCoord[5],
1324 edgeCoord[8]);
1325 subContVol[5].volume_ =
1326 hexahedronVolume(edgeCoord[1],
1327 faceCoord[1],
1328 elementGlobal,
1329 faceCoord[2],
1330 subContVol[5].global,
1331 edgeCoord[9],
1332 faceCoord[5],
1333 edgeCoord[10]);
1334 subContVol[6].volume_ =
1335 hexahedronVolume(edgeCoord[2],
1336 faceCoord[0],
1337 elementGlobal,
1338 faceCoord[3],
1339 subContVol[6].global,
1340 edgeCoord[8],
1341 faceCoord[5],
1342 edgeCoord[11]);
1343 subContVol[7].volume_ =
1344 hexahedronVolume(edgeCoord[3],
1345 faceCoord[3],
1346 elementGlobal,
1347 faceCoord[1],
1348 subContVol[7].global,
1349 edgeCoord[11],
1350 faceCoord[5],
1351 edgeCoord[9]);
1352 break;
1353 default:
1354 throw std::logic_error("Not implemented:VcfvStencil for dim = " + std::to_string(dim) +
1355 ", numVertices = " + std::to_string(numVertices));
1356 }
1357 }
1358 else {
1359 throw std::logic_error("Not implemented:VcfvStencil for dim = " + std::to_string(dim));
1360 }
1361 }
1362#if __GNUC__ || __clang__
1363#pragma GCC diagnostic pop
1364#endif
1365
1366 const GridView& gridView_;
1367 const Mapper& vertexMapper_;
1368
1369 Element element_;
1370
1371#if HAVE_DUNE_LOCALFUNCTIONS
1372 static LocalFiniteElementCache feCache_;
1373#endif // HAVE_DUNE_LOCALFUNCTIONS
1374
1376 LocalPosition elementLocal;
1377
1379 GlobalPosition elementGlobal;
1380
1382 Scalar elementVolume;
1383
1385 std::array<SubControlVolume, maxNC> subContVol{};
1386
1388 std::array<SubControlVolumeFace, maxNE> subContVolFace{};
1389
1391 std::array<BoundaryFace, maxBF> boundaryFace_{};
1392
1393 unsigned numBoundarySegments_{};
1394
1396 std::array<GlobalPosition, maxNE> edgeCoord{};
1397
1399 std::array<GlobalPosition, maxNF> faceCoord{};
1400
1402 unsigned numVertices{};
1403
1405 unsigned numEdges{};
1406
1408 unsigned numFaces{};
1409
1410 Dune::GeometryType geometryType_;
1411};
1412
1413#if HAVE_DUNE_LOCALFUNCTIONS
1414template<class Scalar, class GridView>
1415typename VcfvStencil<Scalar, GridView>::LocalFiniteElementCache
1416VcfvStencil<Scalar, GridView>::feCache_;
1417#endif // HAVE_DUNE_LOCALFUNCTIONS
1418
1419template <class Scalar>
1420std::array<typename VcfvScvGeometries<Scalar, /*dim=*/1, ElementType::cube>::ScvLocalGeometry,
1421 VcfvScvGeometries<Scalar, /*dim=*/1, ElementType::cube>::numScv>
1422VcfvScvGeometries<Scalar, /*dim=*/1, ElementType::cube>::scvGeoms_{};
1423
1424template <class Scalar>
1425std::array<typename VcfvScvGeometries<Scalar, /*dim=*/2, ElementType::simplex>::ScvLocalGeometry,
1426 VcfvScvGeometries<Scalar, /*dim=*/2, ElementType::simplex>::numScv>
1427VcfvScvGeometries<Scalar, /*dim=*/2, ElementType::simplex>::scvGeoms_{};
1428
1429template <class Scalar>
1430std::array<typename VcfvScvGeometries<Scalar, /*dim=*/2, ElementType::cube>::ScvLocalGeometry,
1431 VcfvScvGeometries<Scalar, /*dim=*/2, ElementType::cube>::numScv>
1432VcfvScvGeometries<Scalar, /*dim=*/2, ElementType::cube>::scvGeoms_{};
1433
1434template <class Scalar>
1435std::array<typename VcfvScvGeometries<Scalar, /*dim=*/3, ElementType::simplex>::ScvLocalGeometry,
1436 VcfvScvGeometries<Scalar, /*dim=*/3, ElementType::simplex>::numScv>
1437VcfvScvGeometries<Scalar, /*dim=*/3, ElementType::simplex>::scvGeoms_{};
1438
1439template <class Scalar>
1440std::array<typename VcfvScvGeometries<Scalar, /*dim=*/3, ElementType::cube>::ScvLocalGeometry,
1441 VcfvScvGeometries<Scalar, /*dim=*/3, ElementType::cube>::numScv>
1442VcfvScvGeometries<Scalar, /*dim=*/3, ElementType::cube>::scvGeoms_{};
1443
1444} // namespace Opm
1445
1446#endif
Quadrature geometry for quadrilaterals.
Definition: quadraturegeometries.hh:44
const GlobalPosition & corner(unsigned cornerIdx) const
Return the position of the corner with a given index.
Definition: quadraturegeometries.hh:133
const GlobalPosition & center() const
Returns the center of weight of the polyhedron.
Definition: quadraturegeometries.hh:73
Definition: vcfvstencil.hh:723
const GlobalPosition center() const
Definition: vcfvstencil.hh:725
const ScvLocalGeometry * localGeometry_
Definition: vcfvstencil.hh:737
const Element * element_
Definition: vcfvstencil.hh:738
const GlobalPosition corner(unsigned cornerIdx) const
Definition: vcfvstencil.hh:728
const GlobalPosition global(const LocalPosition &localPos) const
Definition: vcfvstencil.hh:731
const ScvLocalGeometry & localGeometry() const
Definition: vcfvstencil.hh:734
Represents the finite volume geometry of a single element in the VCFV discretization.
Definition: vcfvstencil.hh:452
unsigned numPrimaryDof() const
Definition: vcfvstencil.hh:1081
const SubControlVolume & subControlVolume(unsigned scvIdx) const
Definition: vcfvstencil.hh:1087
VcfvStencil(const GridView &gridView, const Mapper &mapper)
Definition: vcfvstencil.hh:813
unsigned numBoundaryFaces() const
Definition: vcfvstencil.hh:1096
unsigned globalSpaceIndex(unsigned dofIdx) const
Return the global space index given the index of a degree of freedom.
Definition: vcfvstencil.hh:1109
unsigned numDof() const
Definition: vcfvstencil.hh:1078
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > Mapper
exported Mapper type
Definition: vcfvstencil.hh:720
void update(const Element &e)
Definition: vcfvstencil.hh:874
unsigned numInteriorFaces() const
Definition: vcfvstencil.hh:1093
Entity entity(unsigned dofIdx) const
Return the global space index given the index of a degree of freedom.
Definition: vcfvstencil.hh:1120
void updatePrimaryTopology(const Element &element)
Definition: vcfvstencil.hh:866
typename GridView::Traits::template Codim< dim >::Entity Entity
Definition: vcfvstencil.hh:463
const SubControlVolumeFace & interiorFace(unsigned faceIdx) const
Definition: vcfvstencil.hh:1099
const BoundaryFace & boundaryFace(unsigned bfIdx) const
Definition: vcfvstencil.hh:1102
Dune::PartitionType partitionType(unsigned scvIdx) const
Definition: vcfvstencil.hh:1084
void updateTopology(const Element &e)
Update the non-geometric part of the stencil.
Definition: vcfvstencil.hh:838
void updateScvGeometry(const Element &element)
Definition: vcfvstencil.hh:1033
static constexpr int dim
Definition: structuredgridvanguard.hh:68
Definition: blackoilbioeffectsmodules.hh:45
ElementType
The types of reference elements available.
Definition: vcfvstencil.hh:57
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
interior face of a sub control volume
Definition: vcfvstencil.hh:775
Scalar area() const
Definition: vcfvstencil.hh:785
Scalar area_
area of face
Definition: vcfvstencil.hh:807
unsigned short j
Definition: vcfvstencil.hh:795
unsigned short exteriorIndex() const
Definition: vcfvstencil.hh:782
GlobalPosition ipGlobal_
integration point in global coords
Definition: vcfvstencil.hh:801
const GlobalPosition & integrationPos() const
Definition: vcfvstencil.hh:791
unsigned short interiorIndex() const
Definition: vcfvstencil.hh:779
DimVector normal_
normal on face pointing to CV j or outward of the domain with length equal to |scvf|
Definition: vcfvstencil.hh:804
unsigned short i
scvf seperates corner i and j of elem
Definition: vcfvstencil.hh:795
const DimVector & normal() const
Definition: vcfvstencil.hh:776
const LocalPosition & localPos() const
Definition: vcfvstencil.hh:788
LocalPosition ipLocal_
integration point in local coords
Definition: vcfvstencil.hh:798
finite volume intersected with element
Definition: vcfvstencil.hh:742
ScvGeometry geometry_
The geometry of the sub-control volume in local coordinates.
Definition: vcfvstencil.hh:768
const ScvGeometry & geometry() const
Definition: vcfvstencil.hh:755
LocalPosition local
local vertex position
Definition: vcfvstencil.hh:759
GlobalPosition global
global vertex position
Definition: vcfvstencil.hh:762
Scalar volume() const
Definition: vcfvstencil.hh:749
const GlobalPosition center() const
Definition: vcfvstencil.hh:746
Scalar volume_
space occupied by the sub-control volume
Definition: vcfvstencil.hh:765
const GlobalPosition & globalPos() const
Definition: vcfvstencil.hh:743
const ScvLocalGeometry & localGeometry() const
Definition: vcfvstencil.hh:752
Dune::FieldVector< DimVector, maxNC > gradCenter
derivative of shape function at the center of the sub control volume
Definition: vcfvstencil.hh:771