3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
Loading...
Searching...
No Matches
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
46// for g++ > 5.3, this can be replaced by a lambda
47struct hasResize
48{
49 template<class Container>
50 auto operator()(Container&& c)
51 -> decltype(c.resize(1))
52 {}
53};
54} // end namespace Detail
55#endif
56
64template<class GridView, int upwindSchemeOrder>
66{
69 using Scalar = typename GridView::ctype;
72
73
74 using Grid = typename GridView::Grid;
75 static constexpr int dim = Grid::dimension;
76 static constexpr int dimWorld = Grid::dimensionworld;
77
78 // we use geometry traits that use static corner vectors to and a fixed geometry type
79 template <class ct>
80 struct ScvfMLGTraits : public Dune::MultiLinearGeometryTraits<ct>
81 {
82 // we use static vectors to store the corners as we know
83 // the number of corners in advance (2^(dim-1) corners (1<<(dim-1))
84 template< int mydim, int cdim >
86 {
87 using Type = std::array< Dune::FieldVector< ct, cdim >, (1<<(dim-1)) >;
88 };
89 };
90
91 using Geometry = Dune::MultiLinearGeometry<Scalar, dim-1, dimWorld, ScvfMLGTraits<Scalar> >;
93 using GlobalPosition = typename CornerStorage::value_type;
94};
95
101template<class GV,
102 int upwindSchemeOrder,
105: public SubControlVolumeFaceBase<FreeFlowStaggeredSubControlVolumeFace<GV, upwindSchemeOrder, T>, T>
106{
108 using ParentType = SubControlVolumeFaceBase<ThisType, T>;
109 using Geometry = typename T::Geometry;
110 using GridIndexType = typename IndexTraits<GV>::GridIndex;
111 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
112 using CornerStorage = typename T::CornerStorage;
113
114 using PairData = typename T::PairData;
115 using AxisData = typename T::AxisData;
116
117 using Scalar = typename T::Scalar;
118 static const int dim = GV::dimension;
119
120 static constexpr int numPairs = 2 * (dim - 1);
121
122 static constexpr bool useHigherOrder = upwindSchemeOrder > 1;
123
124public:
125 using GlobalPosition = typename T::GlobalPosition;
126
128 using Traits = T;
129
130 // The default constructor
132
134 template <class Intersection, class GeometryHelper>
136 const typename Intersection::Geometry& isGeometry,
137 GridIndexType scvfIndex,
138 const std::vector<GridIndexType>& scvIndices,
139 const GeometryHelper& geometryHelper)
140 : ParentType(),
141 geomType_(isGeometry.type()),
142 area_(isGeometry.volume()),
143 center_(isGeometry.center()),
144 unitOuterNormal_(is.centerUnitOuterNormal()),
145 scvfIndex_(scvfIndex),
146 scvIndices_(scvIndices),
147 boundary_(is.boundary()),
148
149 axisData_(geometryHelper.axisData()),
150 pairData_(std::move(geometryHelper.pairData())),
151 localFaceIdx_(geometryHelper.localFaceIndex()),
152 dirIdx_(geometryHelper.directionIndex()),
153 outerNormalSign_(sign(unitOuterNormal_[directionIndex()])),
154 isGhostFace_(false)
155 {
156 using HasResize = decltype(isValid(Detail::hasResize())(corners_));
157 maybeResizeCornerStorage_(HasResize{}, isGeometry.corners());
158 for (int i = 0; i < isGeometry.corners(); ++i)
159 corners_[i] = isGeometry.corner(i);
160 }
161
163 const GlobalPosition& center() const
164 {
165 return center_;
166 }
167
170 {
171 return center_;
172 }
173
176 {
177 // Return center for now
178 return center_;
179 }
180
182 Scalar area() const
183 {
184 return area_;
185 }
186
188 bool boundary() const
189 {
190 return boundary_;
191 }
192
195 {
196 return unitOuterNormal_;
197 }
198
200 GridIndexType insideScvIdx() const
201 {
202 return scvIndices_[0];
203 }
204
206 // This results in undefined behaviour if boundary is true
207 GridIndexType outsideScvIdx() const
208 {
209 return scvIndices_[1];
210 }
211
213 GridIndexType index() const
214 {
215 return scvfIndex_;
216 }
217
219 const GlobalPosition& corner(unsigned int localIdx) const
220 {
221 assert(localIdx < corners_.size() && "provided index exceeds the number of corners");
222 return corners_[localIdx];
223 }
224
226 const Geometry geometry() const
227 {
228 return Geometry(geomType_, corners_);
229 }
230
232 LocalIndexType localFaceIdx() const
233 {
234 return localFaceIdx_;
235 }
236
238 unsigned int directionIndex() const
239 {
240 return dirIdx_;
241 }
242
245 {
246 return directionSign() > 0;
247 }
248
250 int directionSign() const
251 {
252 return outerNormalSign_;
253 }
254
256 const PairData& pairData(const int idx) const
257 {
258 return pairData_[idx];
259 }
260
262 const std::array<PairData, numPairs>& pairData() const
263 {
264 return pairData_;
265 }
266
268 const AxisData& axisData() const
269 {
270 return axisData_;
271 }
272
274 bool isGhostFace() const
275 {
276 return isGhostFace_;
277 }
278
280 Scalar faceLength(const int localSubFaceIdx) const
281 {
282 if (dim == 3)
283 {
284 if (localSubFaceIdx < 2)
285 return (corner(1) - corner(0)).two_norm();
286 else
287 return (corner(2) - corner(0)).two_norm();
288 }
289 else
290 return (corner(1) - corner(0)).two_norm();
291 }
292
299 bool hasParallelNeighbor(const int localSubFaceIdx, const int parallelDegreeIdx) const
300 {
301 return pairData(localSubFaceIdx).hasParallelNeighbor[parallelDegreeIdx];
302 }
303
309 bool hasOuterLateral(const int localSubFaceIdx) const
310 {
311 return pairData(localSubFaceIdx).hasOuterLateral;
312 }
313
319 template<bool enable = useHigherOrder, std::enable_if_t<enable, int> = 0>
320 bool hasBackwardNeighbor(const int backwardIdx) const
321 {
322 return axisData().hasBackwardNeighbor[backwardIdx];
323 }
324
330 template<bool enable = useHigherOrder, std::enable_if_t<enable, int> = 0>
331 bool hasForwardNeighbor(const int forwardIdx) const
332 {
333 return axisData().hasForwardNeighbor[forwardIdx];
334 }
335
337 GridIndexType dofIndex() const
338 {
339 return axisData().selfDof;
340 }
341
343 GridIndexType dofIndexOpposingFace() const
344 {
345 return axisData().oppositeDof;
346 }
347
349 GridIndexType dofIndexForwardFace() const
350 {
351 return axisData().inAxisForwardDofs[0];
352 }
353
355 GridIndexType dofIndexBackwardFace() const
356 {
357 return axisData().inAxisBackwardDofs[0];
358 }
359
362 {
363 return axisData().selfToOppositeDistance;
364 }
365
372 Scalar parallelDofsDistance(const int localSubFaceIdx, const int parallelDegreeIdx) const
373 {
374 if (parallelDegreeIdx == 0)
375 return (faceLength(localSubFaceIdx) + pairData(localSubFaceIdx).parallelCellWidths[0]) * 0.5;
376 // pairData(localSubFaceIdx).parallelCellWidths[0]) will return 0.0 if the subface perpendicular the scvf lies on a boundary
377 else
378 {
379 assert((parallelDegreeIdx == 1) && "Only the width of the first two parallel cells (indicies 0 and 1) is stored for each scvf.");
380 return (pairData(localSubFaceIdx).parallelCellWidths[0] + pairData(localSubFaceIdx).parallelCellWidths[1]) * 0.5;
381 }
382 }
383
392 {
393 FreeFlowStaggeredSubControlVolumeFace boundaryFace = *this;
394 boundaryFace.center_ = pos;
395 boundaryFace.boundary_ = true;
396 boundaryFace.isGhostFace_ = true;
397 return boundaryFace;
398 }
399
400private:
401 void maybeResizeCornerStorage_(std::true_type /*hasResize*/, std::size_t size)
402 { corners_.resize(size); }
403
404 void maybeResizeCornerStorage_(std::false_type /*hasResize*/, std::size_t size)
405 {}
406
407 Dune::GeometryType geomType_;
408 CornerStorage corners_;
409 Scalar area_;
410 GlobalPosition center_;
411 GlobalPosition unitOuterNormal_;
412 GridIndexType scvfIndex_;
413 std::vector<GridIndexType> scvIndices_;
414 bool boundary_;
415
416 Scalar selfToOppositeDistance_;
417 AxisData axisData_;
418 std::array<PairData, numPairs> pairData_;
419
420 int localFaceIdx_;
421 unsigned int dirIdx_;
422 int outerNormalSign_;
423 bool isGhostFace_;
424};
425
426} // end namespace Dumux
427
428#endif // DUMUX_DISCRETIZATION_STAGGERED_FREE_FLOW_SUBCONTROLVOLUMEFACE_HH
A helper function for class member function introspection.
Defines the index types used for grid and local indices.
Helper class constructing the dual grid finite volume geometries for the free flow staggered discreti...
Base class for a sub control volume face.
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition math.hh:618
make the local view function available whenever we use the grid geometry
Definition adapt.hh:29
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
Definition state.hh:28
const auto hasResize
Definition test_isvalid.cc:22
typename GridView::IndexSet::IndexType GridIndex
Definition indextraits.hh:39
unsigned int LocalIndex
Definition indextraits.hh:40
Detail::AxisData< GridView, upwindSchemeOrder > AxisData
Definition staggeredgeometryhelper.hh:147
Detail::PairData< GridView, upwindSchemeOrder > PairData
Definition staggeredgeometryhelper.hh:146
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:66
static constexpr int dimWorld
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:76
typename ScvfMLGTraits< Scalar >::template CornerStorage< dim-1, dimWorld >::Type CornerStorage
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:92
typename IndexTraits< GridView >::GridIndex GridIndexType
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:67
static constexpr int dim
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:75
typename GridView::Grid Grid
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:74
typename CornerStorage::value_type GlobalPosition
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:93
Dune::MultiLinearGeometry< Scalar, dim-1, dimWorld, ScvfMLGTraits< Scalar > > Geometry
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:91
typename FreeFlowStaggeredGeometryHelper< GridView, upwindSchemeOrder >::AxisData AxisData
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:71
typename GridView::ctype Scalar
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:69
typename FreeFlowStaggeredGeometryHelper< GridView, upwindSchemeOrder >::PairData PairData
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:70
typename IndexTraits< GridView >::LocalIndex LocalIndexType
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:68
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:81
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:86
std::array< Dune::FieldVector< ct, cdim >,(1<<(dim-1)) > Type
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:87
LocalIndexType localFaceIdx() const
The local index of this sub control volume face.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:232
const GlobalPosition & center() const
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:163
GridIndexType dofIndex() const
Returns the dof of the face.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:337
const AxisData & axisData() const
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:268
bool hasOuterLateral(const int localSubFaceIdx) const
Check if the face has an outer normal neighbor.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:309
GridIndexType dofIndexForwardFace() const
Returns the dof the first forward face.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:349
bool normalInPosCoordDir() const
Returns whether the unitNormal of the face points in positive coordinate direction.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:244
const std::array< PairData, numPairs > & pairData() const
Return an array of all pair data.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:262
GridIndexType dofIndexOpposingFace() const
Returns the dof of the opposing face.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:343
const Geometry geometry() const
The geometry of the sub control volume face.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:226
Scalar parallelDofsDistance(const int localSubFaceIdx, const int parallelDegreeIdx) const
Returns the distance between the parallel dofs.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:372
Scalar selfToOppositeDistance() const
Returns the distance between the face and the opposite one.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:361
const GlobalPosition & ipGlobal() const
The integration point for flux evaluations in global coordinates.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:175
unsigned int directionIndex() const
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:238
bool isGhostFace() const
Returns true if the face is a ghost face.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:274
const GlobalPosition & unitOuterNormal() const
The unit outer normal vector.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:194
FreeFlowStaggeredSubControlVolumeFace(const Intersection &is, const typename Intersection::Geometry &isGeometry, GridIndexType scvfIndex, const std::vector< GridIndexType > &scvIndices, const GeometryHelper &geometryHelper)
Constructor with intersection.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:135
int directionSign() const
Returns the sign of the unit outer normal's vector.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:250
FreeFlowStaggeredSubControlVolumeFace makeBoundaryFace(const GlobalPosition &pos) const
Returns a copy of the own scvf whith a user-specified center position. This is needed for retrieving ...
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:391
bool hasForwardNeighbor(const int forwardIdx) const
Check if the face has a forward neighbor.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:331
GridIndexType outsideScvIdx() const
index of the outside sub control volume for spatial param evaluation
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:207
typename FreeFlowStaggeredDefaultScvfGeometryTraits< GridView, upwindSchemeOrder >::GlobalPosition GlobalPosition
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:125
const GlobalPosition & corner(unsigned int localIdx) const
The positions of the corners.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:219
bool hasBackwardNeighbor(const int backwardIdx) const
Check if the face has a backward neighbor.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:320
FreeFlowStaggeredDefaultScvfGeometryTraits< GridView, upwindSchemeOrder > Traits
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:128
bool boundary() const
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:188
const PairData & pairData(const int idx) const
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:256
const GlobalPosition & dofPosition() const
The position of the dof living on the face.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:169
GridIndexType index() const
The global index of this sub control volume face.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:213
GridIndexType dofIndexBackwardFace() const
Returns the dof of the first backward face.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:355
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:280
Scalar area() const
The area of the sub control volume face.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:182
bool hasParallelNeighbor(const int localSubFaceIdx, const int parallelDegreeIdx) const
Check if the face has a parallel neighbor.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:299
GridIndexType insideScvIdx() const
Index of the inside sub control volume for spatial param evaluation.
Definition discretization/staggered/freeflow/subcontrolvolumeface.hh:200
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
Class for a sub control volume face in the staggered method, i.e a part of the boundary of a sub cont...