28#ifndef EWOMS_ECFV_STENCIL_HH
29#define EWOMS_ECFV_STENCIL_HH
31#include <dune/common/fvector.hh>
32#include <dune/common/version.hh>
34#include <dune/geometry/type.hh>
36#include <dune/grid/common/mcmgmapper.hh>
37#include <dune/grid/common/intersectioniterator.hh>
39#include <opm/common/ErrorMacros.hpp>
41#include <opm/input/eclipse/EclipseState/Grid/FaceDir.hpp>
43#include <opm/material/common/ConditionalStorage.hpp>
64template <
class Scalar,
66 bool needFaceIntegrationPos =
true,
67 bool needFaceNormal =
true>
70 enum { dimWorld = GridView::dimensionworld };
72 using CoordScalar =
typename GridView::ctype;
73 using Intersection =
typename GridView::Intersection;
74 using Element =
typename GridView::template Codim<0>::Entity;
76 using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
78 using GlobalPosition = Dune::FieldVector<CoordScalar, dimWorld>;
80 using WorldVector = Dune::FieldVector<Scalar, dimWorld>;
83 using Entity = Element;
84 using Mapper = ElementMapper;
86 using LocalGeometry =
typename Element::Geometry;
96 class SubControlVolume
104 explicit SubControlVolume(
const Element&
element)
108 void update(
const Element&
element)
118 {
return element_.geometry().center(); }
124 {
return element_.geometry().center(); }
130 {
return element_.geometry().volume(); }
136 {
return element_.geometry(); }
142 {
return element_.geometryInFather(); }
151 template <
bool needIntegrationPos,
bool needNormal>
152 class EcfvSubControlVolumeFace
155 EcfvSubControlVolumeFace()
158 EcfvSubControlVolumeFace(
const Intersection& intersection,
unsigned localNeighborIdx)
160 exteriorIdx_ =
static_cast<unsigned short>(localNeighborIdx);
162 if constexpr (needNormal) {
163 (*normal_) = intersection.centerUnitOuterNormal();
165 const auto& geometry = intersection.geometry();
166 if constexpr (needIntegrationPos) {
167 (*integrationPos_) = geometry.center();
169 area_ = geometry.volume();
174 dirId_ = intersection.indexInInside();
196 {
return exteriorIdx_; }
203 {
return *integrationPos_; }
232 using Dir = FaceDir::DirEnum;
234 OPM_THROW(std::runtime_error,
"NNC faces does not have a face id");
247 OPM_THROW(std::runtime_error,
248 "Unexpected face id" + std::to_string(dirId_));
253 ConditionalStorage<needIntegrationPos, GlobalPosition> integrationPos_;
254 ConditionalStorage<needNormal, WorldVector> normal_;
258 unsigned short exteriorIdx_;
261 using SubControlVolumeFace = EcfvSubControlVolumeFace<needFaceIntegrationPos, needFaceNormal>;
262 using BoundaryFace = EcfvSubControlVolumeFace<
true, needFaceNormal>;
266 , elementMapper_(mapper)
269 assert(
int(
gridView.size(0)) ==
int(elementMapper_.size()));
272 void updateTopology(
const Element&
element)
275 subControlVolumes_.clear();
276 subControlVolumes_.emplace_back(
element);
278 elements_.emplace_back(
element);
280 interiorFaces_.clear();
281 boundaryFaces_.clear();
283 for (
const auto& intersection : intersections(gridView_,
element)) {
287 if (intersection.neighbor()) {
288 elements_.emplace_back(intersection.outside());
289 subControlVolumes_.emplace_back(elements_.back());
290 interiorFaces_.emplace_back(intersection, subControlVolumes_.size() - 1);
293 boundaryFaces_.emplace_back(intersection, - 10000);
298 void updatePrimaryTopology(
const Element&
element)
301 subControlVolumes_.clear();
302 subControlVolumes_.emplace_back(
element);
304 elements_.emplace_back(
element);
307 void update(
const Element&
element)
312 void updateCenterGradients()
328 {
return *gridView_; }
335 {
return subControlVolumes_.size(); }
355 assert(dofIdx <
numDof());
357 return static_cast<unsigned>(elementMapper_.index(
element(dofIdx)));
364 {
return elements_[dofIdx]->partitionType(); }
375 assert(dofIdx <
numDof());
377 return elements_[dofIdx];
384 const Entity&
entity(
unsigned dofIdx)
const
392 {
return subControlVolumes_[dofIdx]; }
398 {
return interiorFaces_.size(); }
405 {
return interiorFaces_[faceIdx]; }
411 {
return boundaryFaces_.size(); }
418 {
return boundaryFaces_[bfIdx]; }
421 const GridView& gridView_;
422 const ElementMapper& elementMapper_;
424 std::vector<Element> elements_;
425 std::vector<SubControlVolume> subControlVolumes_;
426 std::vector<SubControlVolumeFace> interiorFaces_;
427 std::vector<BoundaryFace> boundaryFaces_;
const GlobalPosition & integrationPos() const
Returns the global position of the face's integration point.
Definition ecfvstencil.hh:202
unsigned short interiorIndex() const
Returns the local index of the degree of freedom to the face's interior.
Definition ecfvstencil.hh:181
unsigned short exteriorIndex() const
Returns the local index of the degree of freedom to the face's outside.
Definition ecfvstencil.hh:195
const WorldVector & normal() const
Returns the outer unit normal at the face's integration point.
Definition ecfvstencil.hh:209
FaceDir::DirEnum faceDirFromDirId() const
Returns the direction of the face.
Definition ecfvstencil.hh:230
Scalar area() const
Returns the area [m^2] of the face.
Definition ecfvstencil.hh:215
int dirId() const
Returns the direction id of the face w.r.t the cell.
Definition ecfvstencil.hh:224
const LocalGeometry geometry() const
The geometry of the sub-control volume.
Definition ecfvstencil.hh:135
Scalar volume() const
The volume [m^3] occupied by the sub-control volume.
Definition ecfvstencil.hh:129
const LocalGeometry localGeometry() const
Geometry of the sub-control volume relative to parent.
Definition ecfvstencil.hh:141
decltype(auto) center() const
The center of the sub-control volume.
Definition ecfvstencil.hh:123
decltype(auto) globalPos() const
The global position associated with the sub-control volume.
Definition ecfvstencil.hh:117
Represents the stencil (finite volume geometry) of a single element in the ECFV discretization.
Definition ecfvstencil.hh:69
const SubControlVolumeFace & interiorFace(unsigned faceIdx) const
Returns the face object belonging to a given face index in the interior of the domain.
Definition ecfvstencil.hh:404
const BoundaryFace & boundaryFace(unsigned bfIdx) const
Returns the boundary face object belonging to a given boundary face index.
Definition ecfvstencil.hh:417
const Element & element(unsigned dofIdx) const
Return the element given the index of a degree of freedom.
Definition ecfvstencil.hh:373
const Element & element() const
Return the element to which the stencil refers.
Definition ecfvstencil.hh:320
const Entity & entity(unsigned dofIdx) const
Return the entity given the index of a degree of freedom.
Definition ecfvstencil.hh:384
std::size_t numDof() const
Returns the number of degrees of freedom which the current element interacts with.
Definition ecfvstencil.hh:334
unsigned globalSpaceIndex(unsigned dofIdx) const
Return the global space index given the index of a degree of freedom.
Definition ecfvstencil.hh:353
std::size_t numPrimaryDof() const
Returns the number of degrees of freedom which are contained by within the current element.
Definition ecfvstencil.hh:346
const GridView & gridView() const
Return the grid view of the element to which the stencil refers.
Definition ecfvstencil.hh:327
Dune::PartitionType partitionType(unsigned dofIdx) const
Return partition type of a given degree of freedom.
Definition ecfvstencil.hh:363
std::size_t numInteriorFaces() const
Returns the number of interior faces of the stencil.
Definition ecfvstencil.hh:397
std::size_t numBoundaryFaces() const
Returns the number of boundary faces of the stencil.
Definition ecfvstencil.hh:410
const SubControlVolume & subControlVolume(unsigned dofIdx) const
Returns the sub-control volume object belonging to a given degree of freedom.
Definition ecfvstencil.hh:391
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
Quadrature geometry for quadrilaterals.