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/geometry/type.hh>
22#include <dune/geometry/referenceelements.hh>
36template <
typename T,
typename ...Ts>
37using MoleFractionDetector =
decltype(std::declval<T>().moleFraction(std::declval<Ts>()...));
39template<
class T,
typename ...Args>
40static constexpr bool hasMoleFraction()
41{
return Dune::Std::is_detected<MoleFractionDetector, T, Args...>::value; }
43template <
typename T,
typename ...Ts>
44using MassFractionDetector =
decltype(std::declval<T>().massFraction(std::declval<Ts>()...));
46template<
class T,
typename ...Args>
47static constexpr bool hasMassFraction()
48{
return Dune::Std::is_detected<MassFractionDetector, T, Args...>::value; }
57template<
class Gr
idVariables,
class FluxVariables>
60 using GridGeometry =
typename GridVariables::GridGeometry;
62 using FVElementGeometry =
typename GridGeometry::LocalView;
63 using SubControlVolume =
typename GridGeometry::SubControlVolume;
64 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
65 using GridView =
typename GridGeometry::GridView;
66 using Element =
typename GridView::template Codim<0>::Entity;
67 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
68 using ElementFluxVarsCache =
typename GridVariables::GridFluxVariablesCache::LocalView;
69 using VolumeVariables =
typename GridVariables::VolumeVariables;
70 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
71 using FluidSystem =
typename VolumeVariables::FluidSystem;
72 using Scalar =
typename GridVariables::Scalar;
74 using AdvectionType =
typename FluxVariables::AdvectionType;
76 static constexpr int dim = GridView::dimension;
77 static constexpr int dimWorld = GridView::dimensionworld;
82 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
84 using Problem =
typename GridVolumeVariables::Problem;
86 static constexpr bool modelIsCompositional = Detail::hasMoleFraction<typename GridVolumeVariables::VolumeVariables, int, int>() ||
87 Detail::hasMassFraction<typename GridVolumeVariables::VolumeVariables, int, int>();
89 static_assert(VolumeVariables::numFluidPhases() >= 1,
"Velocity output only makes sense for models with fluid phases.");
101 : problem_(gridVariables.curGridVolVars().problem())
102 , gridGeometry_(gridVariables.gridGeometry())
103 , gridVariables_(gridVariables)
106 if constexpr (isBox && dim > 1)
109 cellNum_.assign(gridGeometry_.gridView().size(dim), 0);
111 for (
const auto& element : elements(gridGeometry_.gridView()))
112 for (
unsigned int vIdx = 0; vIdx < element.subEntities(dim); ++vIdx)
113 ++cellNum_[gridGeometry_.vertexMapper().subIndex(element, vIdx, dim)];
120 const Element& element,
121 const FVElementGeometry& fvGeometry,
122 const ElementVolumeVariables& elemVolVars,
123 const ElementFluxVarsCache& elemFluxVarsCache,
126 using Velocity =
typename VelocityVector::value_type;
128 const auto geometry = element.geometry();
129 const Dune::GeometryType geomType = geometry.type();
132 auto upwindTerm = [phaseIdx](
const auto& volVars) {
return volVars.mobility(phaseIdx); };
134 if constexpr (isBox && dim == 1)
136 Velocity tmpVelocity(0.0);
137 tmpVelocity = (geometry.corner(1) - geometry.corner(0));
138 tmpVelocity /= tmpVelocity.two_norm();
140 for (
auto&& scvf : scvfs(fvGeometry))
146 FluxVariables fluxVars;
147 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
151 Scalar localArea = scvfReferenceArea_(geomType, scvf.index());
152 Scalar flux = fluxVars.advectiveFlux(phaseIdx, upwindTerm) / localArea;
153 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
154 flux /= insideVolVars.extrusionFactor() * Extrusion::area(fvGeometry, scvf) / scvf.area();
157 const int eIdxGlobal = gridGeometry_.elementMapper().index(element);
158 velocity[eIdxGlobal] = tmpVelocity;
164 using Dune::referenceElement;
165 const auto refElement = referenceElement(geometry);
166 const auto& localPos = refElement.position(0, 0);
167 const auto jacobianT2 = geometry.jacobianTransposed(localPos);
171 using ScvVelocities = Dune::BlockVector<Velocity>;
172 ScvVelocities scvVelocities(fvGeometry.numScv());
175 for (
auto&& scvf : scvfs(fvGeometry))
181 const auto localPosIP = geometry.local(scvf.ipGlobal());
184 const auto jacobianT1 = geometry.jacobianTransposed(localPosIP);
185 const auto globalNormal = scvf.unitOuterNormal();
186 GlobalPosition localNormal(0);
187 jacobianT1.mv(globalNormal, localNormal);
188 localNormal /= localNormal.two_norm();
191 FluxVariables fluxVars;
192 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
196 Scalar localArea = scvfReferenceArea_(geomType, scvf.index());
197 Scalar flux = fluxVars.advectiveFlux(phaseIdx, upwindTerm) / localArea;
198 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
199 flux /= insideVolVars.extrusionFactor() * Extrusion::area(fvGeometry, scvf) / scvf.area();
202 Velocity tmpVelocity = localNormal;
205 scvVelocities[scvf.insideScvIdx()] += tmpVelocity;
206 scvVelocities[scvf.outsideScvIdx()] += tmpVelocity;
210 for (
auto&& scv : scvs(fvGeometry))
212 int vIdxGlobal = scv.dofIndex();
215 Velocity scvVelocity(0);
217 jacobianT2.mtv(scvVelocities[scv.indexInElement()], scvVelocity);
218 scvVelocity /= geometry.integrationElement(localPos)*cellNum_[vIdxGlobal];
220 velocity[vIdxGlobal] += scvVelocity;
223 else if constexpr (isDiamond)
224 DUNE_THROW(Dune::NotImplemented,
"Velocity output with diamond discretization");
231 const int numScvfsPerFace = isMpfa ? element.template subEntity<1>(0).geometry().corners() : 1;
233 if (fvGeometry.numScvf() != element.subEntities(1)*numScvfsPerFace)
234 DUNE_THROW(Dune::NotImplemented,
"Velocity output for non-conforming grids");
236 if (!geomType.isCube() && !geomType.isSimplex())
237 DUNE_THROW(Dune::NotImplemented,
"Velocity output for other geometry types than cube and simplex");
243 std::vector<bool> handledScvf;
245 handledScvf.resize(element.subEntities(1),
false);
248 std::vector<unsigned int> scvfIndexInInside(fvGeometry.numScvf());
249 int localScvfIdx = 0;
250 for (
const auto& intersection : intersections(gridGeometry_.gridView(), element))
253 if (handledScvf[intersection.indexInInside()])
256 if (intersection.neighbor() || intersection.boundary())
258 for (
int i = 0; i < numScvfsPerFace; ++i)
259 scvfIndexInInside[localScvfIdx++] = intersection.indexInInside();
263 handledScvf[intersection.indexInInside()] =
true;
267 std::vector<Scalar> scvfFluxes(element.subEntities(1), 0.0);
269 for (
auto&& scvf : scvfs(fvGeometry))
271 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
272 if (!scvf.boundary())
274 FluxVariables fluxVars;
275 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
276 scvfFluxes[scvfIndexInInside[localScvfIdx]] += fluxVars.advectiveFlux(phaseIdx, upwindTerm)/insideVolVars.extrusionFactor();
280 auto bcTypes = problem_.boundaryTypes(element, scvf);
281 if (bcTypes.hasOnlyDirichlet())
283 FluxVariables fluxVars;
284 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
285 scvfFluxes[scvfIndexInInside[localScvfIdx]] += fluxVars.advectiveFlux(phaseIdx, upwindTerm)/insideVolVars.extrusionFactor();
300 for (
auto&& scvf : scvfs(fvGeometry))
304 auto bcTypes = problem_.boundaryTypes(element, scvf);
305 if (bcTypes.hasNeumann())
311 if (stationaryVelocityField)
313 const auto flux = AdvectionType::flux(problem_, element, fvGeometry, elemVolVars,
314 scvf, phaseIdx, elemFluxVarsCache);
316 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
317 scvfFluxes[scvfIndexInInside[localScvfIdx]] += flux / insideVolVars.extrusionFactor();
322 const auto neumannFlux = problem_.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
323 using NumEqVector = std::decay_t<
decltype(neumannFlux)>;
324 if (Dune::FloatCmp::eq<NumEqVector, Dune::FloatCmp::CmpStyle::absolute>(neumannFlux,
NumEqVector(0.0), 1e-30))
325 scvfFluxes[scvfIndexInInside[localScvfIdx]] = 0;
329 else if (dim == 1 || geomType.isCube())
331 const auto fIdx = scvfIndexInInside[localScvfIdx];
333 if constexpr (!modelIsCompositional)
336 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
337 const auto eqIdx = VolumeVariables::Indices::conti0EqIdx + phaseIdx;
338 scvfFluxes[fIdx] += neumannFlux[eqIdx] / insideVolVars.density(phaseIdx) * scvf.area() * insideVolVars.extrusionFactor();
344 const auto fIdxOpposite = fIdx%2 ? fIdx-1 : fIdx+1;
345 scvfFluxes[fIdx] = -scvfFluxes[fIdxOpposite];
350 else if (geomType.isSimplex())
351 scvfFluxes[scvfIndexInInside[localScvfIdx]] = 0;
354 DUNE_THROW(Dune::NotImplemented,
"Velocity computation at Neumann boundaries for cell-centered and prism/pyramid");
364 Velocity refVelocity;
367 if (dim == 1 || geomType.isCube())
369 for (
int i = 0; i < dim; i++)
370 refVelocity[i] = 0.5 * (scvfFluxes[2*i + 1] - scvfFluxes[2*i]);
373 else if (geomType.isSimplex())
375 for (
int dimIdx = 0; dimIdx < dim; dimIdx++)
377 refVelocity[dimIdx] = -scvfFluxes[dim - 1 - dimIdx];
378 for (
int fIdx = 0; fIdx < dim + 1; fIdx++)
379 refVelocity[dimIdx] += scvfFluxes[fIdx]/(dim + 1);
384 DUNE_THROW(Dune::NotImplemented,
"Velocity computation for cell-centered and prism/pyramid");
386 Velocity scvVelocity(0);
387 jacobianT2.mtv(refVelocity, scvVelocity);
389 scvVelocity /= geometry.integrationElement(localPos);
391 int eIdxGlobal = gridGeometry_.elementMapper().index(element);
393 velocity[eIdxGlobal] = scvVelocity;
402 static Scalar scvfReferenceArea_(Dune::GeometryType geomType,
int fIdx)
404 if (dim == 1 || geomType.isCube())
406 return 1.0/(1 << (dim-1));
408 else if (geomType.isTriangle())
410 static const Scalar faceToArea[] = {0.372677996249965,
413 return faceToArea[fIdx];
415 else if (geomType.isTetrahedron())
417 static const Scalar faceToArea[] = {0.102062072615966,
423 return faceToArea[fIdx];
425 else if (geomType.isPyramid())
427 static const Scalar faceToArea[] = {0.130437298687488,
435 return faceToArea[fIdx];
437 else if (geomType.isPrism())
439 static const Scalar faceToArea[] = {0.166666666666667,
448 return faceToArea[fIdx];
451 DUNE_THROW(Dune::NotImplemented,
"scvf area for unknown GeometryType");
456 const Problem& problem_;
457 const GridGeometry& gridGeometry_;
458 const GridVariables& gridVariables_;
460 std::vector<int> cellNum_;
Velocity computation for implicit (porous media) models.
Definition: velocity.hh:59
void calculateVelocity(VelocityVector &velocity, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVarsCache &elemFluxVarsCache, int phaseIdx) const
Definition: velocity.hh:119
PorousMediumFlowVelocity(const GridVariables &gridVariables)
Constructor initializes the static data with the initial solution.
Definition: velocity.hh:100
std::vector< Dune::FieldVector< Scalar, dimWorld > > VelocityVector
Definition: velocity.hh:93
static constexpr int numFluidPhases
Definition: velocity.hh:92
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
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