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