24#ifndef DUMUX_DISCRETIZATION_CC_TPFA_DARCYS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_TPFA_DARCYS_LAW_HH
37template<
class TypeTag, DiscretizationMethod discMethod>
38class DarcysLawImplementation;
48template<
class Scalar,
class Gr
idGeometry,
bool isNetwork>
56template <
class TypeTag>
59 GetPropType<TypeTag, Properties::GridGeometry>,
60 (GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimension < GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimensionworld)>
67template<class GridGeometry>
68class TpfaDarcysLawCacheFiller
70 using FVElementGeometry = typename GridGeometry::LocalView;
71 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
72 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
78 template<class FluxVariablesCache, class Problem, class ElementVolumeVariables, class FluxVariablesCacheFiller>
79 static void fill(FluxVariablesCache& scvfFluxVarsCache,
80 const Problem& problem,
81 const Element& element,
82 const FVElementGeometry& fvGeometry,
83 const ElementVolumeVariables& elemVolVars,
84 const SubControlVolumeFace& scvf,
85 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
87 scvfFluxVarsCache.updateAdvection(problem, element, fvGeometry, elemVolVars, scvf);
95template<class AdvectionType, class GridGeometry>
96class TpfaDarcysLawCache
98 using Scalar = typename AdvectionType::Scalar;
99 using FVElementGeometry = typename GridGeometry::LocalView;
100 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
101 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
104 using Filler = TpfaDarcysLawCacheFiller<GridGeometry>;
106 template<class Problem, class ElementVolumeVariables>
107 void updateAdvection(const Problem& problem,
108 const Element& element,
109 const FVElementGeometry& fvGeometry,
110 const ElementVolumeVariables& elemVolVars,
111 const SubControlVolumeFace &scvf)
113 tij_ = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
116 const Scalar& advectionTij() const
127template<class ScalarType, class GridGeometry>
128class CCTpfaDarcysLaw<ScalarType, GridGeometry, false>
130 using ThisType = CCTpfaDarcysLaw<ScalarType, GridGeometry, false>;
131 using FVElementGeometry = typename GridGeometry::LocalView;
132 using SubControlVolume = typename GridGeometry::SubControlVolume;
133 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
134 using GridView = typename GridGeometry::GridView;
135 using Element = typename GridView::template Codim<0>::Entity;
137 static constexpr int dim = GridView::dimension;
138 static constexpr int dimWorld = GridView::dimensionworld;
140 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
144 using Scalar = ScalarType;
147 static const DiscretizationMethod discMethod = DiscretizationMethod::cctpfa;
150 using Cache = TpfaDarcysLawCache<ThisType, GridGeometry>;
153 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
154 static Scalar flux(const Problem& problem,
155 const Element& element,
156 const FVElementGeometry& fvGeometry,
157 const ElementVolumeVariables& elemVolVars,
158 const SubControlVolumeFace& scvf,
160 const ElementFluxVarsCache& elemFluxVarsCache)
162 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
164 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
167 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
168 const auto& insideVolVars = elemVolVars[insideScv];
169 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
174 const auto rho = scvf.boundary() ? outsideVolVars.density(phaseIdx)
175 : (insideVolVars.density(phaseIdx) + outsideVolVars.density(phaseIdx))*0.5;
178 const auto pInside = insideVolVars.pressure(phaseIdx);
179 const auto pOutside = outsideVolVars.pressure(phaseIdx);
181 const auto& tij = fluxVarsCache.advectionTij();
182 const auto& g = problem.spatialParams().gravity(scvf.ipGlobal());
185 const auto alpha_inside = vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g)*insideVolVars.extrusionFactor();
187 Scalar flux = tij*(pInside - pOutside) + rho*scvf.area()*alpha_inside;
190 if (!scvf.boundary())
192 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
193 const auto outsideK = outsideVolVars.permeability();
194 const auto outsideTi = fvGeometry.gridGeometry().isPeriodic()
195 ? computeTpfaTransmissibility(fvGeometry.flipScvf(scvf.index()), outsideScv, outsideK, outsideVolVars.extrusionFactor())
196 : -1.0*computeTpfaTransmissibility(scvf, outsideScv, outsideK, outsideVolVars.extrusionFactor());
197 const auto alpha_outside = vtmv(scvf.unitOuterNormal(), outsideK, g)*outsideVolVars.extrusionFactor();
199 flux -= rho*tij/outsideTi*(alpha_inside - alpha_outside);
207 const auto pInside = insideVolVars.pressure(phaseIdx);
208 const auto pOutside = outsideVolVars.pressure(phaseIdx);
211 return fluxVarsCache.advectionTij()*(pInside - pOutside);
217 template<class Problem, class ElementVolumeVariables>
218 static Scalar calculateTransmissibility(const Problem& problem,
219 const Element& element,
220 const FVElementGeometry& fvGeometry,
221 const ElementVolumeVariables& elemVolVars,
222 const SubControlVolumeFace& scvf)
226 const auto insideScvIdx = scvf.insideScvIdx();
227 const auto& insideScv = fvGeometry.scv(insideScvIdx);
228 const auto& insideVolVars = elemVolVars[insideScvIdx];
230 const Scalar ti = computeTpfaTransmissibility(scvf, insideScv,
231 getPermeability_(problem, insideVolVars, scvf.ipGlobal()),
232 insideVolVars.extrusionFactor());
236 tij = scvf.area()*ti;
241 const auto outsideScvIdx = scvf.outsideScvIdx();
244 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
245 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
246 const Scalar tj = fvGeometry.gridGeometry().isPeriodic()
247 ? computeTpfaTransmissibility(fvGeometry.flipScvf(scvf.index()), outsideScv, getPermeability_(problem, outsideVolVars, scvf.ipGlobal()), outsideVolVars.extrusionFactor())
248 : -1.0*computeTpfaTransmissibility(scvf, outsideScv, getPermeability_(problem, outsideVolVars, scvf.ipGlobal()), outsideVolVars.extrusionFactor());
255 tij = scvf.area()*(ti * tj)/(ti + tj);
262 template<class Problem, class VolumeVariables,
263 std::enable_if_t<!Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
264 static decltype(auto) getPermeability_(const Problem& problem,
265 const VolumeVariables& volVars,
266 const GlobalPosition& scvfIpGlobal)
267 { return volVars.permeability(); }
269 template<class Problem, class VolumeVariables,
270 std::enable_if_t<Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
271 static decltype(auto) getPermeability_(const Problem& problem,
272 const VolumeVariables& volVars,
273 const GlobalPosition& scvfIpGlobal)
274 { return problem.spatialParams().permeabilityAtPos(scvfIpGlobal); }
281template<class ScalarType, class GridGeometry>
282class CCTpfaDarcysLaw<ScalarType, GridGeometry, true>
284 using ThisType = CCTpfaDarcysLaw<ScalarType, GridGeometry, true>;
285 using FVElementGeometry = typename GridGeometry::LocalView;
286 using SubControlVolume = typename GridGeometry::SubControlVolume;
287 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
288 using GridView = typename GridGeometry::GridView;
289 using Element = typename GridView::template Codim<0>::Entity;
291 static constexpr int dim = GridView::dimension;
292 static constexpr int dimWorld = GridView::dimensionworld;
294 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
298 using Scalar = ScalarType;
301 static const DiscretizationMethod discMethod = DiscretizationMethod::cctpfa;
304 using Cache = TpfaDarcysLawCache<ThisType, GridGeometry>;
307 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
308 static Scalar flux(const Problem& problem,
309 const Element& element,
310 const FVElementGeometry& fvGeometry,
311 const ElementVolumeVariables& elemVolVars,
312 const SubControlVolumeFace& scvf,
314 const ElementFluxVarsCache& elemFluxVarsCache)
316 static const bool gravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
318 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
321 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
322 const auto& insideVolVars = elemVolVars[insideScv];
323 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
328 const auto rho = [&]()
332 return outsideVolVars.density(phaseIdx);
335 else if (scvf.numOutsideScvs() == 1)
336 return (insideVolVars.density(phaseIdx) + outsideVolVars.density(phaseIdx))*0.5;
341 Scalar rho(insideVolVars.density(phaseIdx));
342 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
344 const auto outsideScvIdx = scvf.outsideScvIdx(i);
345 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
346 rho += outsideVolVars.density(phaseIdx);
348 return rho/(scvf.numOutsideScvs()+1);
352 const auto& tij = fluxVarsCache.advectionTij();
353 const auto& g = problem.spatialParams().gravity(scvf.ipGlobal());
356 const auto pInside = insideVolVars.pressure(phaseIdx);
357 const auto pOutside = [&]()
360 if (scvf.numOutsideScvs() == 1)
361 return outsideVolVars.pressure(phaseIdx);
367 Scalar sumPTi(tij*pInside);
370 sumPTi += rho*scvf.area()
371 *insideVolVars.extrusionFactor()
372 *vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g);
374 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
376 const auto outsideScvIdx = scvf.outsideScvIdx(i);
377 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
378 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
379 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
380 sumTi += outsideFluxVarsCache.advectionTij();
381 sumPTi += outsideFluxVarsCache.advectionTij()*outsideVolVars.pressure(phaseIdx);
384 sumPTi += rho*scvf.area()
385 *outsideVolVars.extrusionFactor()
386 *vtmv(flippedScvf.unitOuterNormal(), outsideVolVars.permeability(), g);
393 const auto alpha_inside = vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g)*insideVolVars.extrusionFactor();
395 Scalar flux = tij*(pInside - pOutside) + scvf.area()*rho*alpha_inside;
398 if (!scvf.boundary() && scvf.numOutsideScvs() == 1)
400 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
401 const auto& outsideScvf = fvGeometry.flipScvf(scvf.index());
402 const auto outsideK = outsideVolVars.permeability();
403 const auto outsideTi = computeTpfaTransmissibility(outsideScvf, outsideScv, outsideK, outsideVolVars.extrusionFactor());
404 const auto alpha_outside = vtmv(outsideScvf.unitOuterNormal(), outsideK, g)*outsideVolVars.extrusionFactor();
406 flux -= rho*tij/outsideTi*(alpha_inside + alpha_outside);
414 const auto pInside = insideVolVars.pressure(phaseIdx);
415 const auto pOutside = [&]()
418 if (scvf.numOutsideScvs() <= 1)
419 return outsideVolVars.pressure(phaseIdx);
424 const auto& insideFluxVarsCache = elemFluxVarsCache[scvf];
425 Scalar sumTi(insideFluxVarsCache.advectionTij());
426 Scalar sumPTi(insideFluxVarsCache.advectionTij()*pInside);
428 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
430 const auto outsideScvIdx = scvf.outsideScvIdx(i);
431 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
432 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
433 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
434 sumTi += outsideFluxVarsCache.advectionTij();
435 sumPTi += outsideFluxVarsCache.advectionTij()*outsideVolVars.pressure(phaseIdx);
442 return fluxVarsCache.advectionTij()*(pInside - pOutside);
448 template<class Problem, class ElementVolumeVariables>
449 static Scalar calculateTransmissibility(const Problem& problem,
450 const Element& element,
451 const FVElementGeometry& fvGeometry,
452 const ElementVolumeVariables& elemVolVars,
453 const SubControlVolumeFace& scvf)
457 const auto insideScvIdx = scvf.insideScvIdx();
458 const auto& insideScv = fvGeometry.scv(insideScvIdx);
459 const auto& insideVolVars = elemVolVars[insideScvIdx];
461 const Scalar ti = computeTpfaTransmissibility(scvf, insideScv,
462 getPermeability_(problem, insideVolVars, scvf.ipGlobal()),
463 insideVolVars.extrusionFactor());
466 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
467 tij = scvf.area()*ti;
472 const auto outsideScvIdx = scvf.outsideScvIdx();
475 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
476 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
477 const Scalar tj = computeTpfaTransmissibility(fvGeometry.flipScvf(scvf.index()), outsideScv,
478 getPermeability_(problem, outsideVolVars, scvf.ipGlobal()),
479 outsideVolVars.extrusionFactor());
486 tij = scvf.area()*(ti * tj)/(ti + tj);
493 template<class Problem, class VolumeVariables,
494 std::enable_if_t<!Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
495 static decltype(auto) getPermeability_(const Problem& problem,
496 const VolumeVariables& volVars,
497 const GlobalPosition& scvfIpGlobal)
498 { return volVars.permeability(); }
500 template<class Problem, class VolumeVariables,
501 std::enable_if_t<Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
502 static decltype(auto) getPermeability_(const Problem& problem,
503 const VolumeVariables& volVars,
504 const GlobalPosition& scvfIpGlobal)
505 { return problem.spatialParams().permeabilityAtPos(scvfIpGlobal); }
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
forward declaration of the method-specific implementation
Definition: flux/darcyslaw.hh:38
Darcy's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: flux/cctpfa/darcyslaw.hh:49
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...