3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
velocity.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 *****************************************************************************/
25#ifndef DUMUX_POROUSMEDIUMFLOW_VELOCITY_HH
26#define DUMUX_POROUSMEDIUMFLOW_VELOCITY_HH
27
28#include <vector>
29#include <type_traits>
30
31#include <dune/common/fvector.hh>
32#include <dune/common/float_cmp.hh>
33#include <dune/geometry/type.hh>
34
38#include <dumux/flux/traits.hh>
39
40namespace Dumux {
41
42#ifndef DOXYGEN
43namespace Detail {
44// helper structs and functions detecting if the model is compositional
45
46template <typename T, typename ...Ts>
47using MoleFractionDetector = decltype(std::declval<T>().moleFraction(std::declval<Ts>()...));
48
49template<class T, typename ...Args>
50static constexpr bool hasMoleFraction()
51{ return Dune::Std::is_detected<MoleFractionDetector, T, Args...>::value; }
52
53template <typename T, typename ...Ts>
54using MassFractionDetector = decltype(std::declval<T>().massFraction(std::declval<Ts>()...));
55
56template<class T, typename ...Args>
57static constexpr bool hasMassFraction()
58{ return Dune::Std::is_detected<MassFractionDetector, T, Args...>::value; }
59
60} // end namespace Detail
61#endif // DOXYGEN
62
67template<class GridVariables, class FluxVariables>
69{
70 using GridGeometry = typename GridVariables::GridGeometry;
71 using FVElementGeometry = typename GridGeometry::LocalView;
72 using SubControlVolume = typename GridGeometry::SubControlVolume;
73 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
74 using GridView = typename GridGeometry::GridView;
75 using Element = typename GridView::template Codim<0>::Entity;
76 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
77 using ElementFluxVarsCache = typename GridVariables::GridFluxVariablesCache::LocalView;
78 using VolumeVariables = typename GridVariables::VolumeVariables;
79 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
80 using FluidSystem = typename VolumeVariables::FluidSystem;
81 using Scalar = typename GridVariables::Scalar;
82 using FluxTraits = typename Dumux::FluxTraits<FluxVariables>;
83 using AdvectionType = typename FluxVariables::AdvectionType;
84
85 static constexpr int dim = GridView::dimension;
86 static constexpr int dimWorld = GridView::dimensionworld;
87 static constexpr bool isBox = GridGeometry::discMethod == DiscretizationMethod::box;
88 static constexpr int dofCodim = isBox ? dim : 0;
89 static constexpr bool stationaryVelocityField = FluxTraits::hasStationaryVelocityField();
90
91 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
92
93 using Problem = typename GridVolumeVariables::Problem;
94
95 static constexpr bool modelIsCompositional = Detail::hasMoleFraction<typename GridVolumeVariables::VolumeVariables, int, int>() ||
96 Detail::hasMassFraction<typename GridVolumeVariables::VolumeVariables, int, int>();
97
98 static_assert(VolumeVariables::numFluidPhases() >= 1, "Velocity output only makes sense for models with fluid phases.");
99
100public:
101 static constexpr int numFluidPhases = VolumeVariables::numFluidPhases();
102 using VelocityVector = std::vector<Dune::FieldVector<Scalar, dimWorld>>;
103
109 PorousMediumFlowVelocity(const GridVariables& gridVariables)
110 : problem_(gridVariables.curGridVolVars().problem())
111 , gridGeometry_(gridVariables.gridGeometry())
112 , gridVariables_(gridVariables)
113 {
114 // set the number of scvs the vertices are connected to
115 if constexpr (isBox && dim > 1)
116 {
117 // resize to the number of vertices of the grid
118 cellNum_.assign(gridGeometry_.gridView().size(dim), 0);
119
120 for (const auto& element : elements(gridGeometry_.gridView()))
121 for (unsigned int vIdx = 0; vIdx < element.subEntities(dim); ++vIdx)
122 ++cellNum_[gridGeometry_.vertexMapper().subIndex(element, vIdx, dim)];
123 }
124 }
125
129 const Element& element,
130 const FVElementGeometry& fvGeometry,
131 const ElementVolumeVariables& elemVolVars,
132 const ElementFluxVarsCache& elemFluxVarsCache,
133 int phaseIdx) const
134 {
135 using Velocity = typename VelocityVector::value_type;
136
137 const auto geometry = element.geometry();
138 const Dune::GeometryType geomType = geometry.type();
139
140 // the upwind term to be used for the volume flux evaluation
141 auto upwindTerm = [phaseIdx](const auto& volVars) { return volVars.mobility(phaseIdx); };
142
143 if constexpr (isBox && dim == 1)
144 {
145 Velocity tmpVelocity(0.0);
146 tmpVelocity = (geometry.corner(1) - geometry.corner(0));
147 tmpVelocity /= tmpVelocity.two_norm();
148
149 for (auto&& scvf : scvfs(fvGeometry))
150 {
151 if (scvf.boundary())
152 continue;
153
154 // insantiate the flux variables
155 FluxVariables fluxVars;
156 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
157
158 // get the volume flux divided by the area of the
159 // subcontrolvolume face in the reference element
160 Scalar localArea = scvfReferenceArea_(geomType, scvf.index());
161 Scalar flux = fluxVars.advectiveFlux(phaseIdx, upwindTerm) / localArea;
162 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
163 flux /= insideVolVars.extrusionFactor();
164 tmpVelocity *= flux;
165
166 const int eIdxGlobal = gridGeometry_.elementMapper().index(element);
167 velocity[eIdxGlobal] = tmpVelocity;
168 }
169 return;
170 }
171
172 // get the transposed Jacobian of the element mapping
173 const auto refElement = Dune::referenceElement(geometry);
174 const auto& localPos = refElement.position(0, 0);
175 const auto jacobianT2 = geometry.jacobianTransposed(localPos);
176
177 if constexpr (isBox)
178 {
179 using ScvVelocities = Dune::BlockVector<Velocity>;
180 ScvVelocities scvVelocities(fvGeometry.numScv());
181 scvVelocities = 0;
182
183 for (auto&& scvf : scvfs(fvGeometry))
184 {
185 if (scvf.boundary())
186 continue;
187
188 // local position of integration point
189 const auto localPosIP = geometry.local(scvf.ipGlobal());
190
191 // Transformation of the global normal vector to normal vector in the reference element
192 const auto jacobianT1 = geometry.jacobianTransposed(localPosIP);
193 const auto globalNormal = scvf.unitOuterNormal();
194 GlobalPosition localNormal(0);
195 jacobianT1.mv(globalNormal, localNormal);
196 localNormal /= localNormal.two_norm();
197
198 // instantiate the flux variables
199 FluxVariables fluxVars;
200 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
201
202 // get the volume flux divided by the area of the
203 // subcontrolvolume face in the reference element
204 Scalar localArea = scvfReferenceArea_(geomType, scvf.index());
205 Scalar flux = fluxVars.advectiveFlux(phaseIdx, upwindTerm) / localArea;
206 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
207 flux /= insideVolVars.extrusionFactor();
208
209 // transform the volume flux into a velocity vector
210 Velocity tmpVelocity = localNormal;
211 tmpVelocity *= flux;
212
213 scvVelocities[scvf.insideScvIdx()] += tmpVelocity;
214 scvVelocities[scvf.outsideScvIdx()] += tmpVelocity;
215 }
216
217 // transform vertex velocities from local to global coordinates
218 for (auto&& scv : scvs(fvGeometry))
219 {
220 int vIdxGlobal = scv.dofIndex();
221
222 // calculate the subcontrolvolume velocity by the Piola transformation
223 Velocity scvVelocity(0);
224
225 jacobianT2.mtv(scvVelocities[scv.indexInElement()], scvVelocity);
226 scvVelocity /= geometry.integrationElement(localPos)*cellNum_[vIdxGlobal];
227 // add up the wetting phase subcontrolvolume velocities for each vertex
228 velocity[vIdxGlobal] += scvVelocity;
229 }
230 }
231 else
232 {
233 // For the number of scvfs per facet (mpfa) we simply obtain the number of
234 // corners of the first facet as prisms/pyramids are not supported here anyway
235 // -> for prisms/pyramids the number of scvfs would differ from facet to facet
236 static constexpr bool isMpfa = GridGeometry::discMethod == DiscretizationMethod::ccmpfa;
237 const int numScvfsPerFace = isMpfa ? element.template subEntity<1>(0).geometry().corners() : 1;
238
239 if (fvGeometry.numScvf() != element.subEntities(1)*numScvfsPerFace)
240 DUNE_THROW(Dune::NotImplemented, "Velocity output for non-conforming grids");
241
242 if (!geomType.isCube() && !geomType.isSimplex())
243 DUNE_THROW(Dune::NotImplemented, "Velocity output for other geometry types than cube and simplex");
244
245 // first we extract the corner indices for each scv for the CIV method
246 // for network grids there might be multiple intersection with the same geometryInInside
247 // we identify those by the indexInInside for now (assumes conforming grids at branching facets)
248 // here we keep track of them
249 std::vector<bool> handledScvf;
250 if (dim < dimWorld)
251 handledScvf.resize(element.subEntities(1), false);
252
253 // find the local face indices of the scvfs (for conforming meshes)
254 std::vector<unsigned int> scvfIndexInInside(fvGeometry.numScvf());
255 int localScvfIdx = 0;
256 for (const auto& intersection : intersections(gridGeometry_.gridView(), element))
257 {
258 if (dim < dimWorld)
259 if (handledScvf[intersection.indexInInside()])
260 continue;
261
262 if (intersection.neighbor() || intersection.boundary())
263 {
264 for (int i = 0; i < numScvfsPerFace; ++i)
265 scvfIndexInInside[localScvfIdx++] = intersection.indexInInside();
266
267 // for surface and network grids mark that we handled this face
268 if (dim < dimWorld)
269 handledScvf[intersection.indexInInside()] = true;
270 }
271 }
272
273 std::vector<Scalar> scvfFluxes(element.subEntities(1), 0.0);
274 localScvfIdx = 0;
275 for (auto&& scvf : scvfs(fvGeometry))
276 {
277 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
278 if (!scvf.boundary())
279 {
280 FluxVariables fluxVars;
281 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
282 scvfFluxes[scvfIndexInInside[localScvfIdx]] += fluxVars.advectiveFlux(phaseIdx, upwindTerm)/insideVolVars.extrusionFactor();
283 }
284 else
285 {
286 auto bcTypes = problem_.boundaryTypes(element, scvf);
287 if (bcTypes.hasOnlyDirichlet())
288 {
289 FluxVariables fluxVars;
290 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
291 scvfFluxes[scvfIndexInInside[localScvfIdx]] += fluxVars.advectiveFlux(phaseIdx, upwindTerm)/insideVolVars.extrusionFactor();
292 }
293 }
294
295 // increment scvf counter
296 localScvfIdx++;
297 }
298
299 // Correct boundary fluxes in case of Neumann conditions.
300 // In this general setting, it would be very difficult to
301 // calculate correct phase, i.e., volume, fluxes from arbitrary
302 // Neumann conditions. We approximate the Neumann flux by the
303 // flux on the opposite face. For extremely distorted grids this can
304 // lead to unexpected results (but then TPFA also leads to unexpected results).
305 localScvfIdx = 0;
306 for (auto&& scvf : scvfs(fvGeometry))
307 {
308 if (scvf.boundary())
309 {
310 auto bcTypes = problem_.boundaryTypes(element, scvf);
311 if (bcTypes.hasNeumann())
312 {
313 // for stationary velocity fields we can easily compute the correct velocity
314 // this special treatment makes sure that the velocity field is also correct on Neumann boundaries
315 // of tracer problems where the velocity field is given.
316 // (For Dirichlet boundaries no special treatment is necessary.)
317 if (stationaryVelocityField)
318 {
319 const auto flux = AdvectionType::flux(problem_, element, fvGeometry, elemVolVars,
320 scvf, phaseIdx, elemFluxVarsCache);
321
322 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
323 scvfFluxes[scvfIndexInInside[localScvfIdx]] += flux / insideVolVars.extrusionFactor();
324 }
325 else
326 {
327 // check if we have Neumann no flow, we can just use 0
328 const auto neumannFlux = problem_.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
329 using NumEqVector = std::decay_t<decltype(neumannFlux)>;
330 if (Dune::FloatCmp::eq<NumEqVector, Dune::FloatCmp::CmpStyle::absolute>(neumannFlux, NumEqVector(0.0), 1e-30))
331 scvfFluxes[scvfIndexInInside[localScvfIdx]] = 0;
332
333 // otherwise, we try some reconstruction
334 // for cubes
335 else if (dim == 1 || geomType.isCube())
336 {
337 const auto fIdx = scvfIndexInInside[localScvfIdx];
338
339 if constexpr (!modelIsCompositional)
340 {
341 // We assume that the density at the face equals the one at the cell center and reconstruct a volume flux from the Neumann mass flux.
342 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
343 const auto eqIdx = VolumeVariables::Indices::conti0EqIdx + phaseIdx;
344 scvfFluxes[fIdx] += neumannFlux[eqIdx] / insideVolVars.density(phaseIdx) * scvf.area() * insideVolVars.extrusionFactor();
345 }
346 else
347 {
348 // For compositional models, we generally can't reconstruct the volume flow from the Neumann flux (which is a species flux rather
349 // than a phase flux here). Instead, we use the velocity of the opposing face.
350 const auto fIdxOpposite = fIdx%2 ? fIdx-1 : fIdx+1;
351 scvfFluxes[fIdx] = -scvfFluxes[fIdxOpposite];
352 }
353 }
354
355 // for simplices
356 else if (geomType.isSimplex())
357 scvfFluxes[scvfIndexInInside[localScvfIdx]] = 0;
358
359 else
360 DUNE_THROW(Dune::NotImplemented, "Velocity computation at Neumann boundaries for cell-centered and prism/pyramid");
361 }
362 }
363 }
364
365 // increment scvf counter
366 localScvfIdx++;
367 }
368
369
370 Velocity refVelocity;
371 // cubes: On the reference element simply average over opposite fluxes
372 // note that this is equal to a corner velocity interpolation method
373 if (dim == 1 || geomType.isCube())
374 {
375 for (int i = 0; i < dim; i++)
376 refVelocity[i] = 0.5 * (scvfFluxes[2*i + 1] - scvfFluxes[2*i]);
377 }
378 // simplices: Raviart-Thomas-0 interpolation evaluated at the cell center
379 else if (geomType.isSimplex())
380 {
381 for (int dimIdx = 0; dimIdx < dim; dimIdx++)
382 {
383 refVelocity[dimIdx] = -scvfFluxes[dim - 1 - dimIdx];
384 for (int fIdx = 0; fIdx < dim + 1; fIdx++)
385 refVelocity[dimIdx] += scvfFluxes[fIdx]/(dim + 1);
386 }
387 }
388 // 3D prism and pyramids
389 else
390 DUNE_THROW(Dune::NotImplemented, "Velocity computation for cell-centered and prism/pyramid");
391
392 Velocity scvVelocity(0);
393 jacobianT2.mtv(refVelocity, scvVelocity);
394
395 scvVelocity /= geometry.integrationElement(localPos);
396
397 int eIdxGlobal = gridGeometry_.elementMapper().index(element);
398
399 velocity[eIdxGlobal] = scvVelocity;
400
401 } // cell-centered
402 }
403
404private:
405 // The area of a subcontrolvolume face in a reference element.
406 // The 3d non-cube values have been calculated with quadrilateralArea3D
407 // of boxfvelementgeometry.hh.
408 static Scalar scvfReferenceArea_(Dune::GeometryType geomType, int fIdx)
409 {
410 if (dim == 1 || geomType.isCube())
411 {
412 return 1.0/(1 << (dim-1));
413 }
414 else if (geomType.isTriangle())
415 {
416 static const Scalar faceToArea[] = {0.372677996249965,
417 0.372677996249965,
418 0.235702260395516};
419 return faceToArea[fIdx];
420 }
421 else if (geomType.isTetrahedron())
422 {
423 static const Scalar faceToArea[] = {0.102062072615966,
424 0.102062072615966,
425 0.058925565098879,
426 0.102062072615966,
427 0.058925565098879,
428 0.058925565098879};
429 return faceToArea[fIdx];
430 }
431 else if (geomType.isPyramid())
432 {
433 static const Scalar faceToArea[] = {0.130437298687488,
434 0.130437298687488,
435 0.130437298687488,
436 0.130437298687488,
437 0.150923085635624,
438 0.1092906420717,
439 0.1092906420717,
440 0.0781735959970571};
441 return faceToArea[fIdx];
442 }
443 else if (geomType.isPrism())
444 {
445 static const Scalar faceToArea[] = {0.166666666666667,
446 0.166666666666667,
447 0.166666666666667,
448 0.186338998124982,
449 0.186338998124982,
450 0.117851130197758,
451 0.186338998124982,
452 0.186338998124982,
453 0.117851130197758};
454 return faceToArea[fIdx];
455 }
456 else {
457 DUNE_THROW(Dune::NotImplemented, "scvf area for unknown GeometryType");
458 }
459 }
460
461private:
462 const Problem& problem_;
463 const GridGeometry& gridGeometry_;
464 const GridVariables& gridVariables_;
465
466 std::vector<int> cellNum_;
467};
468
469} // end namespace Dumux
470
471#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Element solution classes and factory functions.
The available discretization methods in Dumux.
Definition: adapt.hh:29
Traits of a flux variables type.
Definition: flux/traits.hh:44
static constexpr bool hasStationaryVelocityField()
Definition: flux/traits.hh:45
Velocity computation for implicit (porous media) models.
Definition: velocity.hh:69
void calculateVelocity(VelocityVector &velocity, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVarsCache &elemFluxVarsCache, int phaseIdx) const
Definition: velocity.hh:128
PorousMediumFlowVelocity(const GridVariables &gridVariables)
Constructor initializes the static data with the initial solution.
Definition: velocity.hh:109
std::vector< Dune::FieldVector< Scalar, dimWorld > > VelocityVector
Definition: velocity.hh:102
static constexpr int numFluidPhases
Definition: velocity.hh:101
Defines the flux traits.