25#ifndef DUMUX_POROUSMEDIUMFLOW_VELOCITY_HH
26#define DUMUX_POROUSMEDIUMFLOW_VELOCITY_HH
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>
47template <
typename T,
typename ...Ts>
48using MoleFractionDetector =
decltype(std::declval<T>().moleFraction(std::declval<Ts>()...));
50template<
class T,
typename ...Args>
51static constexpr bool hasMoleFraction()
52{
return Dune::Std::is_detected<MoleFractionDetector, T, Args...>::value; }
54template <
typename T,
typename ...Ts>
55using MassFractionDetector =
decltype(std::declval<T>().massFraction(std::declval<Ts>()...));
57template<
class T,
typename ...Args>
58static constexpr bool hasMassFraction()
59{
return Dune::Std::is_detected<MassFractionDetector, T, Args...>::value; }
68template<
class Gr
idVariables,
class FluxVariables>
71 using GridGeometry =
typename GridVariables::GridGeometry;
72 using FVElementGeometry =
typename GridGeometry::LocalView;
73 using SubControlVolume =
typename GridGeometry::SubControlVolume;
74 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
75 using GridView =
typename GridGeometry::GridView;
76 using Element =
typename GridView::template Codim<0>::Entity;
77 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
78 using ElementFluxVarsCache =
typename GridVariables::GridFluxVariablesCache::LocalView;
79 using VolumeVariables =
typename GridVariables::VolumeVariables;
80 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
81 using FluidSystem =
typename VolumeVariables::FluidSystem;
82 using Scalar =
typename GridVariables::Scalar;
84 using AdvectionType =
typename FluxVariables::AdvectionType;
86 static constexpr int dim = GridView::dimension;
87 static constexpr int dimWorld = GridView::dimensionworld;
89 static constexpr int dofCodim = isBox ? dim : 0;
92 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
94 using Problem =
typename GridVolumeVariables::Problem;
96 static constexpr bool modelIsCompositional = Detail::hasMoleFraction<typename GridVolumeVariables::VolumeVariables, int, int>() ||
97 Detail::hasMassFraction<typename GridVolumeVariables::VolumeVariables, int, int>();
99 static_assert(VolumeVariables::numFluidPhases() >= 1,
"Velocity output only makes sense for models with fluid phases.");
111 : problem_(gridVariables.curGridVolVars().problem())
112 , gridGeometry_(gridVariables.gridGeometry())
113 , gridVariables_(gridVariables)
116 if constexpr (isBox && dim > 1)
119 cellNum_.assign(gridGeometry_.gridView().size(dim), 0);
121 for (
const auto& element : elements(gridGeometry_.gridView()))
122 for (
unsigned int vIdx = 0; vIdx < element.subEntities(dim); ++vIdx)
123 ++cellNum_[gridGeometry_.vertexMapper().subIndex(element, vIdx, dim)];
130 const Element& element,
131 const FVElementGeometry& fvGeometry,
132 const ElementVolumeVariables& elemVolVars,
133 const ElementFluxVarsCache& elemFluxVarsCache,
136 using Velocity =
typename VelocityVector::value_type;
138 const auto geometry = element.geometry();
139 const Dune::GeometryType geomType = geometry.type();
142 auto upwindTerm = [phaseIdx](
const auto& volVars) {
return volVars.mobility(phaseIdx); };
144 if constexpr (isBox && dim == 1)
146 Velocity tmpVelocity(0.0);
147 tmpVelocity = (geometry.corner(1) - geometry.corner(0));
148 tmpVelocity /= tmpVelocity.two_norm();
150 for (
auto&& scvf : scvfs(fvGeometry))
156 FluxVariables fluxVars;
157 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
161 Scalar localArea = scvfReferenceArea_(geomType, scvf.index());
162 Scalar flux = fluxVars.advectiveFlux(phaseIdx, upwindTerm) / localArea;
163 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
164 flux /= insideVolVars.extrusionFactor();
167 const int eIdxGlobal = gridGeometry_.elementMapper().index(element);
168 velocity[eIdxGlobal] = tmpVelocity;
174 const auto refElement = referenceElement(geometry);
175 const auto& localPos = refElement.position(0, 0);
176 const auto jacobianT2 = geometry.jacobianTransposed(localPos);
180 using ScvVelocities = Dune::BlockVector<Velocity>;
181 ScvVelocities scvVelocities(fvGeometry.numScv());
184 for (
auto&& scvf : scvfs(fvGeometry))
190 const auto localPosIP = geometry.local(scvf.ipGlobal());
193 const auto jacobianT1 = geometry.jacobianTransposed(localPosIP);
194 const auto globalNormal = scvf.unitOuterNormal();
195 GlobalPosition localNormal(0);
196 jacobianT1.mv(globalNormal, localNormal);
197 localNormal /= localNormal.two_norm();
200 FluxVariables fluxVars;
201 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
205 Scalar localArea = scvfReferenceArea_(geomType, scvf.index());
206 Scalar flux = fluxVars.advectiveFlux(phaseIdx, upwindTerm) / localArea;
207 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
208 flux /= insideVolVars.extrusionFactor();
211 Velocity tmpVelocity = localNormal;
214 scvVelocities[scvf.insideScvIdx()] += tmpVelocity;
215 scvVelocities[scvf.outsideScvIdx()] += tmpVelocity;
219 for (
auto&& scv : scvs(fvGeometry))
221 int vIdxGlobal = scv.dofIndex();
224 Velocity scvVelocity(0);
226 jacobianT2.mtv(scvVelocities[scv.indexInElement()], scvVelocity);
227 scvVelocity /= geometry.integrationElement(localPos)*cellNum_[vIdxGlobal];
229 velocity[vIdxGlobal] += scvVelocity;
238 const int numScvfsPerFace = isMpfa ? element.template subEntity<1>(0).geometry().corners() : 1;
240 if (fvGeometry.numScvf() != element.subEntities(1)*numScvfsPerFace)
241 DUNE_THROW(Dune::NotImplemented,
"Velocity output for non-conforming grids");
243 if (!geomType.isCube() && !geomType.isSimplex())
244 DUNE_THROW(Dune::NotImplemented,
"Velocity output for other geometry types than cube and simplex");
250 std::vector<bool> handledScvf;
252 handledScvf.resize(element.subEntities(1),
false);
255 std::vector<unsigned int> scvfIndexInInside(fvGeometry.numScvf());
256 int localScvfIdx = 0;
257 for (
const auto& intersection : intersections(gridGeometry_.gridView(), element))
260 if (handledScvf[intersection.indexInInside()])
263 if (intersection.neighbor() || intersection.boundary())
265 for (
int i = 0; i < numScvfsPerFace; ++i)
266 scvfIndexInInside[localScvfIdx++] = intersection.indexInInside();
270 handledScvf[intersection.indexInInside()] =
true;
274 std::vector<Scalar> scvfFluxes(element.subEntities(1), 0.0);
276 for (
auto&& scvf : scvfs(fvGeometry))
278 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
279 if (!scvf.boundary())
281 FluxVariables fluxVars;
282 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
283 scvfFluxes[scvfIndexInInside[localScvfIdx]] += fluxVars.advectiveFlux(phaseIdx, upwindTerm)/insideVolVars.extrusionFactor();
287 auto bcTypes = problem_.boundaryTypes(element, scvf);
288 if (bcTypes.hasOnlyDirichlet())
290 FluxVariables fluxVars;
291 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
292 scvfFluxes[scvfIndexInInside[localScvfIdx]] += fluxVars.advectiveFlux(phaseIdx, upwindTerm)/insideVolVars.extrusionFactor();
307 for (
auto&& scvf : scvfs(fvGeometry))
311 auto bcTypes = problem_.boundaryTypes(element, scvf);
312 if (bcTypes.hasNeumann())
318 if (stationaryVelocityField)
320 const auto flux = AdvectionType::flux(problem_, element, fvGeometry, elemVolVars,
321 scvf, phaseIdx, elemFluxVarsCache);
323 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
324 scvfFluxes[scvfIndexInInside[localScvfIdx]] += flux / insideVolVars.extrusionFactor();
329 const auto neumannFlux = problem_.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
330 using NumEqVector = std::decay_t<
decltype(neumannFlux)>;
331 if (Dune::FloatCmp::eq<NumEqVector, Dune::FloatCmp::CmpStyle::absolute>(neumannFlux,
NumEqVector(0.0), 1e-30))
332 scvfFluxes[scvfIndexInInside[localScvfIdx]] = 0;
336 else if (dim == 1 || geomType.isCube())
338 const auto fIdx = scvfIndexInInside[localScvfIdx];
340 if constexpr (!modelIsCompositional)
343 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
344 const auto eqIdx = VolumeVariables::Indices::conti0EqIdx + phaseIdx;
345 scvfFluxes[fIdx] += neumannFlux[eqIdx] / insideVolVars.density(phaseIdx) * scvf.area() * insideVolVars.extrusionFactor();
351 const auto fIdxOpposite = fIdx%2 ? fIdx-1 : fIdx+1;
352 scvfFluxes[fIdx] = -scvfFluxes[fIdxOpposite];
357 else if (geomType.isSimplex())
358 scvfFluxes[scvfIndexInInside[localScvfIdx]] = 0;
361 DUNE_THROW(Dune::NotImplemented,
"Velocity computation at Neumann boundaries for cell-centered and prism/pyramid");
371 Velocity refVelocity;
374 if (dim == 1 || geomType.isCube())
376 for (
int i = 0; i < dim; i++)
377 refVelocity[i] = 0.5 * (scvfFluxes[2*i + 1] - scvfFluxes[2*i]);
380 else if (geomType.isSimplex())
382 for (
int dimIdx = 0; dimIdx < dim; dimIdx++)
384 refVelocity[dimIdx] = -scvfFluxes[dim - 1 - dimIdx];
385 for (
int fIdx = 0; fIdx < dim + 1; fIdx++)
386 refVelocity[dimIdx] += scvfFluxes[fIdx]/(dim + 1);
391 DUNE_THROW(Dune::NotImplemented,
"Velocity computation for cell-centered and prism/pyramid");
393 Velocity scvVelocity(0);
394 jacobianT2.mtv(refVelocity, scvVelocity);
396 scvVelocity /= geometry.integrationElement(localPos);
398 int eIdxGlobal = gridGeometry_.elementMapper().index(element);
400 velocity[eIdxGlobal] = scvVelocity;
409 static Scalar scvfReferenceArea_(Dune::GeometryType geomType,
int fIdx)
411 if (dim == 1 || geomType.isCube())
413 return 1.0/(1 << (dim-1));
415 else if (geomType.isTriangle())
417 static const Scalar faceToArea[] = {0.372677996249965,
420 return faceToArea[fIdx];
422 else if (geomType.isTetrahedron())
424 static const Scalar faceToArea[] = {0.102062072615966,
430 return faceToArea[fIdx];
432 else if (geomType.isPyramid())
434 static const Scalar faceToArea[] = {0.130437298687488,
442 return faceToArea[fIdx];
444 else if (geomType.isPrism())
446 static const Scalar faceToArea[] = {0.166666666666667,
455 return faceToArea[fIdx];
458 DUNE_THROW(Dune::NotImplemented,
"scvf area for unknown GeometryType");
463 const Problem& problem_;
464 const GridGeometry& gridGeometry_;
465 const GridVariables& gridVariables_;
467 std::vector<int> cellNum_;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Element solution classes and factory functions.
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
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:70
void calculateVelocity(VelocityVector &velocity, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVarsCache &elemFluxVarsCache, int phaseIdx) const
Definition: velocity.hh:129
PorousMediumFlowVelocity(const GridVariables &gridVariables)
Constructor initializes the static data with the initial solution.
Definition: velocity.hh:110
std::vector< Dune::FieldVector< Scalar, dimWorld > > VelocityVector
Definition: velocity.hh:103
static constexpr int numFluidPhases
Definition: velocity.hh:102