3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
discretization/staggered/freeflow/subcontrolvolumeface.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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program 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 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program 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 this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_DISCRETIZATION_STAGGERED_FREE_FLOW_SUBCONTROLVOLUMEFACE_HH
25#define DUMUX_DISCRETIZATION_STAGGERED_FREE_FLOW_SUBCONTROLVOLUMEFACE_HH
26
27#include <array>
28#include <utility>
29#include <dune/common/fvector.hh>
30#include <dune/geometry/type.hh>
31#include <dune/geometry/multilineargeometry.hh>
32
38
39#include <typeinfo>
40
41namespace Dumux {
42
43#ifndef DOXYGEN
44namespace Detail {
45// helper struct detecting if the container class storing the scvf's corners has a resize function
46struct HasResize
47{
48 template<class Container> auto operator()(Container&& c)
49 -> decltype(c.resize(1)) {}
50};
51
52template<class C>
53constexpr bool hasResize()
54{ return decltype(isValid(HasResize{})(std::declval<C>()))::value; }
55
56} // end namespace Detail
57#endif
58
66template<class GridView, int upwindSchemeOrder>
68{
71 using Scalar = typename GridView::ctype;
75
76 using Grid = typename GridView::Grid;
77 static constexpr int dim = Grid::dimension;
78 static constexpr int dimWorld = Grid::dimensionworld;
79
80 // we use geometry traits that use static corner vectors to and a fixed geometry type
81 template <class ct>
82 struct ScvfMLGTraits : public Dune::MultiLinearGeometryTraits<ct>
83 {
84 // we use static vectors to store the corners as we know
85 // the number of corners in advance (2^(dim-1) corners (1<<(dim-1))
86 template< int mydim, int cdim >
88 {
89 using Type = std::array< Dune::FieldVector< ct, cdim >, (1<<(dim-1)) >;
90 };
91 };
92
93 using Geometry = Dune::MultiLinearGeometry<Scalar, dim-1, dimWorld, ScvfMLGTraits<Scalar> >;
95 using GlobalPosition = typename CornerStorage::value_type;
96};
97
103template<class SubControlVolumeFace>
104SubControlVolumeFace makeStaggeredBoundaryFace(const SubControlVolumeFace& scvf,
105 const typename SubControlVolumeFace::GlobalPosition& newCenter)
106{
107 auto bf = scvf;
108 bf.setCenter(newCenter);
109 bf.setBoundary(true);
110 bf.setIsGhostFace(true);
111 return bf;
112}
113
119template<class GV,
120 int upwindSchemeOrder,
121 class T = FreeFlowStaggeredDefaultScvfGeometryTraits<GV, upwindSchemeOrder>>
123: public SubControlVolumeFaceBase<FreeFlowStaggeredSubControlVolumeFace<GV, upwindSchemeOrder, T>, T>
124{
127 using Geometry = typename T::Geometry;
128 using GridIndexType = typename IndexTraits<GV>::GridIndex;
129 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
130 using CornerStorage = typename T::CornerStorage;
131
132 using PairData = typename T::PairData;
133 using AxisData = typename T::AxisData;
134
135 using Scalar = typename T::Scalar;
136 static const int dim = GV::dimension;
137
138 static constexpr int numPairs = 2 * (dim - 1);
139
140 static constexpr bool useHigherOrder = upwindSchemeOrder > 1;
141
142public:
143 using GlobalPosition = typename T::GlobalPosition;
144
145 static constexpr int numCornersPerFace = 2 * (dim - 1);
146
148 using Traits = T;
149
150 // The default constructor
152
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)
160 : ParentType(),
161 geomType_(isGeometry.type()),
162 area_(isGeometry.volume()),
163 center_(isGeometry.center()),
164 unitOuterNormal_(is.centerUnitOuterNormal()),
165 scvfIndex_(scvfIndex),
166 scvIndices_(scvIndices),
167 boundary_(is.boundary()),
168 axisData_(geometryHelper.axisData()),
169 pairData_(std::move(geometryHelper.pairData())),
170 localFaceIdx_(geometryHelper.localFaceIndex()),
171 dirIdx_(geometryHelper.directionIndex()),
172 outerNormalSign_(sign(unitOuterNormal_[directionIndex()])),
173 isGhostFace_(false)
174 {
175 if constexpr (Detail::hasResize<CornerStorage>())
176 corners_.resize(isGeometry.corners());
177
178 for (int i = 0; i < isGeometry.corners(); ++i)
179 corners_[i] = isGeometry.corner(i);
180 }
181
183 const GlobalPosition& center() const
184 {
185 return center_;
186 }
187
190 {
191 return center_;
192 }
193
196 {
197 // Return center for now
198 return center_;
199 }
200
202 Scalar area() const
203 {
204 return area_;
205 }
206
208 bool boundary() const
209 {
210 return boundary_;
211 }
212
215 {
216 return unitOuterNormal_;
217 }
218
220 GridIndexType insideScvIdx() const
221 {
222 return scvIndices_[0];
223 }
224
226 GridIndexType outsideScvIdx() const
227 {
228 return scvIndices_[1];
229 }
230
232 GridIndexType index() const
233 {
234 return scvfIndex_;
235 }
236
238 const GlobalPosition& corner(unsigned int localIdx) const
239 {
240 assert(localIdx < corners_.size() && "provided index exceeds the number of corners");
241 return corners_[localIdx];
242 }
243
245 const Geometry geometry() const
246 {
247 return Geometry(geomType_, corners_);
248 }
249
251 LocalIndexType localFaceIdx() const
252 {
253 return localFaceIdx_;
254 }
255
257 unsigned int directionIndex() const
258 {
259 return dirIdx_;
260 }
261
264 {
265 return directionSign() > 0;
266 }
267
269 int directionSign() const
270 {
271 return outerNormalSign_;
272 }
273
275 const PairData& pairData(const int idx) const
276 {
277 return pairData_[idx];
278 }
279
281 const std::array<PairData, numPairs>& pairData() const
282 {
283 return pairData_;
284 }
285
287 const AxisData& axisData() const
288 {
289 return axisData_;
290 }
291
293 bool isGhostFace() const
294 {
295 return isGhostFace_;
296 }
297
299 Scalar faceLength(const int localSubFaceIdx) const
300 {
301 if (dim == 3)
302 {
303 if (localSubFaceIdx < 2)
304 return (corner(1) - corner(0)).two_norm();
305 else
306 return (corner(2) - corner(0)).two_norm();
307 }
308 else
309 return (corner(1) - corner(0)).two_norm();
310 }
311
318 bool hasParallelNeighbor(const int localSubFaceIdx, const int parallelDegreeIdx) const
319 {
320 return pairData(localSubFaceIdx).hasParallelNeighbor[parallelDegreeIdx];
321 }
322
341 bool hasHalfParallelNeighbor(const int localSubFaceIdx) const
342 {
343 return pairData(localSubFaceIdx).hasHalfParallelNeighbor;
344 }
345
365 bool hasCornerParallelNeighbor(const int localSubFaceIdx) const
366 {
367 return pairData(localSubFaceIdx).hasCornerParallelNeighbor;
368 }
369
375 bool hasOuterLateral(const int localSubFaceIdx) const
376 {
377 return pairData(localSubFaceIdx).hasOuterLateral;
378 }
379
385 template<bool enable = useHigherOrder, std::enable_if_t<enable, int> = 0>
386 bool hasBackwardNeighbor(const int backwardIdx) const
387 {
388 return axisData().hasBackwardNeighbor[backwardIdx];
389 }
390
396 template<bool enable = useHigherOrder, std::enable_if_t<enable, int> = 0>
397 bool hasForwardNeighbor(const int forwardIdx) const
398 {
399 return axisData().hasForwardNeighbor[forwardIdx];
400 }
401
403 GridIndexType dofIndex() const
404 {
405 return axisData().selfDof;
406 }
407
409 GridIndexType dofIndexOpposingFace() const
410 {
411 return axisData().oppositeDof;
412 }
413
415 GridIndexType dofIndexForwardFace() const
416 {
417 return axisData().inAxisForwardDofs[0];
418 }
419
421 GridIndexType dofIndexBackwardFace() const
422 {
423 return axisData().inAxisBackwardDofs[0];
424 }
425
428 {
429 return axisData().selfToOppositeDistance;
430 }
431
438 Scalar parallelDofsDistance(const int localSubFaceIdx, const int parallelDegreeIdx) const
439 {
440 if (parallelDegreeIdx == 0)
441 return (faceLength(localSubFaceIdx) + pairData(localSubFaceIdx).parallelCellWidths[0]) * 0.5;
442 // pairData(localSubFaceIdx).parallelCellWidths[0]) will return 0.0 if the subface perpendicular the scvf lies on a boundary
443 else
444 {
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;
447 }
448 }
449
452 { center_ = center; }
453
455 void setBoundary(bool boundaryFlag)
456 { boundary_ = boundaryFlag; }
457
459 void setIsGhostFace(bool isGhostFaceFlag)
460 { isGhostFace_ = isGhostFaceFlag; }
461
462private:
463 Dune::GeometryType geomType_;
464 CornerStorage corners_;
465 Scalar area_;
466 GlobalPosition center_;
467 GlobalPosition unitOuterNormal_;
468 GridIndexType scvfIndex_;
469 std::vector<GridIndexType> scvIndices_;
470 bool boundary_;
471
472 Scalar selfToOppositeDistance_;
473 AxisData axisData_;
474 std::array<PairData, numPairs> pairData_;
475
476 int localFaceIdx_;
477 unsigned int dirIdx_;
478 int outerNormalSign_;
479 bool isGhostFace_;
480};
481
482} // end namespace Dumux
483
484#endif
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
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 &center)
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