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>
48template <
typename T,
typename ...Ts>
49using MoleFractionDetector =
decltype(std::declval<T>().moleFraction(std::declval<Ts>()...));
51template<
class T,
typename ...Args>
52static constexpr bool hasMoleFraction()
53{
return Dune::Std::is_detected<MoleFractionDetector, T, Args...>::value; }
55template <
typename T,
typename ...Ts>
56using MassFractionDetector =
decltype(std::declval<T>().massFraction(std::declval<Ts>()...));
58template<
class T,
typename ...Args>
59static constexpr bool hasMassFraction()
60{
return Dune::Std::is_detected<MassFractionDetector, T, Args...>::value; }
69template<
class Gr
idVariables,
class FluxVariables>
72 using GridGeometry =
typename GridVariables::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;
86 using AdvectionType =
typename FluxVariables::AdvectionType;
88 static constexpr int dim = GridView::dimension;
89 static constexpr int dimWorld = GridView::dimensionworld;
91 static constexpr int dofCodim = isBox ? dim : 0;
94 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
96 using Problem =
typename GridVolumeVariables::Problem;
98 static constexpr bool modelIsCompositional = Detail::hasMoleFraction<typename GridVolumeVariables::VolumeVariables, int, int>() ||
99 Detail::hasMassFraction<typename GridVolumeVariables::VolumeVariables, int, int>();
101 static_assert(VolumeVariables::numFluidPhases() >= 1,
"Velocity output only makes sense for models with fluid phases.");
113 : problem_(gridVariables.curGridVolVars().problem())
114 , gridGeometry_(gridVariables.gridGeometry())
115 , gridVariables_(gridVariables)
118 if constexpr (isBox && dim > 1)
121 cellNum_.assign(gridGeometry_.gridView().size(dim), 0);
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)];
132 const Element& element,
133 const FVElementGeometry& fvGeometry,
134 const ElementVolumeVariables& elemVolVars,
135 const ElementFluxVarsCache& elemFluxVarsCache,
138 using Velocity =
typename VelocityVector::value_type;
140 const auto geometry = element.geometry();
141 const Dune::GeometryType geomType = geometry.type();
144 auto upwindTerm = [phaseIdx](
const auto& volVars) {
return volVars.mobility(phaseIdx); };
146 if constexpr (isBox && dim == 1)
148 Velocity tmpVelocity(0.0);
149 tmpVelocity = (geometry.corner(1) - geometry.corner(0));
150 tmpVelocity /= tmpVelocity.two_norm();
152 for (
auto&& scvf : scvfs(fvGeometry))
158 FluxVariables fluxVars;
159 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
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();
169 const int eIdxGlobal = gridGeometry_.elementMapper().index(element);
170 velocity[eIdxGlobal] = tmpVelocity;
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);
183 using ScvVelocities = Dune::BlockVector<Velocity>;
184 ScvVelocities scvVelocities(fvGeometry.numScv());
187 for (
auto&& scvf : scvfs(fvGeometry))
193 const auto localPosIP = geometry.local(scvf.ipGlobal());
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();
203 FluxVariables fluxVars;
204 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
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();
214 Velocity tmpVelocity = localNormal;
217 scvVelocities[scvf.insideScvIdx()] += tmpVelocity;
218 scvVelocities[scvf.outsideScvIdx()] += tmpVelocity;
222 for (
auto&& scv : scvs(fvGeometry))
224 int vIdxGlobal = scv.dofIndex();
227 Velocity scvVelocity(0);
229 jacobianT2.mtv(scvVelocities[scv.indexInElement()], scvVelocity);
230 scvVelocity /= geometry.integrationElement(localPos)*cellNum_[vIdxGlobal];
232 velocity[vIdxGlobal] += scvVelocity;
241 const int numScvfsPerFace = isMpfa ? element.template subEntity<1>(0).geometry().corners() : 1;
243 if (fvGeometry.numScvf() != element.subEntities(1)*numScvfsPerFace)
244 DUNE_THROW(Dune::NotImplemented,
"Velocity output for non-conforming grids");
246 if (!geomType.isCube() && !geomType.isSimplex())
247 DUNE_THROW(Dune::NotImplemented,
"Velocity output for other geometry types than cube and simplex");
253 std::vector<bool> handledScvf;
255 handledScvf.resize(element.subEntities(1),
false);
258 std::vector<unsigned int> scvfIndexInInside(fvGeometry.numScvf());
259 int localScvfIdx = 0;
260 for (
const auto& intersection : intersections(gridGeometry_.gridView(), element))
263 if (handledScvf[intersection.indexInInside()])
266 if (intersection.neighbor() || intersection.boundary())
268 for (
int i = 0; i < numScvfsPerFace; ++i)
269 scvfIndexInInside[localScvfIdx++] = intersection.indexInInside();
273 handledScvf[intersection.indexInInside()] =
true;
277 std::vector<Scalar> scvfFluxes(element.subEntities(1), 0.0);
279 for (
auto&& scvf : scvfs(fvGeometry))
281 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
282 if (!scvf.boundary())
284 FluxVariables fluxVars;
285 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
286 scvfFluxes[scvfIndexInInside[localScvfIdx]] += fluxVars.advectiveFlux(phaseIdx, upwindTerm)/insideVolVars.extrusionFactor();
290 auto bcTypes = problem_.boundaryTypes(element, scvf);
291 if (bcTypes.hasOnlyDirichlet())
293 FluxVariables fluxVars;
294 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
295 scvfFluxes[scvfIndexInInside[localScvfIdx]] += fluxVars.advectiveFlux(phaseIdx, upwindTerm)/insideVolVars.extrusionFactor();
310 for (
auto&& scvf : scvfs(fvGeometry))
314 auto bcTypes = problem_.boundaryTypes(element, scvf);
315 if (bcTypes.hasNeumann())
321 if (stationaryVelocityField)
323 const auto flux = AdvectionType::flux(problem_, element, fvGeometry, elemVolVars,
324 scvf, phaseIdx, elemFluxVarsCache);
326 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
327 scvfFluxes[scvfIndexInInside[localScvfIdx]] += flux / insideVolVars.extrusionFactor();
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;
339 else if (dim == 1 || geomType.isCube())
341 const auto fIdx = scvfIndexInInside[localScvfIdx];
343 if constexpr (!modelIsCompositional)
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();
354 const auto fIdxOpposite = fIdx%2 ? fIdx-1 : fIdx+1;
355 scvfFluxes[fIdx] = -scvfFluxes[fIdxOpposite];
360 else if (geomType.isSimplex())
361 scvfFluxes[scvfIndexInInside[localScvfIdx]] = 0;
364 DUNE_THROW(Dune::NotImplemented,
"Velocity computation at Neumann boundaries for cell-centered and prism/pyramid");
374 Velocity refVelocity;
377 if (dim == 1 || geomType.isCube())
379 for (
int i = 0; i < dim; i++)
380 refVelocity[i] = 0.5 * (scvfFluxes[2*i + 1] - scvfFluxes[2*i]);
383 else if (geomType.isSimplex())
385 for (
int dimIdx = 0; dimIdx < dim; dimIdx++)
387 refVelocity[dimIdx] = -scvfFluxes[dim - 1 - dimIdx];
388 for (
int fIdx = 0; fIdx < dim + 1; fIdx++)
389 refVelocity[dimIdx] += scvfFluxes[fIdx]/(dim + 1);
394 DUNE_THROW(Dune::NotImplemented,
"Velocity computation for cell-centered and prism/pyramid");
396 Velocity scvVelocity(0);
397 jacobianT2.mtv(refVelocity, scvVelocity);
399 scvVelocity /= geometry.integrationElement(localPos);
401 int eIdxGlobal = gridGeometry_.elementMapper().index(element);
403 velocity[eIdxGlobal] = scvVelocity;
412 static Scalar scvfReferenceArea_(Dune::GeometryType geomType,
int fIdx)
414 if (dim == 1 || geomType.isCube())
416 return 1.0/(1 << (dim-1));
418 else if (geomType.isTriangle())
420 static const Scalar faceToArea[] = {0.372677996249965,
423 return faceToArea[fIdx];
425 else if (geomType.isTetrahedron())
427 static const Scalar faceToArea[] = {0.102062072615966,
433 return faceToArea[fIdx];
435 else if (geomType.isPyramid())
437 static const Scalar faceToArea[] = {0.130437298687488,
445 return faceToArea[fIdx];
447 else if (geomType.isPrism())
449 static const Scalar faceToArea[] = {0.166666666666667,
458 return faceToArea[fIdx];
461 DUNE_THROW(Dune::NotImplemented,
"scvf area for unknown GeometryType");
466 const Problem& problem_;
467 const GridGeometry& gridGeometry_;
468 const GridVariables& gridVariables_;
470 std::vector<int> cellNum_;
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
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