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