12#ifndef DUMUX_DISCRETIZATION_CC_TPFA_DARCYS_LAW_HH
13#define DUMUX_DISCRETIZATION_CC_TPFA_DARCYS_LAW_HH
26template<
class TypeTag,
class DiscretizationMethod>
27class DarcysLawImplementation;
37template<
class Scalar,
class Gr
idGeometry,
bool isNetwork>
45template <
class TypeTag>
48 GetPropType<TypeTag, Properties::GridGeometry>,
49 (GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimension < GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimensionworld)>
56template<class GridGeometry>
57class TpfaDarcysLawCacheFiller
59 using FVElementGeometry = typename GridGeometry::LocalView;
60 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
61 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
67 template<class FluxVariablesCache, class Problem, class ElementVolumeVariables, class FluxVariablesCacheFiller>
68 static void fill(FluxVariablesCache& scvfFluxVarsCache,
69 const Problem& problem,
70 const Element& element,
71 const FVElementGeometry& fvGeometry,
72 const ElementVolumeVariables& elemVolVars,
73 const SubControlVolumeFace& scvf,
74 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
76 scvfFluxVarsCache.updateAdvection(problem, element, fvGeometry, elemVolVars, scvf);
84template<class AdvectionType, class GridGeometry>
85class TpfaDarcysLawCache
87 using Scalar = typename AdvectionType::Scalar;
88 using FVElementGeometry = typename GridGeometry::LocalView;
89 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
90 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
93 using Filler = TpfaDarcysLawCacheFiller<GridGeometry>;
95 template<class Problem, class ElementVolumeVariables>
96 void updateAdvection(const Problem& problem,
97 const Element& element,
98 const FVElementGeometry& fvGeometry,
99 const ElementVolumeVariables& elemVolVars,
100 const SubControlVolumeFace &scvf)
102 tij_ = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
105 const Scalar& advectionTij() const
116template<class ScalarType, class GridGeometry>
117class CCTpfaDarcysLaw<ScalarType, GridGeometry, false>
119 using ThisType = CCTpfaDarcysLaw<ScalarType, GridGeometry, false>;
120 using FVElementGeometry = typename GridGeometry::LocalView;
121 using SubControlVolume = typename GridGeometry::SubControlVolume;
122 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
123 using Extrusion = Extrusion_t<GridGeometry>;
124 using GridView = typename GridGeometry::GridView;
125 using Element = typename GridView::template Codim<0>::Entity;
127 static constexpr int dim = GridView::dimension;
128 static constexpr int dimWorld = GridView::dimensionworld;
130 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
134 using Scalar = ScalarType;
136 using DiscretizationMethod = DiscretizationMethods::CCTpfa;
138 static constexpr DiscretizationMethod discMethod{};
141 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*Extrusion::area(fvGeometry, scvf)*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, fvGeometry.flipScvf(scvf.index()), outsideScv, outsideK, outsideVolVars.extrusionFactor())
196 : -1.0*computeTpfaTransmissibility(fvGeometry, 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(fvGeometry, scvf, insideScv,
231 getPermeability_(problem, insideVolVars, scvf.ipGlobal()),
232 insideVolVars.extrusionFactor());
236 tij = Extrusion::area(fvGeometry, scvf)*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, fvGeometry.flipScvf(scvf.index()), outsideScv, getPermeability_(problem, outsideVolVars, scvf.ipGlobal()), outsideVolVars.extrusionFactor())
248 : -1.0*computeTpfaTransmissibility(fvGeometry, scvf, outsideScv, getPermeability_(problem, outsideVolVars, scvf.ipGlobal()), outsideVolVars.extrusionFactor());
255 tij = Extrusion::area(fvGeometry, scvf)*(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 Extrusion = Extrusion_t<GridGeometry>;
289 using GridView = typename GridGeometry::GridView;
290 using Element = typename GridView::template Codim<0>::Entity;
292 static constexpr int dim = GridView::dimension;
293 static constexpr int dimWorld = GridView::dimensionworld;
295 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
299 using Scalar = ScalarType;
301 using DiscretizationMethod = DiscretizationMethods::CCTpfa;
303 static constexpr DiscretizationMethod discMethod{};
306 using Cache = TpfaDarcysLawCache<ThisType, GridGeometry>;
318 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
319 static Scalar flux(const Problem& problem,
320 const Element& element,
321 const FVElementGeometry& fvGeometry,
322 const ElementVolumeVariables& elemVolVars,
323 const SubControlVolumeFace& scvf,
325 const ElementFluxVarsCache& elemFluxVarsCache)
327 static const bool gravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
329 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
332 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
333 const auto& insideVolVars = elemVolVars[insideScv];
334 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
339 const auto rho = [&]()
343 return outsideVolVars.density(phaseIdx);
346 else if (scvf.numOutsideScvs() == 1)
347 return (insideVolVars.density(phaseIdx) + outsideVolVars.density(phaseIdx))*0.5;
352 Scalar rho(insideVolVars.density(phaseIdx));
353 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
355 const auto outsideScvIdx = scvf.outsideScvIdx(i);
356 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
357 rho += outsideVolVars.density(phaseIdx);
359 return rho/(scvf.numOutsideScvs()+1);
363 const auto& tij = fluxVarsCache.advectionTij();
364 const auto& g = problem.spatialParams().gravity(scvf.ipGlobal());
367 const auto pInside = insideVolVars.pressure(phaseIdx);
368 const auto pOutside = [&]()
371 if (scvf.numOutsideScvs() == 1)
372 return outsideVolVars.pressure(phaseIdx);
378 Scalar sumPTi(tij*pInside);
381 sumPTi += rho*Extrusion::area(fvGeometry, scvf)
382 *insideVolVars.extrusionFactor()
383 *vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g);
385 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
387 const auto outsideScvIdx = scvf.outsideScvIdx(i);
388 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
389 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
390 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
391 sumTi += outsideFluxVarsCache.advectionTij();
392 sumPTi += outsideFluxVarsCache.advectionTij()*outsideVolVars.pressure(phaseIdx);
395 sumPTi += rho*Extrusion::area(fvGeometry, scvf)
396 *outsideVolVars.extrusionFactor()
397 *vtmv(flippedScvf.unitOuterNormal(), outsideVolVars.permeability(), g);
404 const auto alpha_inside = vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g)*insideVolVars.extrusionFactor();
406 Scalar flux = tij*(pInside - pOutside) + Extrusion::area(fvGeometry, scvf)*rho*alpha_inside;
409 if (!scvf.boundary() && scvf.numOutsideScvs() == 1)
411 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
412 const auto& outsideScvf = fvGeometry.flipScvf(scvf.index());
413 const auto outsideK = outsideVolVars.permeability();
414 const auto outsideTi = computeTpfaTransmissibility(fvGeometry, outsideScvf, outsideScv, outsideK, outsideVolVars.extrusionFactor());
415 const auto alpha_outside = vtmv(outsideScvf.unitOuterNormal(), outsideK, g)*outsideVolVars.extrusionFactor();
417 flux -= rho*tij/outsideTi*(alpha_inside + alpha_outside);
425 const auto pInside = insideVolVars.pressure(phaseIdx);
426 const auto pOutside = [&]()
429 if (scvf.numOutsideScvs() <= 1)
430 return outsideVolVars.pressure(phaseIdx);
435 const auto& insideFluxVarsCache = elemFluxVarsCache[scvf];
436 Scalar sumTi(insideFluxVarsCache.advectionTij());
437 Scalar sumPTi(insideFluxVarsCache.advectionTij()*pInside);
439 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
441 const auto outsideScvIdx = scvf.outsideScvIdx(i);
442 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
443 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
444 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
445 sumTi += outsideFluxVarsCache.advectionTij();
446 sumPTi += outsideFluxVarsCache.advectionTij()*outsideVolVars.pressure(phaseIdx);
453 return fluxVarsCache.advectionTij()*(pInside - pOutside);
459 template<class Problem, class ElementVolumeVariables>
460 static Scalar calculateTransmissibility(const Problem& problem,
461 const Element& element,
462 const FVElementGeometry& fvGeometry,
463 const ElementVolumeVariables& elemVolVars,
464 const SubControlVolumeFace& scvf)
468 const auto insideScvIdx = scvf.insideScvIdx();
469 const auto& insideScv = fvGeometry.scv(insideScvIdx);
470 const auto& insideVolVars = elemVolVars[insideScvIdx];
472 const Scalar ti = computeTpfaTransmissibility(fvGeometry, scvf, insideScv,
473 getPermeability_(problem, insideVolVars, scvf.ipGlobal()),
474 insideVolVars.extrusionFactor());
477 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
478 tij = Extrusion::area(fvGeometry, scvf)*ti;
483 const auto outsideScvIdx = scvf.outsideScvIdx();
486 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
487 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
488 const Scalar tj = computeTpfaTransmissibility(fvGeometry, fvGeometry.flipScvf(scvf.index()), outsideScv,
489 getPermeability_(problem, outsideVolVars, scvf.ipGlobal()),
490 outsideVolVars.extrusionFactor());
497 tij = Extrusion::area(fvGeometry, scvf)*(ti * tj)/(ti + tj);
504 template<class Problem, class VolumeVariables,
505 std::enable_if_t<!Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
506 static decltype(auto) getPermeability_(const Problem& problem,
507 const VolumeVariables& volVars,
508 const GlobalPosition& scvfIpGlobal)
509 { return volVars.permeability(); }
511 template<class Problem, class VolumeVariables,
512 std::enable_if_t<Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
513 static decltype(auto) getPermeability_(const Problem& problem,
514 const VolumeVariables& volVars,
515 const GlobalPosition& scvfIpGlobal)
516 { return problem.spatialParams().permeabilityAtPos(scvfIpGlobal); }
Darcy's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: flux/cctpfa/darcyslaw.hh:38
forward declaration of the method-specific implementation
Definition: flux/ccmpfa/darcyslaw.hh:27
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
Define some often used mathematical functions.
The available discretization methods in Dumux.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...