29#ifndef DUMUX_POROUSMEDIUMFLOW_BOXDFM_FV_ELEMENT_GEOMETRY_HH
30#define DUMUX_POROUSMEDIUMFLOW_BOXDFM_FV_ELEMENT_GEOMETRY_HH
35#include <dune/geometry/type.hh>
36#include <dune/localfunctions/lagrange/pqkfactory.hh>
52template<
class GG,
bool enableGr
idGeometryCache>
59 using GridView =
typename GG::GridView;
60 static constexpr int dim = GridView::dimension;
61 static constexpr int dimWorld = GridView::dimensionworld;
62 using GridIndexType =
typename GridView::IndexSet::IndexType;
63 using CoordScalar =
typename GridView::ctype;
64 using FeLocalBasis =
typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
67 using Element =
typename GridView::template Codim<0>::Entity;
77 static constexpr std::size_t maxNumElementScvs = (1<<dim)*3;
81 : gridGeometryPtr_(&gridGeometry) {}
85 {
return gridGeometry().scvs(eIdx_)[scvIdx]; }
89 {
return gridGeometry().scvfs(eIdx_)[scvfIdx]; }
99 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
102 const auto& g = fvGeometry.gridGeometry();
103 using Iter =
typename std::vector<SubControlVolume>::const_iterator;
104 return Dune::IteratorRange<Iter>(g.scvs(fvGeometry.eIdx_).begin(), g.scvs(fvGeometry.eIdx_).end());
115 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
118 const auto& g = fvGeometry.gridGeometry();
119 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
120 return Dune::IteratorRange<Iter>(g.scvfs(fvGeometry.eIdx_).begin(), g.scvfs(fvGeometry.eIdx_).end());
125 {
return gridGeometry().feCache().get(element_->type()).localBasis(); }
129 {
return gridGeometry().scvs(eIdx_).size(); }
133 {
return gridGeometry().scvfs(eIdx_).size(); }
142 this->bindElement(element);
143 return std::move(*
this);
149 { this->bindElement(element); }
158 this->bindElement(element);
159 return std::move(*
this);
171 eIdx_ = gridGeometry().elementMapper().index(element);
176 {
return *gridGeometryPtr_; }
180 {
return static_cast<bool>(element_); }
184 {
return *element_; }
187 const GridGeometry* gridGeometryPtr_;
189 std::optional<Element> element_;
197 using GridView =
typename GG::GridView;
198 static constexpr int dim = GridView::dimension;
199 static constexpr int dimWorld = GridView::dimensionworld;
201 using GridIndexType =
typename GridView::IndexSet::IndexType;
203 using CoordScalar =
typename GridView::ctype;
204 using FeLocalBasis =
typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
206 typename GG::SubControlVolume,
207 typename GG::SubControlVolumeFace>;
210 using Element =
typename GridView::template Codim<0>::Entity;
219 static constexpr std::size_t maxNumElementScvs = (1<<dim)*3;
223 : gridGeometryPtr_(&gridGeometry) {}
227 {
return scvs_[scvIdx]; }
231 {
return scvfs_[scvfIdx]; }
241 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
244 using Iter =
typename std::vector<SubControlVolume>::const_iterator;
245 return Dune::IteratorRange<Iter>(fvGeometry.scvs_.begin(), fvGeometry.scvs_.end());
256 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
259 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
260 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
265 {
return gridGeometry().feCache().get(element_->type()).localBasis(); }
269 {
return scvs_.size(); }
273 {
return scvfs_.size(); }
282 this->bindElement(element);
283 return std::move(*
this);
293 { this->bindElement(element); }
302 this->bindElement(element);
303 return std::move(*
this);
313 eIdx_ = gridGeometry().elementMapper().index(element);
314 makeElementGeometries_();
319 {
return *gridGeometryPtr_; }
323 {
return static_cast<bool>(element_); }
327 {
return *element_; }
331 void makeElementGeometries_()
334 const auto& element = *element_;
335 const auto elementGeometry = element.geometry();
336 const auto refElement = referenceElement(element);
339 GeometryHelper geometryHelper(elementGeometry);
342 scvs_.resize(elementGeometry.corners());
343 using LocalIndexType =
typename SubControlVolumeFace::Traits::LocalIndexType;
344 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
347 const auto dofIdxGlobal = gridGeometry().vertexMapper().subIndex(element, scvLocalIdx, dim);
350 scvs_[scvLocalIdx] = SubControlVolume(geometryHelper,
357 const auto numInnerScvf = (dim==1) ? 1 : element.subEntities(dim-1);
358 scvfs_.resize(numInnerScvf);
360 unsigned int scvfLocalIdx = 0;
361 for (; scvfLocalIdx < numInnerScvf; ++scvfLocalIdx)
364 std::vector<LocalIndexType> localScvIndices({
static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
365 static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
367 scvfs_[scvfLocalIdx] = SubControlVolumeFace(geometryHelper,
371 std::move(localScvIndices));
385 LocalIndexType scvLocalIdx =
element.subEntities(dim);
386 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
389 const auto& isGeometry = intersection.geometry();
391 const auto idxInInside = intersection.indexInInside();
393 std::vector<GridIndexType> isVertexIndices(
numCorners);
394 for (
unsigned int vIdxLocal = 0; vIdxLocal <
numCorners; ++vIdxLocal)
395 isVertexIndices[vIdxLocal] = gridGeometry().vertexMapper().subIndex(element,
396 refElement.subEntity(idxInInside, 1, vIdxLocal, dim),
399 if (intersection.boundary())
401 for (
unsigned int isScvfLocalIdx = 0; isScvfLocalIdx <
numCorners; ++isScvfLocalIdx)
404 const LocalIndexType insideScvIdx =
static_cast<LocalIndexType
>(refElement.subEntity(idxInInside, 1, isScvfLocalIdx, dim));
405 std::vector<LocalIndexType> localScvIndices = {insideScvIdx, insideScvIdx};
407 scvfs_.emplace_back(geometryHelper,
412 std::move(localScvIndices));
420 if (this->gridGeometry().isOnFracture(element, intersection))
423 const auto curNumScvs = scvs_.size();
425 for (
unsigned int vIdxLocal = 0; vIdxLocal <
numCorners; ++vIdxLocal)
426 scvs_.emplace_back(geometryHelper,
430 static_cast<LocalIndexType
>(refElement.subEntity(idxInInside, 1, vIdxLocal, dim)),
434 isVertexIndices[vIdxLocal]);
439 const auto& faceRefElement = referenceElement(isGeometry);
440 for (
unsigned int edgeIdx = 0; edgeIdx < faceRefElement.size(1); ++edgeIdx)
443 std::vector<LocalIndexType> localScvIndices({
static_cast<LocalIndexType
>(faceRefElement.subEntity(edgeIdx, 1, 0, dim-1)),
444 static_cast<LocalIndexType
>(faceRefElement.subEntity(edgeIdx, 1, 1, dim-1))});
447 std::for_each( localScvIndices.begin(),
448 localScvIndices.end(),
449 [curNumScvs] (
auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
452 scvfs_.emplace_back(geometryHelper,
457 std::move(localScvIndices),
458 intersection.boundary());
466 std::vector<LocalIndexType> localScvIndices({0, 1});
469 std::for_each( localScvIndices.begin(),
470 localScvIndices.end(),
471 [curNumScvs] (
auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
474 scvfs_.emplace_back(geometryHelper,
479 std::move(localScvIndices),
480 intersection.boundary());
487 std::optional<Element> element_;
491 const GridGeometry* gridGeometryPtr_;
494 std::vector<SubControlVolume> scvs_;
495 std::vector<SubControlVolumeFace> scvfs_;
Class providing iterators over sub control volumes and sub control volume faces of an element.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:232
Base class for the finite volume geometry vector for box discrete fracture model.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:53
typename GG::SubControlVolumeFace SubControlVolumeFace
Export type of subcontrol volume face.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:71
void bindElement(const Element &element) &
Binding of an element, has to be called before using the fvgeometries.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:168
std::size_t numScv() const
The total number of sub control volumes.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:128
const SubControlVolumeFace & scvf(std::size_t scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:88
GG GridGeometry
Export type of finite volume grid geometry.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:73
const Element & element() const
The bound element.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:183
void bind(const Element &element) &
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:148
BoxDfmFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:80
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:132
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:124
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:67
typename GG::SubControlVolume SubControlVolume
Export type of subcontrol volume.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:69
const SubControlVolume & scv(std::size_t scvIdx) const
Get a sub control volume with a local scv index.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:84
BoxDfmFVElementGeometry 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: porousmediumflow/boxdfm/fvelementgeometry.hh:140
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxDfmFVElementGeometry &fvGeometry)
Iterator range for sub control volumes.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:100
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:175
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:179
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxDfmFVElementGeometry &fvGeometry)
Iterator range for sub control volumes faces.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:116
BoxDfmFVElementGeometry 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: porousmediumflow/boxdfm/fvelementgeometry.hh:156
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:264
const SubControlVolumeFace & scvf(std::size_t scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:230
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:318
BoxDfmFVElementGeometry 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: porousmediumflow/boxdfm/fvelementgeometry.hh:300
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:210
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:272
void bind(const Element &element) &
Binding of an element, has to be called before using the fvgeometries Prepares all the volume variabl...
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:292
BoxDfmFVElementGeometry 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: porousmediumflow/boxdfm/fvelementgeometry.hh:280
const Element & element() const
The bound element.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:326
BoxDfmFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:222
std::size_t numScv() const
The total number of sub control volumes.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:268
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:322
const SubControlVolume & scv(std::size_t scvIdx) const
Get a sub control volume with a local scv index.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:226
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxDfmFVElementGeometry &fvGeometry)
Iterator range for sub control volumes.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:242
GG GridGeometry
Export type of finite volume grid geometry.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:216
typename GG::SubControlVolumeFace SubControlVolumeFace
Export type of subcontrol volume face.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:214
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxDfmFVElementGeometry &fvGeometry)
Iterator range for sub control volumes faces.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:257
typename GG::SubControlVolume SubControlVolume
Export type of subcontrol volume.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:212
void bindElement(const Element &element) &
Binding of an element, has to be called before using the fvgeometries Prepares all the volume variabl...
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:310
Create sub control volumes and sub control volume face geometries.
Definition: porousmediumflow/boxdfm/geometryhelper.hh:59
Helper class constructing the dual grid finite volume geometries for the box discrete fracture model.