3.2-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
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{
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 GridIndexType outsideScvIdx() const
207 {
208 return scvIndices_[1];
209 }
210
212 GridIndexType index() const
213 {
214 return scvfIndex_;
215 }
216
218 const GlobalPosition& corner(unsigned int localIdx) const
219 {
220 assert(localIdx < corners_.size() && "provided index exceeds the number of corners");
221 return corners_[localIdx];
222 }
223
225 const Geometry geometry() const
226 {
227 return Geometry(geomType_, corners_);
228 }
229
231 LocalIndexType localFaceIdx() const
232 {
233 return localFaceIdx_;
234 }
235
237 unsigned int directionIndex() const
238 {
239 return dirIdx_;
240 }
241
244 {
245 return directionSign() > 0;
246 }
247
249 int directionSign() const
250 {
251 return outerNormalSign_;
252 }
253
255 const PairData& pairData(const int idx) const
256 {
257 return pairData_[idx];
258 }
259
261 const std::array<PairData, numPairs>& pairData() const
262 {
263 return pairData_;
264 }
265
267 const AxisData& axisData() const
268 {
269 return axisData_;
270 }
271
273 bool isGhostFace() const
274 {
275 return isGhostFace_;
276 }
277
279 Scalar faceLength(const int localSubFaceIdx) const
280 {
281 if (dim == 3)
282 {
283 if (localSubFaceIdx < 2)
284 return (corner(1) - corner(0)).two_norm();
285 else
286 return (corner(2) - corner(0)).two_norm();
287 }
288 else
289 return (corner(1) - corner(0)).two_norm();
290 }
291
298 bool hasParallelNeighbor(const int localSubFaceIdx, const int parallelDegreeIdx) const
299 {
300 return pairData(localSubFaceIdx).hasParallelNeighbor[parallelDegreeIdx];
301 }
302
308 bool hasOuterLateral(const int localSubFaceIdx) const
309 {
310 return pairData(localSubFaceIdx).hasOuterLateral;
311 }
312
318 template<bool enable = useHigherOrder, std::enable_if_t<enable, int> = 0>
319 bool hasBackwardNeighbor(const int backwardIdx) const
320 {
321 return axisData().hasBackwardNeighbor[backwardIdx];
322 }
323
329 template<bool enable = useHigherOrder, std::enable_if_t<enable, int> = 0>
330 bool hasForwardNeighbor(const int forwardIdx) const
331 {
332 return axisData().hasForwardNeighbor[forwardIdx];
333 }
334
336 GridIndexType dofIndex() const
337 {
338 return axisData().selfDof;
339 }
340
342 GridIndexType dofIndexOpposingFace() const
343 {
344 return axisData().oppositeDof;
345 }
346
348 GridIndexType dofIndexForwardFace() const
349 {
350 return axisData().inAxisForwardDofs[0];
351 }
352
354 GridIndexType dofIndexBackwardFace() const
355 {
356 return axisData().inAxisBackwardDofs[0];
357 }
358
361 {
362 return axisData().selfToOppositeDistance;
363 }
364
371 Scalar parallelDofsDistance(const int localSubFaceIdx, const int parallelDegreeIdx) const
372 {
373 if (parallelDegreeIdx == 0)
374 return (faceLength(localSubFaceIdx) + pairData(localSubFaceIdx).parallelCellWidths[0]) * 0.5;
375 // pairData(localSubFaceIdx).parallelCellWidths[0]) will return 0.0 if the subface perpendicular the scvf lies on a boundary
376 else
377 {
378 assert((parallelDegreeIdx == 1) && "Only the width of the first two parallel cells (indicies 0 and 1) is stored for each scvf.");
379 return (pairData(localSubFaceIdx).parallelCellWidths[0] + pairData(localSubFaceIdx).parallelCellWidths[1]) * 0.5;
380 }
381 }
382
391 {
392 FreeFlowStaggeredSubControlVolumeFace boundaryFace = *this;
393 boundaryFace.center_ = pos;
394 boundaryFace.boundary_ = true;
395 boundaryFace.isGhostFace_ = true;
396 return boundaryFace;
397 }
398
399private:
400 void maybeResizeCornerStorage_(std::true_type /*hasResize*/, std::size_t size)
401 { corners_.resize(size); }
402
403 void maybeResizeCornerStorage_(std::false_type /*hasResize*/, std::size_t size)
404 {}
405
406 Dune::GeometryType geomType_;
407 CornerStorage corners_;
408 Scalar area_;
409 GlobalPosition center_;
410 GlobalPosition unitOuterNormal_;
411 GridIndexType scvfIndex_;
412 std::vector<GridIndexType> scvIndices_;
413 bool boundary_;
414
415 Scalar selfToOppositeDistance_;
416 AxisData axisData_;
417 std::array<PairData, numPairs> pairData_;
418
419 int localFaceIdx_;
420 unsigned int dirIdx_;
421 int outerNormalSign_;
422 bool isGhostFace_;
423};
424
425} // end namespace Dumux
426
427#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.
Base class for a sub control volume face.
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:618
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
Struture 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
Parallel Data stored per sub face.
Definition: staggeredgeometryhelper.hh:48
In Axis Data stored per sub face.
Definition: staggeredgeometryhelper.hh:77
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
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:106
LocalIndexType localFaceIdx() const
The local index of this sub control volume face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:231
const GlobalPosition & center() const
The center of the sub control volume face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:163
GridIndexType dofIndex() const
Returns the dof of the face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:336
const AxisData & axisData() const
Return an array of all pair data.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:267
bool hasOuterLateral(const int localSubFaceIdx) const
Check if the face has an outer normal neighbor.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:308
GridIndexType dofIndexForwardFace() const
Returns the dof the first forward face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:348
bool normalInPosCoordDir() const
Returns whether the unitNormal of the face points in positive coordinate direction.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:243
const std::array< PairData, numPairs > & pairData() const
Return an array of all pair data.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:261
GridIndexType dofIndexOpposingFace() const
Returns the dof of the opposing face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:342
const Geometry geometry() const
The geometry of the sub control volume face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:225
Scalar parallelDofsDistance(const int localSubFaceIdx, const int parallelDegreeIdx) const
Returns the distance between the parallel dofs.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:371
Scalar selfToOppositeDistance() const
Returns the distance between the face and the opposite one.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:360
const GlobalPosition & ipGlobal() const
The integration point for flux evaluations in global coordinates.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:175
unsigned int directionIndex() const
Returns the dirction index of the facet (0 = x, 1 = y, 2 = z)
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:237
bool isGhostFace() const
Returns true if the face is a ghost face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:273
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:249
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:390
bool hasForwardNeighbor(const int forwardIdx) const
Check if the face has a forward neighbor.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:330
GridIndexType outsideScvIdx() const
index of the outside sub control volume for spatial param evaluation
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:206
typename T::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:218
bool hasBackwardNeighbor(const int backwardIdx) const
Check if the face has a backward neighbor.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:319
T Traits
State the traits public and thus export all types.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:128
bool boundary() const
Returns bolean if the sub control volume face is on the boundary.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:188
const PairData & pairData(const int idx) const
Returns the data for one sub face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:255
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:212
GridIndexType dofIndexBackwardFace() const
Returns the dof of the first backward face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:354
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:279
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:298
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