version 3.11-dev
discretization/box/fvelementgeometry.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-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_DISCRETIZATION_BOX_FV_ELEMENT_GEOMETRY_HH
15#define DUMUX_DISCRETIZATION_BOX_FV_ELEMENT_GEOMETRY_HH
16
17#include <optional>
18#include <utility>
19#include <unordered_map>
20#include <array>
21#include <vector>
22
23#include <dune/geometry/type.hh>
24#include <dune/localfunctions/lagrange/pqkfactory.hh>
25
32
33namespace Dumux {
34
43template<class GG, bool enableGridGeometryCache>
45
47template<class GG>
48class BoxFVElementGeometry<GG, true>
49{
50 using GridView = typename GG::GridView;
51 static constexpr int dim = GridView::dimension;
52 static constexpr int dimWorld = GridView::dimensionworld;
53 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
54 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
55 using CoordScalar = typename GridView::ctype;
56 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
57 using GGCache = typename GG::Cache;
58 using GeometryHelper = typename GGCache::GeometryHelper;
59
61 typename GridView::template Codim<0>::Entity::Geometry::LocalCoordinate,
62 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate
63 >;
64
65public:
67 using Element = typename GridView::template Codim<0>::Entity;
69 using SubControlVolume = typename GG::SubControlVolume;
71 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
73 using GridGeometry = GG;
75 using ScvQuadratureRule = typename GG::ScvQuadratureRule;
77 using ScvfQuadratureRule = typename GG::ScvfQuadratureRule;
79 static constexpr std::size_t maxNumElementScvs = (1<<dim);
80
85 BoxFVElementGeometry(const GGCache& ggCache)
86 : ggCache_(&ggCache) {}
87
89 const SubControlVolume& scv(LocalIndexType scvIdx) const
90 {
91 return ggCache_->scvs(eIdx_)[scvIdx];
92 }
93
95 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
96 {
97 return ggCache_->scvfs(eIdx_)[scvfIdx];
98 }
99
105 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
106 scvs(const BoxFVElementGeometry& fvGeometry)
107 {
108 using Iter = typename std::vector<SubControlVolume>::const_iterator;
109 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
110 return Dune::IteratorRange<Iter>(s.begin(), s.end());
111 }
112
118 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
119 scvfs(const BoxFVElementGeometry& fvGeometry)
120 {
121 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
122 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
123 return Dune::IteratorRange<Iter>(s.begin(), s.end());
124 }
125
127 template<class Intersection>
128 friend inline auto localDofs(const BoxFVElementGeometry& fvGeometry, const Intersection& intersection)
129 {
130 return Dune::transformedRangeView(
131 Dune::range(Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).size(intersection.indexInInside(), 1, dim)),
132 [&](const auto i)
133 {
134 auto localDofIdx = Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).subEntity(intersection.indexInInside(), 1, i, dim);
135 return CVFE::LocalDof
136 {
137 static_cast<LocalIndexType>(localDofIdx),
138 fvGeometry.gridGeometry().dofMapper().subIndex(fvGeometry.element(), localDofIdx, dim),
139 static_cast<GridIndexType>(fvGeometry.elementIndex())
140 };
141 }
142 );
143 }
144
146 const FeLocalBasis& feLocalBasis() const
147 {
148 return gridGeometry().feCache().get(element_->type()).localBasis();
149 }
150
152 std::size_t numLocalDofs() const
153 {
154 return numScv();
155 }
156
158 std::size_t numScv() const
159 {
160 return ggCache_->scvs(eIdx_).size();
161 }
162
164 std::size_t numScvf() const
165 {
166 return ggCache_->scvfs(eIdx_).size();
167 }
168
175 {
176 this->bindElement(element);
177 return std::move(*this);
178 }
179
183 void bind(const Element& element) &
184 { this->bindElement(element); }
185
192 {
193 this->bindElement(element);
194 return std::move(*this);
195 }
196
200 void bindElement(const Element& element) &
201 {
202 element_ = element;
203 elementGeometry_.emplace(element.geometry());
204 // cache element index
205 eIdx_ = gridGeometry().elementMapper().index(element);
206 }
207
209 bool isBound() const
210 { return static_cast<bool>(element_); }
211
213 const Element& element() const
214 { return *element_; }
215
217 const typename Element::Geometry& elementGeometry() const
218 { return *elementGeometry_; }
219
221 GridIndexType elementIndex() const
222 { return eIdx_; }
223
225 std::size_t intersectionIndex(const SubControlVolumeFace& scvf) const
226 {
227 const auto localScvfIdx = scvf.index() - GeometryHelper::numInteriorScvf(element().type());
228 return ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx][0];
229 }
230
233 { return ggCache_->gridGeometry(); }
234
236 bool hasBoundaryScvf() const
237 { return ggCache_->hasBoundaryScvf(eIdx_); }
238
240 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
241 {
242 assert(isBound());
243 return { Dune::GeometryTypes::cube(dim), GeometryHelper(*elementGeometry_).getScvCorners(scv.indexInElement()) };
244 }
245
247 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
248 {
249 assert(isBound());
250 const GeometryHelper geometryHelper(*elementGeometry_);
251 if (scvf.boundary())
252 {
253 const auto localBoundaryIndex = scvf.index() - geometryHelper.numInteriorScvf();
254 const auto& key = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localBoundaryIndex];
255 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getBoundaryScvfCorners(key[0], key[1]) };
256 }
257 else
258 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(scvf.index()) };
259 }
260
262 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const SubControlVolume& scv)
263 {
264 const auto type = fvGeometry.element().type();
265 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(scv.localDofIndex());
266
267 return CVFE::LocalDofInterpolationPointData{ GeometryHelper::localDofPosition(type, localKey), scv.dofPosition(), scv.localDofIndex() };
268 }
269
271 template<class LocalDof>
272 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const LocalDof& localDof)
273 {
274 const auto type = fvGeometry.element().type();
275 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(localDof.index());
276 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
277
278 return CVFE::LocalDofInterpolationPointData{ localPos, fvGeometry.elementGeometry().global(localPos), localDof.index() };
279 }
280
282 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const typename Element::Geometry::GlobalCoordinate& globalPos)
283 {
284 // Create ipData that does not automatically calculate the local position but only if it is called
286 [&] (const typename Element::Geometry::GlobalCoordinate& pos)
287 { return fvGeometry.elementGeometry().local(pos); },
288 globalPos
289 };
290 }
291
293 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const SubControlVolumeFace& scvf)
294 {
296 { scvf.unitOuterNormal(), scvf.index(), fvGeometry.elementGeometry().local(scvf.ipGlobal()), scvf.ipGlobal() };
297 }
298
299private:
300 const GGCache* ggCache_;
301 GridIndexType eIdx_;
302
303 std::optional<Element> element_;
304 std::optional<typename Element::Geometry> elementGeometry_;
305};
306
308template<class GG>
309class BoxFVElementGeometry<GG, false>
310{
311 using GridView = typename GG::GridView;
312 static constexpr int dim = GridView::dimension;
313 static constexpr int dimWorld = GridView::dimensionworld;
314 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
315 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
316 using CoordScalar = typename GridView::ctype;
317 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
318 using GGCache = typename GG::Cache;
319 using GeometryHelper = typename GGCache::GeometryHelper;
320
322 typename GridView::template Codim<0>::Entity::Geometry::LocalCoordinate,
323 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate
324 >;
325
326public:
328 using Element = typename GridView::template Codim<0>::Entity;
330 using SubControlVolume = typename GG::SubControlVolume;
332 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
334 using GridGeometry = GG;
336 using ScvQuadratureRule = typename GG::ScvQuadratureRule;
338 using ScvfQuadratureRule = typename GG::ScvfQuadratureRule;
340 static constexpr std::size_t maxNumElementScvs = (1<<dim);
341
346 BoxFVElementGeometry(const GGCache& ggCache)
347 : ggCache_(&ggCache) {}
348
350 const SubControlVolume& scv(LocalIndexType scvIdx) const
351 {
352 return scvs_[scvIdx];
353 }
354
356 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
357 {
358 return scvfs_[scvfIdx];
359 }
360
366 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
367 scvs(const BoxFVElementGeometry& fvGeometry)
368 {
369 using Iter = typename std::vector<SubControlVolume>::const_iterator;
370 return Dune::IteratorRange<Iter>(fvGeometry.scvs_.begin(), fvGeometry.scvs_.end());
371 }
372
378 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
379 scvfs(const BoxFVElementGeometry& fvGeometry)
380 {
381 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
382 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
383 }
384
386 template<class Intersection>
387 friend inline auto localDofs(const BoxFVElementGeometry& fvGeometry, const Intersection& intersection)
388 {
389 return Dune::transformedRangeView(
390 Dune::range(Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).size(intersection.indexInInside(), 1, dim)),
391 [&](const auto i)
392 {
393 auto localDofIdx = Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).subEntity(intersection.indexInInside(), 1, i, dim);
394 return CVFE::LocalDof
395 {
396 static_cast<LocalIndexType>(localDofIdx),
397 fvGeometry.gridGeometry().dofMapper().subIndex(fvGeometry.element(), localDofIdx, dim),
398 static_cast<GridIndexType>(fvGeometry.elementIndex())
399 };
400 }
401 );
402 }
403
404
406 const FeLocalBasis& feLocalBasis() const
407 {
408 return gridGeometry().feCache().get(element_->type()).localBasis();
409 }
410
412 std::size_t numLocalDofs() const
413 {
414 return numScv();
415 }
416
418 std::size_t numScv() const
419 {
420 return scvs_.size();
421 }
422
424 std::size_t numScvf() const
425 {
426 return scvfs_.size();
427 }
428
435 {
436 this->bindElement(element);
437 return std::move(*this);
438 }
439
443 void bind(const Element& element) &
444 { this->bindElement(element); }
445
452 {
453 this->bindElement(element);
454 return std::move(*this);
455 }
456
460 void bindElement(const Element& element) &
461 {
462 element_ = element;
463 eIdx_ = gridGeometry().elementMapper().index(element);
464 elementGeometry_.emplace(element.geometry());
465 makeElementGeometries_();
466 }
467
469 bool isBound() const
470 { return static_cast<bool>(element_); }
471
473 const Element& element() const
474 { return *element_; }
475
477 const typename Element::Geometry& elementGeometry() const
478 { return *elementGeometry_; }
479
481 GridIndexType elementIndex() const
482 { return eIdx_; }
483
485 std::size_t intersectionIndex(const SubControlVolumeFace& scvf) const
486 {
487 const auto localScvfIdx = scvf.index() - GeometryHelper::numInteriorScvf(element().type());
488 return scvfBoundaryGeometryKeys_[localScvfIdx][0];
489 }
490
493 { return ggCache_->gridGeometry(); }
494
496 bool hasBoundaryScvf() const
497 { return hasBoundaryScvf_; }
498
500 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
501 {
502 assert(isBound());
503 return { Dune::GeometryTypes::cube(dim), GeometryHelper(*elementGeometry_).getScvCorners(scv.indexInElement()) };
504 }
505
507 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
508 {
509 assert(isBound());
510 const GeometryHelper geometryHelper(*elementGeometry_);
511 if (scvf.boundary())
512 {
513 const auto localBoundaryIndex = scvf.index() - geometryHelper.numInteriorScvf();
514 const auto& key = scvfBoundaryGeometryKeys_[localBoundaryIndex];
515 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getBoundaryScvfCorners(key[0], key[1]) };
516 }
517 else
518 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(scvf.index()) };
519 }
520
522 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const SubControlVolume& scv)
523 {
524 const auto type = fvGeometry.element().type();
525 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(scv.localDofIndex());
526
527 return CVFE::LocalDofInterpolationPointData{ GeometryHelper::localDofPosition(type, localKey), scv.dofPosition(), scv.localDofIndex() };
528 }
529
531 template<class LocalDof>
532 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const LocalDof& localDof)
533 {
534 const auto type = fvGeometry.element().type();
535 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(localDof.index());
536 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
537
538 return CVFE::LocalDofInterpolationPointData{ localPos, fvGeometry.elementGeometry().global(localPos), localDof.index() };
539 }
540
542 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const typename Element::Geometry::GlobalCoordinate& globalPos)
543 {
544 // Create ipData that does not automatically calculate the local position but only if it is called
546 [&] (const typename Element::Geometry::GlobalCoordinate& pos) { return fvGeometry.elementGeometry().local(pos); },
547 globalPos
548 };
549 }
550
552 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const SubControlVolumeFace& scvf)
553 {
555 { scvf.unitOuterNormal(), scvf.index(), fvGeometry.elementGeometry().local(scvf.ipGlobal()), scvf.ipGlobal() };
556 }
557
558private:
559 void makeElementGeometries_()
560 {
561 hasBoundaryScvf_ = false;
562
563 // get the element geometry
564 const auto& element = *element_;
565 const auto& elementGeometry = *elementGeometry_;
566 const auto refElement = referenceElement(elementGeometry);
567
568 // get the sub control volume geometries of this element
569 GeometryHelper geometryHelper(elementGeometry);
570
571 // construct the sub control volumes
572 scvs_.resize(elementGeometry.corners());
573 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
574 {
575 // get associated dof index
576 const auto dofIdxGlobal = gridGeometry().vertexMapper().subIndex(element, scvLocalIdx, dim);
577
578 // add scv to the local container
579 scvs_[scvLocalIdx] = SubControlVolume(
580 geometryHelper.getScvCorners(scvLocalIdx),
581 scvLocalIdx,
582 eIdx_,
583 dofIdxGlobal
584 );
585 }
586
587 // construct the sub control volume faces
588 const auto numInnerScvf = geometryHelper.numInteriorScvf();
589 scvfs_.resize(numInnerScvf);
590 scvfBoundaryGeometryKeys_.clear();
591
592 LocalIndexType scvfLocalIdx = 0;
593 for (; scvfLocalIdx < numInnerScvf; ++scvfLocalIdx)
594 {
595 // find the local scv indices this scvf is connected to
596 std::array<LocalIndexType, 2> localScvIndices{{
597 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
598 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))
599 }};
600
601 const auto& corners = geometryHelper.getScvfCorners(scvfLocalIdx);
602 scvfs_[scvfLocalIdx] = SubControlVolumeFace(
603 corners,
604 geometryHelper.normal(corners, localScvIndices),
605 element,
606 scvfLocalIdx,
607 std::move(localScvIndices)
608 );
609 }
610
611 // construct the sub control volume faces on the domain boundary
612 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
613 {
614 if (intersection.boundary() && !intersection.neighbor())
615 {
616 const auto isGeometry = intersection.geometry();
617 hasBoundaryScvf_ = true;
618
619 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < isGeometry.corners(); ++isScvfLocalIdx)
620 {
621 // find the scv this scvf is connected to
622 const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(refElement.subEntity(intersection.indexInInside(), 1, isScvfLocalIdx, dim));
623 std::array<LocalIndexType, 2> localScvIndices{{insideScvIdx, insideScvIdx}};
624
625 scvfs_.emplace_back(
626 geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), isScvfLocalIdx),
627 intersection.centerUnitOuterNormal(),
628 intersection,
629 isScvfLocalIdx,
630 scvfLocalIdx,
631 std::move(localScvIndices)
632 );
633
634 scvfBoundaryGeometryKeys_.emplace_back(std::array<LocalIndexType, 2>{{
635 static_cast<LocalIndexType>(intersection.indexInInside()),
636 static_cast<LocalIndexType>(isScvfLocalIdx)
637 }});
638
639 // increment local counter
640 scvfLocalIdx++;
641 }
642 }
643 }
644 }
645
647 GridIndexType eIdx_;
648 std::optional<Element> element_;
649 std::optional<typename Element::Geometry> elementGeometry_;
650
652 const GGCache* ggCache_;
653
655 std::vector<SubControlVolume> scvs_;
656 std::vector<SubControlVolumeFace> scvfs_;
657 std::vector<std::array<LocalIndexType, 2>> scvfBoundaryGeometryKeys_;
658
659 bool hasBoundaryScvf_ = false;
660};
661
662} // end namespace Dumux
663
664#endif
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/box/fvelementgeometry.hh:469
friend auto localDofs(const BoxFVElementGeometry &fvGeometry, const Intersection &intersection)
an iterator over all local dofs related to an intersection
Definition: discretization/box/fvelementgeometry.hh:387
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/box/fvelementgeometry.hh:424
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/box/fvelementgeometry.hh:334
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/box/fvelementgeometry.hh:507
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/box/fvelementgeometry.hh:496
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:367
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvelementgeometry.hh:418
GridIndexType elementIndex() const
The bound element's index in the grid view.
Definition: discretization/box/fvelementgeometry.hh:481
void bindElement(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:460
std::size_t intersectionIndex(const SubControlVolumeFace &scvf) const
The intersection index the scvf belongs to.
Definition: discretization/box/fvelementgeometry.hh:485
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/box/fvelementgeometry.hh:356
BoxFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/box/fvelementgeometry.hh:346
typename GG::ScvfQuadratureRule ScvfQuadratureRule
the quadrature rule type for scvfs
Definition: discretization/box/fvelementgeometry.hh:338
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/box/fvelementgeometry.hh:500
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/box/fvelementgeometry.hh:492
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/box/fvelementgeometry.hh:328
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:379
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition: discretization/box/fvelementgeometry.hh:522
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/box/fvelementgeometry.hh:406
const Element & element() const
The bound element.
Definition: discretization/box/fvelementgeometry.hh:473
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/box/fvelementgeometry.hh:350
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const typename Element::Geometry::GlobalCoordinate &globalPos)
Interpolation point data for a global position.
Definition: discretization/box/fvelementgeometry.hh:542
BoxFVElementGeometry bind(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/box/fvelementgeometry.hh:434
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf)
Interpolation point data for scvf.
Definition: discretization/box/fvelementgeometry.hh:552
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition: discretization/box/fvelementgeometry.hh:532
typename GG::ScvQuadratureRule ScvQuadratureRule
the quadrature rule type for scvs
Definition: discretization/box/fvelementgeometry.hh:336
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition: discretization/box/fvelementgeometry.hh:412
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/box/fvelementgeometry.hh:332
BoxFVElementGeometry bindElement(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/box/fvelementgeometry.hh:451
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition: discretization/box/fvelementgeometry.hh:477
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/box/fvelementgeometry.hh:330
void bind(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:443
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/box/fvelementgeometry.hh:209
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/box/fvelementgeometry.hh:146
friend auto localDofs(const BoxFVElementGeometry &fvGeometry, const Intersection &intersection)
an iterator over all local dofs related to an intersection
Definition: discretization/box/fvelementgeometry.hh:128
BoxFVElementGeometry bindElement(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/box/fvelementgeometry.hh:191
const Element & element() const
The bound element.
Definition: discretization/box/fvelementgeometry.hh:213
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/box/fvelementgeometry.hh:240
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:106
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvelementgeometry.hh:158
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/box/fvelementgeometry.hh:89
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition: discretization/box/fvelementgeometry.hh:152
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/box/fvelementgeometry.hh:164
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/box/fvelementgeometry.hh:69
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/box/fvelementgeometry.hh:232
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/box/fvelementgeometry.hh:236
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:119
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition: discretization/box/fvelementgeometry.hh:262
typename GG::ScvQuadratureRule ScvQuadratureRule
export the scv interpolation point data type
Definition: discretization/box/fvelementgeometry.hh:75
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition: discretization/box/fvelementgeometry.hh:217
void bind(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:183
typename GG::ScvfQuadratureRule ScvfQuadratureRule
the quadrature rule type for scvfs
Definition: discretization/box/fvelementgeometry.hh:77
std::size_t intersectionIndex(const SubControlVolumeFace &scvf) const
The intersection index the scvf belongs to.
Definition: discretization/box/fvelementgeometry.hh:225
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/box/fvelementgeometry.hh:247
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const typename Element::Geometry::GlobalCoordinate &globalPos)
Interpolation point data for a global position.
Definition: discretization/box/fvelementgeometry.hh:282
GridIndexType elementIndex() const
The bound element's index in the grid view.
Definition: discretization/box/fvelementgeometry.hh:221
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf)
Interpolation point data for scvf.
Definition: discretization/box/fvelementgeometry.hh:293
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition: discretization/box/fvelementgeometry.hh:272
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/box/fvelementgeometry.hh:71
void bindElement(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:200
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/box/fvelementgeometry.hh:95
BoxFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/box/fvelementgeometry.hh:85
BoxFVElementGeometry bind(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/box/fvelementgeometry.hh:174
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/box/fvelementgeometry.hh:67
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/box/fvelementgeometry.hh:73
Base class for the finite volume geometry vector for box models This builds up the sub control volume...
Definition: discretization/box/fvelementgeometry.hh:44
An interpolation point related to a face of an element.
Definition: cvfe/interpolationpointdata.hh:103
const GlobalPosition & unitOuterNormal() const
The unit outer normal vector at the quadrature point.
Definition: cvfe/interpolationpointdata.hh:117
An interpolation point related to an element that includes global and local positions.
Definition: cvfe/interpolationpointdata.hh:25
An interpolation point related to a global position of an element, giving its local positions by a ma...
Definition: cvfe/interpolationpointdata.hh:76
An interpolation point related to a localDof of an element, giving its global and local positions.
Definition: cvfe/interpolationpointdata.hh:54
LocalIndex localDofIndex() const
The local index of the corresponding dof.
Definition: cvfe/interpolationpointdata.hh:63
Classes representing interpolation point data for control-volume finite element schemes.
Defines the index types used for grid and local indices.
Class representing dofs on elements for control-volume finite element schemes.
Definition: adapt.hh:17
Quadrature rules over sub-control volumes and sub-control volume faces.
Class providing iterators over sub control volumes and sub control volume faces of an element.
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:27
unsigned int LocalIndex
Definition: indextraits.hh:28