24#ifndef DUMUX_DISCRETIZATION_CC_TPFA_FACET_COUPLING_DARCYS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_TPFA_FACET_COUPLING_DARCYS_LAW_HH
30#include <dune/common/float_cmp.hh>
44template<
class ScalarType,
class Gr
idGeometry,
bool isNetwork>
53template<
class AdvectionType,
class Gr
idGeometry,
bool isNetwork>
62template<
class ScalarType,
class Gr
idGeometry>
65 int(GridGeometry::GridView::dimensionworld) ) >;
71template<
class AdvectionType,
class Gr
idGeometry>
74 using Scalar =
typename AdvectionType::Scalar;
75 using FVElementGeometry =
typename GridGeometry::LocalView;
76 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
77 using Element =
typename GridGeometry::GridView::template Codim<0>::Entity;
86 static constexpr int insideTijIdx = 0;
87 static constexpr int outsideTijIdx = 1;
88 static constexpr int facetTijIdx = 2;
97 template<
class Problem,
class ElementVolumeVariables >
99 const Element& element,
100 const FVElementGeometry& fvGeometry,
101 const ElementVolumeVariables& elemVolVars,
102 const SubControlVolumeFace &scvf)
104 tij_ = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, g_);
112 {
return tij_[insideTijIdx]; }
116 {
return tij_[insideTijIdx]; }
120 {
return tij_[outsideTijIdx]; }
124 {
return tij_[facetTijIdx]; }
131 std::array<Scalar, 3> tij_;
132 GravityCoefficients g_;
139template<
class ScalarType,
class Gr
idGeometry>
145 using FVElementGeometry =
typename GridGeometry::LocalView;
146 using SubControlVolume =
typename GridGeometry::SubControlVolume;
147 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
150 using GridView =
typename GridGeometry::GridView;
151 using Element =
typename GridView::template Codim<0>::Entity;
152 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
163 using TijContainer =
typename Cache::AdvectionTransmissibilityContainer;
167 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache >
169 const Element& element,
170 const FVElementGeometry& fvGeometry,
171 const ElementVolumeVariables& elemVolVars,
172 const SubControlVolumeFace& scvf,
174 const ElementFluxVarsCache& elemFluxVarsCache)
176 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
177 return TpfaDarcysLaw::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
180 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
181 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
182 const auto pInside = insideVolVars.pressure(phaseIdx);
183 const auto pFacet = facetVolVars.pressure(phaseIdx);
186 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
187 Scalar flux = fluxVarsCache.advectionTijInside()*pInside + fluxVarsCache.advectionTijFacet()*pFacet;
190 static const Scalar gravity = getParamFromGroup<bool>(problem.paramGroup(),
"Problem.EnableGravity");
194 const auto& g = problem.spatialParams().gravity(scvf.ipGlobal());
195 const auto rho = 0.5*(insideVolVars.density(phaseIdx) + facetVolVars.density(phaseIdx));
196 const auto rhoTimesArea = rho*Extrusion::area(fvGeometry, scvf);
197 const auto alpha_inside = rhoTimesArea*insideVolVars.extrusionFactor()
198 *
vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g);
200 flux += alpha_inside;
201 if (!scvf.boundary())
203 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
206 if ( problem.interiorBoundaryTypes(element, scvf).hasOnlyNeumann() )
208 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(),
"FacetCoupling.Xi", 1.0);
209 const auto alpha_facet = rhoTimesArea*insideVolVars.extrusionFactor()
210 *
vtmv(scvf.unitOuterNormal(), facetVolVars.permeability(), g);
211 const auto alpha_outside = rhoTimesArea*outsideVolVars.extrusionFactor()
212 *
vtmv(scvf.unitOuterNormal(), outsideVolVars.permeability(), g);
214 flux -= fluxVarsCache.gravityCoefficients()[0]*(xi*alpha_inside - alpha_facet + (1.0 - xi)*alpha_outside);
215 flux += fluxVarsCache.gravityCoefficients()[1]*(xi*alpha_outside - alpha_facet + (1.0 - xi)*alpha_inside);
219 flux += fluxVarsCache.advectionTijOutside()*outsideVolVars.pressure(phaseIdx);
225 return scvf.boundary() ? flux
226 : flux + fluxVarsCache.advectionTijOutside()*elemVolVars[scvf.outsideScvIdx()].pressure(phaseIdx);
231 template<
class Problem,
class ElementVolumeVariables >
233 const Element& element,
234 const FVElementGeometry& fvGeometry,
235 const ElementVolumeVariables& elemVolVars,
236 const SubControlVolumeFace& scvf)
238 typename Cache::GravityCoefficients g;
239 return calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, g);
244 template<
class Problem,
class ElementVolumeVariables >
246 const Element& element,
247 const FVElementGeometry& fvGeometry,
248 const ElementVolumeVariables& elemVolVars,
249 const SubControlVolumeFace& scvf,
250 typename Cache::GravityCoefficients& g)
253 if (!problem.couplingManager().isCoupled(element, scvf))
256 tij[Cache::insideTijIdx] = TpfaDarcysLaw::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
261 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(),
"FacetCoupling.Xi", 1.0);
263 const auto insideScvIdx = scvf.insideScvIdx();
264 const auto& insideScv = fvGeometry.scv(insideScvIdx);
265 const auto& insideVolVars = elemVolVars[insideScvIdx];
266 const auto wIn = Extrusion::area(fvGeometry, scvf)
268 insideVolVars.permeability(),
269 insideVolVars.extrusionFactor());
272 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
275 if (iBcTypes.hasOnlyNeumann())
277 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
278 const auto wFacet = 2.0*Extrusion::area(fvGeometry, scvf)*insideVolVars.extrusionFactor()
279 /facetVolVars.extrusionFactor()
280 *
vtmv(scvf.unitOuterNormal(), facetVolVars.permeability(), scvf.unitOuterNormal());
289 if (!scvf.boundary())
291 const auto outsideScvIdx = scvf.outsideScvIdx();
292 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
293 const auto wOut = -1.0*Extrusion::area(fvGeometry, scvf)
295 outsideVolVars.permeability(),
296 outsideVolVars.extrusionFactor());
298 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
305 const Scalar xiMinusOne = (xi - 1.0);
306 const Scalar a01 = xiMinusOne*wOut;
307 const Scalar a11 = xi*wOut + wFacet;
308 const Scalar detA = (xi*wIn + wFacet)*a11 - xiMinusOne*wIn*a01;
309 g[0] = wIn*a11/detA; g[1] = -wIn*a01/detA;
312 const Scalar factor = wIn * wFacet / ( wIn * wOut * ( 2.0 * xi - 1.0 ) + wFacet * ( xi * ( wIn + wOut ) + wFacet ) );
313 tij[Cache::insideTijIdx] = factor * ( wOut * xi + wFacet );
314 tij[Cache::outsideTijIdx] = factor * ( wOut * ( 1.0 - xi ) );
315 tij[Cache::facetTijIdx] = factor * ( - wOut - wFacet );
319 g[0] = wIn/(wIn+wFacet); g[1] = 0.0;
320 tij[Cache::insideTijIdx] = wFacet*g[0];
321 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
322 tij[Cache::outsideTijIdx] = 0.0;
328 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
329 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
330 tij[Cache::outsideTijIdx] = 0.0;
333 else if (iBcTypes.hasOnlyDirichlet())
335 tij[Cache::insideTijIdx] = wIn;
336 tij[Cache::outsideTijIdx] = 0.0;
337 tij[Cache::facetTijIdx] = -wIn;
340 DUNE_THROW(Dune::NotImplemented,
"Interior boundary types other than pure Dirichlet or Neumann");
350template<
class AdvectionType,
class Gr
idGeometry>
353 using Scalar =
typename AdvectionType::Scalar;
354 using FVElementGeometry =
typename GridGeometry::LocalView;
355 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
356 using Element =
typename GridGeometry::GridView::template Codim<0>::Entity;
365 static constexpr int insideTijIdx = 0;
366 static constexpr int facetTijIdx = 1;
375 template<
class Problem,
class ElementVolumeVariables >
377 const Element& element,
378 const FVElementGeometry& fvGeometry,
379 const ElementVolumeVariables& elemVolVars,
380 const SubControlVolumeFace &scvf)
382 tij_ = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, g_);
390 {
return tij_[insideTijIdx]; }
394 {
return tij_[insideTijIdx]; }
398 {
return tij_[facetTijIdx]; }
405 AdvectionTransmissibilityContainer tij_;
406 GravityCoefficients g_;
413template<
class ScalarType,
class Gr
idGeometry>
419 using FVElementGeometry =
typename GridGeometry::LocalView;
420 using SubControlVolume =
typename GridGeometry::SubControlVolume;
421 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
424 using GridView =
typename GridGeometry::GridView;
425 using Element =
typename GridView::template Codim<0>::Entity;
426 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
437 using TijContainer =
typename Cache::AdvectionTransmissibilityContainer;
440 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache >
442 const Element& element,
443 const FVElementGeometry& fvGeometry,
444 const ElementVolumeVariables& elemVolVars,
445 const SubControlVolumeFace& scvf,
447 const ElementFluxVarsCache& elemFluxVarsCache)
449 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
450 return TpfaDarcysLaw::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
455 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(),
"FacetCoupling.Xi", 1.0);
456 if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
457 DUNE_THROW(Dune::InvalidStateException,
"Xi != 1.0 cannot be used on surface grids");
460 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
461 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
462 const auto pInside = insideVolVars.pressure(phaseIdx);
463 const auto pFacet = facetVolVars.pressure(phaseIdx);
466 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
467 Scalar flux = fluxVarsCache.advectionTijInside()*pInside + fluxVarsCache.advectionTijFacet()*pFacet;
469 static const Scalar gravity = getParamFromGroup<bool>(problem.paramGroup(),
"Problem.EnableGravity");
473 const auto& g = problem.spatialParams().gravity(scvf.ipGlobal());
474 const auto rho = 0.5*(insideVolVars.density(phaseIdx) + facetVolVars.density(phaseIdx));
475 const auto rhoTimesArea = rho*Extrusion::area(fvGeometry, scvf);
476 const auto alpha_inside = rhoTimesArea*insideVolVars.extrusionFactor()
477 *
vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g);
479 flux += alpha_inside;
482 if ( !scvf.boundary() && problem.interiorBoundaryTypes(element, scvf).hasOnlyNeumann() )
484 const auto alpha_facet = rhoTimesArea*insideVolVars.extrusionFactor()
485 *
vtmv(scvf.unitOuterNormal(), facetVolVars.permeability(), g);
487 flux -= fluxVarsCache.gravityCoefficients()[0]*(alpha_inside - alpha_facet);
496 template<
class Problem,
class ElementVolumeVariables >
498 const Element& element,
499 const FVElementGeometry& fvGeometry,
500 const ElementVolumeVariables& elemVolVars,
501 const SubControlVolumeFace& scvf)
503 typename Cache::GravityCoefficients g;
504 return calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, g);
509 template<
class Problem,
class ElementVolumeVariables >
511 const Element& element,
512 const FVElementGeometry& fvGeometry,
513 const ElementVolumeVariables& elemVolVars,
514 const SubControlVolumeFace& scvf,
515 typename Cache::GravityCoefficients& g)
518 if (!problem.couplingManager().isCoupled(element, scvf))
521 tij[Cache::insideTijIdx] = TpfaDarcysLaw::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
526 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(),
"FacetCoupling.Xi", 1.0);
531 if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
532 DUNE_THROW(Dune::InvalidStateException,
"Xi != 1.0 cannot be used on surface grids");
534 const auto area = Extrusion::area(fvGeometry, scvf);
535 const auto insideScvIdx = scvf.insideScvIdx();
536 const auto& insideScv = fvGeometry.scv(insideScvIdx);
537 const auto& insideVolVars = elemVolVars[insideScvIdx];
538 const auto wIn = area*
computeTpfaTransmissibility(fvGeometry, scvf, insideScv, insideVolVars.permeability(), insideVolVars.extrusionFactor());
541 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
544 if (iBcTypes.hasOnlyNeumann())
546 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
551 const auto wFacet = 2.0*area*insideVolVars.extrusionFactor()
552 /sqrt(facetVolVars.extrusionFactor())
553 *
vtmv(scvf.unitOuterNormal(), facetVolVars.permeability(), scvf.unitOuterNormal());
556 g[0] = wIn/(wIn+wFacet);
557 tij[Cache::insideTijIdx] = wFacet*g[0];
558 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
560 else if (iBcTypes.hasOnlyDirichlet())
562 tij[Cache::insideTijIdx] = wIn;
563 tij[Cache::facetTijIdx] = -wIn;
566 DUNE_THROW(Dune::NotImplemented,
"Interior boundary types other than pure Dirichlet or Neumann");
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.
Tensor::field_type computeTpfaTransmissibility(const SubControlVolumeFace &scvf, const SubControlVolume &scv, const Tensor &T, typename SubControlVolume::Traits::Scalar extrusionFactor)
Free function to evaluate the Tpfa transmissibility associated with the flux (in the form of flux = T...
Definition: tpfa/computetransmissibility.hh:48
Dune::DenseMatrix< MAT >::value_type vtmv(const Dune::DenseVector< V1 > &v1, const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V2 > &v2)
Evaluates the scalar product of a vector v2, projected by a matrix M, with a vector v1.
Definition: math.hh:863
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
Class that fills the cache corresponding to tpfa Darcy's Law.
Definition: flux/cctpfa/darcyslaw.hh:70
Specialization of the CCTpfaDarcysLaw grids where dim=dimWorld.
Definition: flux/cctpfa/darcyslaw.hh:130
Specialization of the CCTpfaDarcysLaw grids where dim < dimWorld (network/surface grids)
Definition: flux/cctpfa/darcyslaw.hh:295
Forward declaration of the implementation.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:45
The cache corresponding to tpfa Darcy's Law with facet coupling.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:54
const GravityCoefficients & gravityCoefficients() const
return the coefficients for the computation of gravity at the scvf
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:127
std::array< Scalar, 3 > AdvectionTransmissibilityContainer
Export transmissibility storage type.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:91
Scalar advectionTijInside() const
returns the transmissibility associated with the inside cell
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:115
std::array< Scalar, 2 > GravityCoefficients
Export the type used for the gravity coefficients.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:94
void updateAdvection(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
update subject to a given problem
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:98
Scalar advectionTijOutside() const
returns the transmissibility associated with the outside cell
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:119
Scalar advectionTijFacet() const
returns the transmissibility associated with the outside cell
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:123
Scalar advectionTij() const
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:111
Specialization of CCTpfaFacetCouplingDarcysLawImpl for dim=dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:141
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, int phaseIdx, const ElementFluxVarsCache &elemFluxVarsCache)
Compute the advective flux.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:168
ScalarType Scalar
state the scalar type of the law
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:156
static TijContainer calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:232
static TijContainer calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, typename Cache::GravityCoefficients &g)
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:245
typename Cache::AdvectionTransmissibilityContainer TijContainer
export the type used to store transmissibilities
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:163
std::array< Scalar, 1 > GravityCoefficients
Export the type used for the gravity coefficients.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:372
std::array< Scalar, 2 > AdvectionTransmissibilityContainer
Export transmissibility storage type.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:369
Scalar advectionTij() const
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:389
Scalar advectionTijInside() const
returns the transmissibility associated with the inside cell
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:393
void updateAdvection(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
update subject to a given problem
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:376
const GravityCoefficients & gravityCoefficients() const
return the coefficients for the computation of gravity at the scvf
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:401
Scalar advectionTijFacet() const
returns the transmissibility associated with the outside cell
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:397
Specialization of CCTpfaFacetCouplingDarcysLawImpl for dim<dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:415
static TijContainer calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:497
typename Cache::AdvectionTransmissibilityContainer TijContainer
export the type used to store transmissibilities
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:437
ScalarType Scalar
state the scalar type of the law
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:430
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, int phaseIdx, const ElementFluxVarsCache &elemFluxVarsCache)
Compute the advective flux.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:441
static TijContainer calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, typename Cache::GravityCoefficients &g)
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:510
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...
Darcy's law for cell-centered finite volume schemes with two-point flux approximation.