24#ifndef DUMUX_DISCRETIZATION_CC_TPFA_DARCYS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_TPFA_DARCYS_LAW_HH
38template<
class TypeTag,
class DiscretizationMethod>
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;
148 using DiscretizationMethod = DiscretizationMethods::CCTpfa;
150 static constexpr DiscretizationMethod discMethod{};
153 using Cache = TpfaDarcysLawCache<ThisType, GridGeometry>;
165 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
166 static Scalar flux(const Problem& problem,
167 const Element& element,
168 const FVElementGeometry& fvGeometry,
169 const ElementVolumeVariables& elemVolVars,
170 const SubControlVolumeFace& scvf,
172 const ElementFluxVarsCache& elemFluxVarsCache)
174 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
176 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
179 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
180 const auto& insideVolVars = elemVolVars[insideScv];
181 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
186 const auto rho = scvf.boundary() ? outsideVolVars.density(phaseIdx)
187 : (insideVolVars.density(phaseIdx) + outsideVolVars.density(phaseIdx))*0.5;
190 const auto pInside = insideVolVars.pressure(phaseIdx);
191 const auto pOutside = outsideVolVars.pressure(phaseIdx);
193 const auto& tij = fluxVarsCache.advectionTij();
194 const auto& g = problem.spatialParams().gravity(scvf.ipGlobal());
197 const auto alpha_inside = vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g)*insideVolVars.extrusionFactor();
199 Scalar flux = tij*(pInside - pOutside) + rho*Extrusion::area(fvGeometry, scvf)*alpha_inside;
202 if (!scvf.boundary())
204 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
205 const auto outsideK = outsideVolVars.permeability();
206 const auto outsideTi = fvGeometry.gridGeometry().isPeriodic()
207 ? computeTpfaTransmissibility(fvGeometry, fvGeometry.flipScvf(scvf.index()), outsideScv, outsideK, outsideVolVars.extrusionFactor())
208 : -1.0*computeTpfaTransmissibility(fvGeometry, scvf, outsideScv, outsideK, outsideVolVars.extrusionFactor());
209 const auto alpha_outside = vtmv(scvf.unitOuterNormal(), outsideK, g)*outsideVolVars.extrusionFactor();
211 flux -= rho*tij/outsideTi*(alpha_inside - alpha_outside);
219 const auto pInside = insideVolVars.pressure(phaseIdx);
220 const auto pOutside = outsideVolVars.pressure(phaseIdx);
223 return fluxVarsCache.advectionTij()*(pInside - pOutside);
229 template<class Problem, class ElementVolumeVariables>
230 static Scalar calculateTransmissibility(const Problem& problem,
231 const Element& element,
232 const FVElementGeometry& fvGeometry,
233 const ElementVolumeVariables& elemVolVars,
234 const SubControlVolumeFace& scvf)
238 const auto insideScvIdx = scvf.insideScvIdx();
239 const auto& insideScv = fvGeometry.scv(insideScvIdx);
240 const auto& insideVolVars = elemVolVars[insideScvIdx];
242 const Scalar ti = computeTpfaTransmissibility(fvGeometry, scvf, insideScv,
243 getPermeability_(problem, insideVolVars, scvf.ipGlobal()),
244 insideVolVars.extrusionFactor());
248 tij = Extrusion::area(fvGeometry, scvf)*ti;
253 const auto outsideScvIdx = scvf.outsideScvIdx();
256 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
257 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
258 const Scalar tj = fvGeometry.gridGeometry().isPeriodic()
259 ? computeTpfaTransmissibility(fvGeometry, fvGeometry.flipScvf(scvf.index()), outsideScv, getPermeability_(problem, outsideVolVars, scvf.ipGlobal()), outsideVolVars.extrusionFactor())
260 : -1.0*computeTpfaTransmissibility(fvGeometry, scvf, outsideScv, getPermeability_(problem, outsideVolVars, scvf.ipGlobal()), outsideVolVars.extrusionFactor());
267 tij = Extrusion::area(fvGeometry, scvf)*(ti * tj)/(ti + tj);
274 template<class Problem, class VolumeVariables,
275 std::enable_if_t<!Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
276 static decltype(auto) getPermeability_(const Problem& problem,
277 const VolumeVariables& volVars,
278 const GlobalPosition& scvfIpGlobal)
279 { return volVars.permeability(); }
281 template<class Problem, class VolumeVariables,
282 std::enable_if_t<Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
283 static decltype(auto) getPermeability_(const Problem& problem,
284 const VolumeVariables& volVars,
285 const GlobalPosition& scvfIpGlobal)
286 { return problem.spatialParams().permeabilityAtPos(scvfIpGlobal); }
293template<class ScalarType, class GridGeometry>
294class CCTpfaDarcysLaw<ScalarType, GridGeometry, true>
296 using ThisType = CCTpfaDarcysLaw<ScalarType, GridGeometry, true>;
297 using FVElementGeometry = typename GridGeometry::LocalView;
298 using SubControlVolume = typename GridGeometry::SubControlVolume;
299 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
300 using Extrusion = Extrusion_t<GridGeometry>;
301 using GridView = typename GridGeometry::GridView;
302 using Element = typename GridView::template Codim<0>::Entity;
304 static constexpr int dim = GridView::dimension;
305 static constexpr int dimWorld = GridView::dimensionworld;
307 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
311 using Scalar = ScalarType;
313 using DiscretizationMethod = DiscretizationMethods::CCTpfa;
315 static constexpr DiscretizationMethod discMethod{};
318 using Cache = TpfaDarcysLawCache<ThisType, GridGeometry>;
330 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
331 static Scalar flux(const Problem& problem,
332 const Element& element,
333 const FVElementGeometry& fvGeometry,
334 const ElementVolumeVariables& elemVolVars,
335 const SubControlVolumeFace& scvf,
337 const ElementFluxVarsCache& elemFluxVarsCache)
339 static const bool gravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
341 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
344 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
345 const auto& insideVolVars = elemVolVars[insideScv];
346 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
351 const auto rho = [&]()
355 return outsideVolVars.density(phaseIdx);
358 else if (scvf.numOutsideScvs() == 1)
359 return (insideVolVars.density(phaseIdx) + outsideVolVars.density(phaseIdx))*0.5;
364 Scalar rho(insideVolVars.density(phaseIdx));
365 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
367 const auto outsideScvIdx = scvf.outsideScvIdx(i);
368 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
369 rho += outsideVolVars.density(phaseIdx);
371 return rho/(scvf.numOutsideScvs()+1);
375 const auto& tij = fluxVarsCache.advectionTij();
376 const auto& g = problem.spatialParams().gravity(scvf.ipGlobal());
379 const auto pInside = insideVolVars.pressure(phaseIdx);
380 const auto pOutside = [&]()
383 if (scvf.numOutsideScvs() == 1)
384 return outsideVolVars.pressure(phaseIdx);
390 Scalar sumPTi(tij*pInside);
393 sumPTi += rho*Extrusion::area(fvGeometry, scvf)
394 *insideVolVars.extrusionFactor()
395 *vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g);
397 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
399 const auto outsideScvIdx = scvf.outsideScvIdx(i);
400 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
401 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
402 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
403 sumTi += outsideFluxVarsCache.advectionTij();
404 sumPTi += outsideFluxVarsCache.advectionTij()*outsideVolVars.pressure(phaseIdx);
407 sumPTi += rho*Extrusion::area(fvGeometry, scvf)
408 *outsideVolVars.extrusionFactor()
409 *vtmv(flippedScvf.unitOuterNormal(), outsideVolVars.permeability(), g);
416 const auto alpha_inside = vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g)*insideVolVars.extrusionFactor();
418 Scalar flux = tij*(pInside - pOutside) + Extrusion::area(fvGeometry, scvf)*rho*alpha_inside;
421 if (!scvf.boundary() && scvf.numOutsideScvs() == 1)
423 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
424 const auto& outsideScvf = fvGeometry.flipScvf(scvf.index());
425 const auto outsideK = outsideVolVars.permeability();
426 const auto outsideTi = computeTpfaTransmissibility(fvGeometry, outsideScvf, outsideScv, outsideK, outsideVolVars.extrusionFactor());
427 const auto alpha_outside = vtmv(outsideScvf.unitOuterNormal(), outsideK, g)*outsideVolVars.extrusionFactor();
429 flux -= rho*tij/outsideTi*(alpha_inside + alpha_outside);
437 const auto pInside = insideVolVars.pressure(phaseIdx);
438 const auto pOutside = [&]()
441 if (scvf.numOutsideScvs() <= 1)
442 return outsideVolVars.pressure(phaseIdx);
447 const auto& insideFluxVarsCache = elemFluxVarsCache[scvf];
448 Scalar sumTi(insideFluxVarsCache.advectionTij());
449 Scalar sumPTi(insideFluxVarsCache.advectionTij()*pInside);
451 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
453 const auto outsideScvIdx = scvf.outsideScvIdx(i);
454 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
455 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
456 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
457 sumTi += outsideFluxVarsCache.advectionTij();
458 sumPTi += outsideFluxVarsCache.advectionTij()*outsideVolVars.pressure(phaseIdx);
465 return fluxVarsCache.advectionTij()*(pInside - pOutside);
471 template<class Problem, class ElementVolumeVariables>
472 static Scalar calculateTransmissibility(const Problem& problem,
473 const Element& element,
474 const FVElementGeometry& fvGeometry,
475 const ElementVolumeVariables& elemVolVars,
476 const SubControlVolumeFace& scvf)
480 const auto insideScvIdx = scvf.insideScvIdx();
481 const auto& insideScv = fvGeometry.scv(insideScvIdx);
482 const auto& insideVolVars = elemVolVars[insideScvIdx];
484 const Scalar ti = computeTpfaTransmissibility(fvGeometry, scvf, insideScv,
485 getPermeability_(problem, insideVolVars, scvf.ipGlobal()),
486 insideVolVars.extrusionFactor());
489 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
490 tij = Extrusion::area(fvGeometry, scvf)*ti;
495 const auto outsideScvIdx = scvf.outsideScvIdx();
498 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
499 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
500 const Scalar tj = computeTpfaTransmissibility(fvGeometry, fvGeometry.flipScvf(scvf.index()), outsideScv,
501 getPermeability_(problem, outsideVolVars, scvf.ipGlobal()),
502 outsideVolVars.extrusionFactor());
509 tij = Extrusion::area(fvGeometry, scvf)*(ti * tj)/(ti + tj);
516 template<class Problem, class VolumeVariables,
517 std::enable_if_t<!Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
518 static decltype(auto) getPermeability_(const Problem& problem,
519 const VolumeVariables& volVars,
520 const GlobalPosition& scvfIpGlobal)
521 { return volVars.permeability(); }
523 template<class Problem, class VolumeVariables,
524 std::enable_if_t<Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
525 static decltype(auto) getPermeability_(const Problem& problem,
526 const VolumeVariables& volVars,
527 const GlobalPosition& scvfIpGlobal)
528 { return problem.spatialParams().permeabilityAtPos(scvfIpGlobal); }
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Define some often used mathematical functions.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
forward declaration of the method-specific implementation
Definition: flux/box/darcyslaw.hh:44
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...