24#ifndef DUMUX_DISCRETIZATION_CC_TPFA_DARCYS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_TPFA_DARCYS_LAW_HH
38template<
class TypeTag, DiscretizationMethod discMethod>
39class DarcysLawImplementation;
49template<
class Scalar,
class Gr
idGeometry,
bool isNetwork>
57template <
class TypeTag>
60 GetPropType<TypeTag, Properties::GridGeometry>,
61 (GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimension < GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimensionworld)>
68template<class GridGeometry>
69class TpfaDarcysLawCacheFiller
71 using FVElementGeometry = typename GridGeometry::LocalView;
72 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
73 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
79 template<class FluxVariablesCache, class Problem, class ElementVolumeVariables, class FluxVariablesCacheFiller>
80 static void fill(FluxVariablesCache& scvfFluxVarsCache,
81 const Problem& problem,
82 const Element& element,
83 const FVElementGeometry& fvGeometry,
84 const ElementVolumeVariables& elemVolVars,
85 const SubControlVolumeFace& scvf,
86 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
88 scvfFluxVarsCache.updateAdvection(problem, element, fvGeometry, elemVolVars, scvf);
96template<class AdvectionType, class GridGeometry>
97class TpfaDarcysLawCache
99 using Scalar = typename AdvectionType::Scalar;
100 using FVElementGeometry = typename GridGeometry::LocalView;
101 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
102 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
105 using Filler = TpfaDarcysLawCacheFiller<GridGeometry>;
107 template<class Problem, class ElementVolumeVariables>
108 void updateAdvection(const Problem& problem,
109 const Element& element,
110 const FVElementGeometry& fvGeometry,
111 const ElementVolumeVariables& elemVolVars,
112 const SubControlVolumeFace &scvf)
114 tij_ = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
117 const Scalar& advectionTij() const
128template<class ScalarType, class GridGeometry>
129class CCTpfaDarcysLaw<ScalarType, GridGeometry, false>
131 using ThisType = CCTpfaDarcysLaw<ScalarType, GridGeometry, false>;
132 using FVElementGeometry = typename GridGeometry::LocalView;
133 using SubControlVolume = typename GridGeometry::SubControlVolume;
134 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
135 using Extrusion = Extrusion_t<GridGeometry>;
136 using GridView = typename GridGeometry::GridView;
137 using Element = typename GridView::template Codim<0>::Entity;
139 static constexpr int dim = GridView::dimension;
140 static constexpr int dimWorld = GridView::dimensionworld;
142 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
146 using Scalar = ScalarType;
149 static const DiscretizationMethod discMethod = DiscretizationMethod::cctpfa;
152 using Cache = TpfaDarcysLawCache<ThisType, GridGeometry>;
164 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
165 static Scalar flux(const Problem& problem,
166 const Element& element,
167 const FVElementGeometry& fvGeometry,
168 const ElementVolumeVariables& elemVolVars,
169 const SubControlVolumeFace& scvf,
171 const ElementFluxVarsCache& elemFluxVarsCache)
173 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
175 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
178 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
179 const auto& insideVolVars = elemVolVars[insideScv];
180 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
185 const auto rho = scvf.boundary() ? outsideVolVars.density(phaseIdx)
186 : (insideVolVars.density(phaseIdx) + outsideVolVars.density(phaseIdx))*0.5;
189 const auto pInside = insideVolVars.pressure(phaseIdx);
190 const auto pOutside = outsideVolVars.pressure(phaseIdx);
192 const auto& tij = fluxVarsCache.advectionTij();
193 const auto& g = problem.spatialParams().gravity(scvf.ipGlobal());
196 const auto alpha_inside = vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g)*insideVolVars.extrusionFactor();
198 Scalar flux = tij*(pInside - pOutside) + rho*Extrusion::area(scvf)*alpha_inside;
201 if (!scvf.boundary())
203 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
204 const auto outsideK = outsideVolVars.permeability();
205 const auto outsideTi = fvGeometry.gridGeometry().isPeriodic()
206 ? computeTpfaTransmissibility(fvGeometry.flipScvf(scvf.index()), outsideScv, outsideK, outsideVolVars.extrusionFactor())
207 : -1.0*computeTpfaTransmissibility(scvf, outsideScv, outsideK, outsideVolVars.extrusionFactor());
208 const auto alpha_outside = vtmv(scvf.unitOuterNormal(), outsideK, g)*outsideVolVars.extrusionFactor();
210 flux -= rho*tij/outsideTi*(alpha_inside - alpha_outside);
218 const auto pInside = insideVolVars.pressure(phaseIdx);
219 const auto pOutside = outsideVolVars.pressure(phaseIdx);
222 return fluxVarsCache.advectionTij()*(pInside - pOutside);
228 template<class Problem, class ElementVolumeVariables>
229 static Scalar calculateTransmissibility(const Problem& problem,
230 const Element& element,
231 const FVElementGeometry& fvGeometry,
232 const ElementVolumeVariables& elemVolVars,
233 const SubControlVolumeFace& scvf)
237 const auto insideScvIdx = scvf.insideScvIdx();
238 const auto& insideScv = fvGeometry.scv(insideScvIdx);
239 const auto& insideVolVars = elemVolVars[insideScvIdx];
241 const Scalar ti = computeTpfaTransmissibility(scvf, insideScv,
242 getPermeability_(problem, insideVolVars, scvf.ipGlobal()),
243 insideVolVars.extrusionFactor());
247 tij = Extrusion::area(scvf)*ti;
252 const auto outsideScvIdx = scvf.outsideScvIdx();
255 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
256 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
257 const Scalar tj = fvGeometry.gridGeometry().isPeriodic()
258 ? computeTpfaTransmissibility(fvGeometry.flipScvf(scvf.index()), outsideScv, getPermeability_(problem, outsideVolVars, scvf.ipGlobal()), outsideVolVars.extrusionFactor())
259 : -1.0*computeTpfaTransmissibility(scvf, outsideScv, getPermeability_(problem, outsideVolVars, scvf.ipGlobal()), outsideVolVars.extrusionFactor());
266 tij = Extrusion::area(scvf)*(ti * tj)/(ti + tj);
273 template<class Problem, class VolumeVariables,
274 std::enable_if_t<!Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
275 static decltype(auto) getPermeability_(const Problem& problem,
276 const VolumeVariables& volVars,
277 const GlobalPosition& scvfIpGlobal)
278 { return volVars.permeability(); }
280 template<class Problem, class VolumeVariables,
281 std::enable_if_t<Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
282 static decltype(auto) getPermeability_(const Problem& problem,
283 const VolumeVariables& volVars,
284 const GlobalPosition& scvfIpGlobal)
285 { return problem.spatialParams().permeabilityAtPos(scvfIpGlobal); }
292template<class ScalarType, class GridGeometry>
293class CCTpfaDarcysLaw<ScalarType, GridGeometry, true>
295 using ThisType = CCTpfaDarcysLaw<ScalarType, GridGeometry, true>;
296 using FVElementGeometry = typename GridGeometry::LocalView;
297 using SubControlVolume = typename GridGeometry::SubControlVolume;
298 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
299 using Extrusion = Extrusion_t<GridGeometry>;
300 using GridView = typename GridGeometry::GridView;
301 using Element = typename GridView::template Codim<0>::Entity;
303 static constexpr int dim = GridView::dimension;
304 static constexpr int dimWorld = GridView::dimensionworld;
306 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
310 using Scalar = ScalarType;
313 static const DiscretizationMethod discMethod = DiscretizationMethod::cctpfa;
316 using Cache = TpfaDarcysLawCache<ThisType, GridGeometry>;
328 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
329 static Scalar flux(const Problem& problem,
330 const Element& element,
331 const FVElementGeometry& fvGeometry,
332 const ElementVolumeVariables& elemVolVars,
333 const SubControlVolumeFace& scvf,
335 const ElementFluxVarsCache& elemFluxVarsCache)
337 static const bool gravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
339 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
342 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
343 const auto& insideVolVars = elemVolVars[insideScv];
344 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
349 const auto rho = [&]()
353 return outsideVolVars.density(phaseIdx);
356 else if (scvf.numOutsideScvs() == 1)
357 return (insideVolVars.density(phaseIdx) + outsideVolVars.density(phaseIdx))*0.5;
362 Scalar rho(insideVolVars.density(phaseIdx));
363 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
365 const auto outsideScvIdx = scvf.outsideScvIdx(i);
366 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
367 rho += outsideVolVars.density(phaseIdx);
369 return rho/(scvf.numOutsideScvs()+1);
373 const auto& tij = fluxVarsCache.advectionTij();
374 const auto& g = problem.spatialParams().gravity(scvf.ipGlobal());
377 const auto pInside = insideVolVars.pressure(phaseIdx);
378 const auto pOutside = [&]()
381 if (scvf.numOutsideScvs() == 1)
382 return outsideVolVars.pressure(phaseIdx);
388 Scalar sumPTi(tij*pInside);
391 sumPTi += rho*Extrusion::area(scvf)
392 *insideVolVars.extrusionFactor()
393 *vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g);
395 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
397 const auto outsideScvIdx = scvf.outsideScvIdx(i);
398 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
399 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
400 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
401 sumTi += outsideFluxVarsCache.advectionTij();
402 sumPTi += outsideFluxVarsCache.advectionTij()*outsideVolVars.pressure(phaseIdx);
405 sumPTi += rho*Extrusion::area(scvf)
406 *outsideVolVars.extrusionFactor()
407 *vtmv(flippedScvf.unitOuterNormal(), outsideVolVars.permeability(), g);
414 const auto alpha_inside = vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g)*insideVolVars.extrusionFactor();
416 Scalar flux = tij*(pInside - pOutside) + Extrusion::area(scvf)*rho*alpha_inside;
419 if (!scvf.boundary() && scvf.numOutsideScvs() == 1)
421 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
422 const auto& outsideScvf = fvGeometry.flipScvf(scvf.index());
423 const auto outsideK = outsideVolVars.permeability();
424 const auto outsideTi = computeTpfaTransmissibility(outsideScvf, outsideScv, outsideK, outsideVolVars.extrusionFactor());
425 const auto alpha_outside = vtmv(outsideScvf.unitOuterNormal(), outsideK, g)*outsideVolVars.extrusionFactor();
427 flux -= rho*tij/outsideTi*(alpha_inside + alpha_outside);
435 const auto pInside = insideVolVars.pressure(phaseIdx);
436 const auto pOutside = [&]()
439 if (scvf.numOutsideScvs() <= 1)
440 return outsideVolVars.pressure(phaseIdx);
445 const auto& insideFluxVarsCache = elemFluxVarsCache[scvf];
446 Scalar sumTi(insideFluxVarsCache.advectionTij());
447 Scalar sumPTi(insideFluxVarsCache.advectionTij()*pInside);
449 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
451 const auto outsideScvIdx = scvf.outsideScvIdx(i);
452 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
453 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
454 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
455 sumTi += outsideFluxVarsCache.advectionTij();
456 sumPTi += outsideFluxVarsCache.advectionTij()*outsideVolVars.pressure(phaseIdx);
463 return fluxVarsCache.advectionTij()*(pInside - pOutside);
469 template<class Problem, class ElementVolumeVariables>
470 static Scalar calculateTransmissibility(const Problem& problem,
471 const Element& element,
472 const FVElementGeometry& fvGeometry,
473 const ElementVolumeVariables& elemVolVars,
474 const SubControlVolumeFace& scvf)
478 const auto insideScvIdx = scvf.insideScvIdx();
479 const auto& insideScv = fvGeometry.scv(insideScvIdx);
480 const auto& insideVolVars = elemVolVars[insideScvIdx];
482 const Scalar ti = computeTpfaTransmissibility(scvf, insideScv,
483 getPermeability_(problem, insideVolVars, scvf.ipGlobal()),
484 insideVolVars.extrusionFactor());
487 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
488 tij = Extrusion::area(scvf)*ti;
493 const auto outsideScvIdx = scvf.outsideScvIdx();
496 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
497 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
498 const Scalar tj = computeTpfaTransmissibility(fvGeometry.flipScvf(scvf.index()), outsideScv,
499 getPermeability_(problem, outsideVolVars, scvf.ipGlobal()),
500 outsideVolVars.extrusionFactor());
507 tij = Extrusion::area(scvf)*(ti * tj)/(ti + tj);
514 template<class Problem, class VolumeVariables,
515 std::enable_if_t<!Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
516 static decltype(auto) getPermeability_(const Problem& problem,
517 const VolumeVariables& volVars,
518 const GlobalPosition& scvfIpGlobal)
519 { return volVars.permeability(); }
521 template<class Problem, class VolumeVariables,
522 std::enable_if_t<Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
523 static decltype(auto) getPermeability_(const Problem& problem,
524 const VolumeVariables& volVars,
525 const GlobalPosition& scvfIpGlobal)
526 { 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.
Helper classes to compute the integration elements.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
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:50
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...