24#ifndef DUMUX_DISCRETIZATION_STAGGERED_FREE_FLOW_SUBCONTROLVOLUMEFACE_HH
25#define DUMUX_DISCRETIZATION_STAGGERED_FREE_FLOW_SUBCONTROLVOLUMEFACE_HH
29#include <dune/common/fvector.hh>
30#include <dune/geometry/type.hh>
31#include <dune/geometry/multilineargeometry.hh>
48 template<
class Container>
auto operator()(Container&& c)
49 ->
decltype(c.resize(1)) {}
53constexpr bool hasResize()
54{
return decltype(
isValid(HasResize{})(std::declval<C>()))::value; }
66template<
class Gr
idView,
int upwindSchemeOrder>
71 using Scalar =
typename GridView::ctype;
76 using Grid =
typename GridView::Grid;
77 static constexpr int dim = Grid::dimension;
78 static constexpr int dimWorld = Grid::dimensionworld;
86 template<
int mydim,
int cdim >
89 using Type = std::array< Dune::FieldVector< ct, cdim >, (1<<(
dim-1)) >;
103template<
class SubControlVolumeFace>
105 const typename SubControlVolumeFace::GlobalPosition& newCenter)
108 bf.setCenter(newCenter);
109 bf.setBoundary(
true);
110 bf.setIsGhostFace(
true);
120 int upwindSchemeOrder,
121 class T = FreeFlowStaggeredDefaultScvfGeometryTraits<GV, upwindSchemeOrder>>
127 using Geometry =
typename T::Geometry;
130 using CornerStorage =
typename T::CornerStorage;
132 using PairData =
typename T::PairData;
133 using AxisData =
typename T::AxisData;
135 using Scalar =
typename T::Scalar;
136 static const int dim = GV::dimension;
138 static constexpr int numPairs = 2 * (dim - 1);
140 static constexpr bool useHigherOrder = upwindSchemeOrder > 1;
154 template <
class Intersection>
156 const typename Intersection::Geometry& isGeometry,
157 GridIndexType scvfIndex,
158 const std::vector<GridIndexType>& scvIndices,
159 const typename T::GeometryHelper& geometryHelper)
161 geomType_(isGeometry.type()),
162 area_(isGeometry.
volume()),
163 center_(isGeometry.
center()),
164 unitOuterNormal_(is.centerUnitOuterNormal()),
165 scvfIndex_(scvfIndex),
166 scvIndices_(scvIndices),
168 axisData_(geometryHelper.
axisData()),
169 pairData_(std::move(geometryHelper.
pairData())),
170 localFaceIdx_(geometryHelper.localFaceIndex()),
175 if constexpr (Detail::hasResize<CornerStorage>())
176 corners_.resize(isGeometry.corners());
178 for (
int i = 0; i < isGeometry.corners(); ++i)
179 corners_[i] = isGeometry.corner(i);
216 return unitOuterNormal_;
222 return scvIndices_[0];
228 return scvIndices_[1];
240 assert(localIdx < corners_.size() &&
"provided index exceeds the number of corners");
241 return corners_[localIdx];
247 return Geometry(geomType_, corners_);
253 return localFaceIdx_;
271 return outerNormalSign_;
277 return pairData_[idx];
281 const std::array<PairData, numPairs>&
pairData()
const
303 if (localSubFaceIdx < 2)
320 return pairData(localSubFaceIdx).hasParallelNeighbor[parallelDegreeIdx];
343 return pairData(localSubFaceIdx).hasHalfParallelNeighbor;
367 return pairData(localSubFaceIdx).hasCornerParallelNeighbor;
377 return pairData(localSubFaceIdx).hasOuterLateral;
385 template<
bool enable = useHigherOrder, std::enable_if_t<enable,
int> = 0>
388 return axisData().hasBackwardNeighbor[backwardIdx];
396 template<
bool enable = useHigherOrder, std::enable_if_t<enable,
int> = 0>
399 return axisData().hasForwardNeighbor[forwardIdx];
417 return axisData().inAxisForwardDofs[0];
423 return axisData().inAxisBackwardDofs[0];
429 return axisData().selfToOppositeDistance;
440 if (parallelDegreeIdx == 0)
441 return (
faceLength(localSubFaceIdx) +
pairData(localSubFaceIdx).parallelCellWidths[0]) * 0.5;
445 assert((parallelDegreeIdx == 1) &&
"Only the width of the first two parallel cells (indices 0 and 1) is stored for each scvf.");
446 return (
pairData(localSubFaceIdx).parallelCellWidths[0] +
pairData(localSubFaceIdx).parallelCellWidths[1]) * 0.5;
456 { boundary_ = boundaryFlag; }
460 { isGhostFace_ = isGhostFaceFlag; }
463 Dune::GeometryType geomType_;
464 CornerStorage corners_;
468 GridIndexType scvfIndex_;
469 std::vector<GridIndexType> scvIndices_;
472 Scalar selfToOppositeDistance_;
474 std::array<PairData, numPairs> pairData_;
477 unsigned int dirIdx_;
478 int outerNormalSign_;
A helper function for class member function introspection.
Defines the index types used for grid and local indices.
Base class for a sub control volume face.
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:171
SubControlVolumeFace makeStaggeredBoundaryFace(const SubControlVolumeFace &scvf, const typename SubControlVolumeFace::GlobalPosition &newCenter)
Helper function to turn a given cell scvface into a fake boundary face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:104
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:641
constexpr auto isValid(const Expression &t)
A function that creates a test functor to do class member introspection at compile time.
Definition: isvalid.hh:93
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
Structure to define the index types used for grid and local indices.
Definition: indextraits.hh:38
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:39
unsigned int LocalIndex
Definition: indextraits.hh:40
Helper class constructing the dual grid finite volume geometries for the free flow staggered discreti...
Definition: staggeredgeometryhelper.hh:146
Detail::AxisData< GridView, upwindSchemeOrder > AxisData
Definition: staggeredgeometryhelper.hh:166
Detail::PairData< GridView, upwindSchemeOrder > PairData
Definition: staggeredgeometryhelper.hh:165
Default traits class to be used for the sub-control volume faces for the free-flow staggered finite v...
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:68
static constexpr int dimWorld
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:78
typename ScvfMLGTraits< Scalar >::template CornerStorage< dim-1, dimWorld >::Type CornerStorage
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:94
typename IndexTraits< GridView >::GridIndex GridIndexType
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:69
static constexpr int dim
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:77
typename GridView::Grid Grid
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:76
typename GeometryHelper::AxisData AxisData
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:74
typename CornerStorage::value_type GlobalPosition
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:95
typename GeometryHelper::PairData PairData
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:73
Dune::MultiLinearGeometry< Scalar, dim-1, dimWorld, ScvfMLGTraits< Scalar > > Geometry
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:93
typename GridView::ctype Scalar
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:71
typename IndexTraits< GridView >::LocalIndex LocalIndexType
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:70
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:83
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:88
std::array< Dune::FieldVector< ct, cdim >,(1<<(dim-1)) > Type
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:89
Class for a sub control volume face in the staggered method, i.e a part of the boundary of a sub cont...
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:124
LocalIndexType localFaceIdx() const
The local index of this sub control volume face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:251
const GlobalPosition & center() const
The center of the sub control volume face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:183
GridIndexType dofIndex() const
Returns the dof of the face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:403
void setBoundary(bool boundaryFlag)
set the boundary flag
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:455
const AxisData & axisData() const
Return an array of all pair data.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:287
bool hasOuterLateral(const int localSubFaceIdx) const
Check if the face has an outer normal neighbor.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:375
GridIndexType dofIndexForwardFace() const
Returns the dof the first forward face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:415
bool normalInPosCoordDir() const
Returns whether the unitNormal of the face points in positive coordinate direction.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:263
FreeFlowStaggeredSubControlVolumeFace()=default
const std::array< PairData, numPairs > & pairData() const
Return an array of all pair data.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:281
GridIndexType dofIndexOpposingFace() const
Returns the dof of the opposing face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:409
const Geometry geometry() const
The geometry of the sub control volume face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:245
Scalar parallelDofsDistance(const int localSubFaceIdx, const int parallelDegreeIdx) const
Returns the distance between the parallel dofs.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:438
Scalar selfToOppositeDistance() const
Returns the distance between the face and the opposite one.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:427
const GlobalPosition & ipGlobal() const
The integration point for flux evaluations in global coordinates.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:195
unsigned int directionIndex() const
Returns the direction index of the facet (0 = x, 1 = y, 2 = z)
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:257
bool isGhostFace() const
Returns true if the face is a ghost face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:293
const GlobalPosition & unitOuterNormal() const
The unit outer normal vector.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:214
bool hasHalfParallelNeighbor(const int localSubFaceIdx) const
Check if the face has a half parallel neighbor.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:341
void setIsGhostFace(bool isGhostFaceFlag)
set the ghost face flag
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:459
int directionSign() const
Returns the sign of the unit outer normal's vector.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:269
bool hasForwardNeighbor(const int forwardIdx) const
Check if the face has a forward neighbor.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:397
void setCenter(const GlobalPosition ¢er)
set the center to a different position
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:451
GridIndexType outsideScvIdx() const
index of the outside sub control volume for spatial param evaluation
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:226
static constexpr int numCornersPerFace
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:145
FreeFlowStaggeredSubControlVolumeFace(const Intersection &is, const typename Intersection::Geometry &isGeometry, GridIndexType scvfIndex, const std::vector< GridIndexType > &scvIndices, const typename T::GeometryHelper &geometryHelper)
Constructor with intersection.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:155
typename T::GlobalPosition GlobalPosition
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:143
const GlobalPosition & corner(unsigned int localIdx) const
The positions of the corners.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:238
bool hasBackwardNeighbor(const int backwardIdx) const
Check if the face has a backward neighbor.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:386
T Traits
State the traits public and thus export all types.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:148
bool boundary() const
Returns boolean if the sub control volume face is on the boundary.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:208
const PairData & pairData(const int idx) const
Returns the data for one sub face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:275
const GlobalPosition & dofPosition() const
The position of the dof living on the face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:189
GridIndexType index() const
The global index of this sub control volume face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:232
GridIndexType dofIndexBackwardFace() const
Returns the dof of the first backward face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:421
Scalar faceLength(const int localSubFaceIdx) const
Returns the length of the face in a certain direction (adaptation of area() for 3d)
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:299
Scalar area() const
The area of the sub control volume face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:202
bool hasParallelNeighbor(const int localSubFaceIdx, const int parallelDegreeIdx) const
Check if the face has a parallel neighbor.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:318
bool hasCornerParallelNeighbor(const int localSubFaceIdx) const
Check if the face has a corner parallel neighbor.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:365
GridIndexType insideScvIdx() const
Index of the inside sub control volume for spatial param evaluation.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:220
Base class for a sub control volume face, i.e a part of the boundary of a sub control volume we compu...
Definition: subcontrolvolumefacebase.hh:41