13#ifndef DUMUX_POROUSMEDIUMFLOW_VELOCITY_HH
14#define DUMUX_POROUSMEDIUMFLOW_VELOCITY_HH
19#include <dune/common/fvector.hh>
20#include <dune/common/float_cmp.hh>
21#include <dune/common/reservedvector.hh>
22#include <dune/istl/bvector.hh>
23#include <dune/geometry/type.hh>
24#include <dune/geometry/referenceelements.hh>
38template <
typename T,
typename ...Ts>
39using MoleFractionDetector =
decltype(std::declval<T>().moleFraction(std::declval<Ts>()...));
41template<
class T,
typename ...Args>
42static constexpr bool hasMoleFraction()
43{
return Dune::Std::is_detected<MoleFractionDetector, T, Args...>::value; }
45template <
typename T,
typename ...Ts>
46using MassFractionDetector =
decltype(std::declval<T>().massFraction(std::declval<Ts>()...));
48template<
class T,
typename ...Args>
49static constexpr bool hasMassFraction()
50{
return Dune::Std::is_detected<MassFractionDetector, T, Args...>::value; }
59template<
class Gr
idVariables,
class FluxVariables>
62 using GridGeometry =
typename GridVariables::GridGeometry;
64 using FVElementGeometry =
typename GridGeometry::LocalView;
65 using SubControlVolume =
typename GridGeometry::SubControlVolume;
66 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
67 using GridView =
typename GridGeometry::GridView;
68 using Element =
typename GridView::template Codim<0>::Entity;
69 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
70 using ElementFluxVarsCache =
typename GridVariables::GridFluxVariablesCache::LocalView;
71 using VolumeVariables =
typename GridVariables::VolumeVariables;
72 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
73 using FluidSystem =
typename VolumeVariables::FluidSystem;
74 using Scalar =
typename GridVariables::Scalar;
76 using AdvectionType =
typename FluxVariables::AdvectionType;
78 static constexpr int dim = GridView::dimension;
79 static constexpr int dimWorld = GridView::dimensionworld;
84 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
86 using Problem =
typename GridVolumeVariables::Problem;
88 static constexpr bool modelIsCompositional = Detail::hasMoleFraction<typename GridVolumeVariables::VolumeVariables, int, int>() ||
89 Detail::hasMassFraction<typename GridVolumeVariables::VolumeVariables, int, int>();
91 static_assert(VolumeVariables::numFluidPhases() >= 1,
"Velocity output only makes sense for models with fluid phases.");
93 using Velocity = Dune::FieldVector<Scalar, dimWorld>;
94 using ReferenceElementVelocity = Dune::FieldVector<Scalar, dim>;
105 : problem_(gridVariables.curGridVolVars().problem())
106 , gridGeometry_(gridVariables.gridGeometry())
107 , gridVariables_(gridVariables)
110 if constexpr (isBox && dim > 1)
113 cellNum_.assign(gridGeometry_.gridView().size(dim), 0);
115 for (
const auto& element : elements(gridGeometry_.gridView()))
116 for (
unsigned int vIdx = 0; vIdx < element.subEntities(dim); ++vIdx)
117 ++cellNum_[gridGeometry_.vertexMapper().subIndex(element, vIdx, dim)];
124 const Element& element,
125 const FVElementGeometry& fvGeometry,
126 const ElementVolumeVariables& elemVolVars,
127 const ElementFluxVarsCache& elemFluxVarsCache,
130 const auto geometry = element.geometry();
131 const Dune::GeometryType geomType = geometry.type();
134 auto upwindTerm = [phaseIdx](
const auto& volVars) {
return volVars.mobility(phaseIdx); };
136 if constexpr (isBox && dim == 1)
138 Velocity tmpVelocity(0.0);
139 tmpVelocity = (geometry.corner(1) - geometry.corner(0));
140 tmpVelocity /= tmpVelocity.two_norm();
142 for (
auto&& scvf : scvfs(fvGeometry))
148 FluxVariables fluxVars;
149 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
153 Scalar localArea = scvfReferenceArea_(geomType, scvf.index());
154 Scalar flux = fluxVars.advectiveFlux(phaseIdx, upwindTerm) / localArea;
155 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
156 flux /= insideVolVars.extrusionFactor() * Extrusion::area(fvGeometry, scvf) / scvf.area();
159 const int eIdxGlobal = gridGeometry_.elementMapper().index(element);
160 velocity[eIdxGlobal] = tmpVelocity;
166 using Dune::referenceElement;
167 const auto localPos = referenceElement(geometry).position(0, 0);
168 const auto jacobianT2 = geometry.jacobianTransposed(localPos);
172 using ScvVelocities = Dune::BlockVector<Velocity>;
173 ScvVelocities scvVelocities(fvGeometry.numScv());
176 for (
auto&& scvf : scvfs(fvGeometry))
182 const auto localPosIP = geometry.local(scvf.ipGlobal());
185 const auto jacobianT1 = geometry.jacobianTransposed(localPosIP);
186 const auto globalNormal = scvf.unitOuterNormal();
187 GlobalPosition localNormal(0);
188 jacobianT1.mv(globalNormal, localNormal);
189 localNormal /= localNormal.two_norm();
192 FluxVariables fluxVars;
193 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
197 Scalar localArea = scvfReferenceArea_(geomType, scvf.index());
198 Scalar flux = fluxVars.advectiveFlux(phaseIdx, upwindTerm) / localArea;
199 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
200 flux /= insideVolVars.extrusionFactor() * Extrusion::area(fvGeometry, scvf) / scvf.area();
203 Velocity tmpVelocity = localNormal;
206 scvVelocities[scvf.insideScvIdx()] += tmpVelocity;
207 scvVelocities[scvf.outsideScvIdx()] += tmpVelocity;
211 for (
auto&& scv :
scvs(fvGeometry))
213 int vIdxGlobal = scv.dofIndex();
216 Velocity scvVelocity(0);
218 jacobianT2.mtv(scvVelocities[scv.indexInElement()], scvVelocity);
219 scvVelocity /= geometry.integrationElement(localPos)*cellNum_[vIdxGlobal];
221 velocity[vIdxGlobal] += scvVelocity;
224 else if constexpr (isDiamond)
225 DUNE_THROW(Dune::NotImplemented,
"Velocity output with diamond discretization");
232 const int numScvfsPerFace = isMpfa ? element.template subEntity<1>(0).geometry().corners() : 1;
234 if (fvGeometry.numScvf() != element.subEntities(1)*numScvfsPerFace)
235 DUNE_THROW(Dune::NotImplemented,
"Velocity output for non-conforming grids");
237 if (!geomType.isCube() && !geomType.isSimplex())
238 DUNE_THROW(Dune::NotImplemented,
"Velocity output for other geometry types than cube and simplex");
244 std::vector<bool> handledScvf;
246 handledScvf.resize(element.subEntities(1),
false);
249 std::vector<unsigned int> scvfIndexInInside(fvGeometry.numScvf());
250 int localScvfIdx = 0;
251 for (
const auto& intersection : intersections(gridGeometry_.gridView(), element))
254 if (handledScvf[intersection.indexInInside()])
257 if (intersection.neighbor() || intersection.boundary())
259 for (
int i = 0; i < numScvfsPerFace; ++i)
260 scvfIndexInInside[localScvfIdx++] = intersection.indexInInside();
264 handledScvf[intersection.indexInInside()] =
true;
268 Dune::ReservedVector<Scalar, 2*dim> scvfFluxes(element.subEntities(1), 0.0);
270 for (
auto&& scvf : scvfs(fvGeometry))
272 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
273 if (!scvf.boundary())
275 FluxVariables fluxVars;
276 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
277 scvfFluxes[scvfIndexInInside[localScvfIdx]] += fluxVars.advectiveFlux(phaseIdx, upwindTerm)/insideVolVars.extrusionFactor();
281 auto bcTypes = problem_.boundaryTypes(element, scvf);
282 if (bcTypes.hasOnlyDirichlet())
284 FluxVariables fluxVars;
285 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
286 scvfFluxes[scvfIndexInInside[localScvfIdx]] += fluxVars.advectiveFlux(phaseIdx, upwindTerm)/insideVolVars.extrusionFactor();
301 for (
auto&& scvf : scvfs(fvGeometry))
305 auto bcTypes = problem_.boundaryTypes(element, scvf);
306 if (bcTypes.hasNeumann())
312 if (stationaryVelocityField)
314 const auto flux = AdvectionType::flux(problem_, element, fvGeometry, elemVolVars,
315 scvf, phaseIdx, elemFluxVarsCache);
317 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
318 scvfFluxes[scvfIndexInInside[localScvfIdx]] += flux / insideVolVars.extrusionFactor();
323 const auto neumannFlux = problem_.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
324 using NumEqVector = std::decay_t<
decltype(neumannFlux)>;
325 if (Dune::FloatCmp::eq<NumEqVector, Dune::FloatCmp::CmpStyle::absolute>(neumannFlux,
NumEqVector(0.0), 1e-30))
326 scvfFluxes[scvfIndexInInside[localScvfIdx]] = 0;
330 else if (dim == 1 || geomType.isCube())
332 const auto fIdx = scvfIndexInInside[localScvfIdx];
334 if constexpr (!modelIsCompositional)
337 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
338 const auto eqIdx = VolumeVariables::Indices::conti0EqIdx + phaseIdx;
339 scvfFluxes[fIdx] += neumannFlux[eqIdx] / insideVolVars.density(phaseIdx) * scvf.area();
345 const auto fIdxOpposite = fIdx%2 ? fIdx-1 : fIdx+1;
346 scvfFluxes[fIdx] = -scvfFluxes[fIdxOpposite];
351 else if (geomType.isSimplex())
352 scvfFluxes[scvfIndexInInside[localScvfIdx]] = 0;
355 DUNE_THROW(Dune::NotImplemented,
"Velocity computation at Neumann boundaries for cell-centered and prism/pyramid");
365 ReferenceElementVelocity refVelocity;
368 if (dim == 1 || geomType.isCube())
370 assert(scvfFluxes.size() == 2*dim);
371 for (
int i = 0; i < dim; ++i)
372 refVelocity[i] = 0.5 * (scvfFluxes[2*i + 1] - scvfFluxes[2*i]);
376 else if (geomType.isSimplex())
378 assert(scvfFluxes.size() == dim+1);
379 for (
int i = 0; i < dim; ++i)
381 refVelocity[i] = -scvfFluxes[dim - 1 - i];
382 for (
int fIdx = 0; fIdx < dim + 1; ++fIdx)
383 refVelocity[i] += scvfFluxes[fIdx]/(dim + 1);
389 DUNE_THROW(Dune::NotImplemented,
"Velocity computation for cell-centered and prism/pyramid");
391 Velocity scvVelocity(0);
392 jacobianT2.mtv(refVelocity, scvVelocity);
394 scvVelocity /= geometry.integrationElement(localPos);
396 int eIdxGlobal = gridGeometry_.elementMapper().index(element);
398 velocity[eIdxGlobal] = scvVelocity;
407 static Scalar scvfReferenceArea_(Dune::GeometryType geomType,
int fIdx)
409 if (dim == 1 || geomType.isCube())
411 return 1.0/(1 << (dim-1));
413 else if (geomType.isTriangle())
415 static const Scalar faceToArea[] = {0.372677996249965,
418 return faceToArea[fIdx];
420 else if (geomType.isTetrahedron())
422 static const Scalar faceToArea[] = {0.102062072615966,
428 return faceToArea[fIdx];
430 else if (geomType.isPyramid())
432 static const Scalar faceToArea[] = {0.130437298687488,
440 return faceToArea[fIdx];
442 else if (geomType.isPrism())
444 static const Scalar faceToArea[] = {0.166666666666667,
453 return faceToArea[fIdx];
456 DUNE_THROW(Dune::NotImplemented,
"scvf area for unknown GeometryType");
461 const Problem& problem_;
462 const GridGeometry& gridGeometry_;
463 const GridVariables& gridVariables_;
465 std::vector<int> cellNum_;
Velocity computation for implicit (porous media) models.
Definition: velocity.hh:61
void calculateVelocity(VelocityVector &velocity, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVarsCache &elemFluxVarsCache, int phaseIdx) const
Definition: velocity.hh:123
PorousMediumFlowVelocity(const GridVariables &gridVariables)
Constructor initializes the static data with the initial solution.
Definition: velocity.hh:104
std::vector< Velocity > VelocityVector
Definition: velocity.hh:97
static constexpr int numFluidPhases
Definition: velocity.hh:96
Element solution classes and factory functions.
Helper classes to compute the integration elements.
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:34
The available discretization methods in Dumux.
constexpr CCMpfa ccmpfa
Definition: method.hh:146
constexpr FCDiamond fcdiamond
Definition: method.hh:152
constexpr Box box
Definition: method.hh:147
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
std::ranges::range auto scvs(const FVElementGeometry &fvGeometry, const LocalDof &localDof)
Definition: localdof.hh:79
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Traits of a flux variables type.
Definition: flux/traits.hh:32
static constexpr bool hasStationaryVelocityField()
Definition: flux/traits.hh:33