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