24#ifndef DUMUX_DISCRETIZATION_CC_TPFA_FACET_COUPLING_DARCYS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_TPFA_FACET_COUPLING_DARCYS_LAW_HH
29#include <dune/common/fmatrix.hh>
30#include <dune/common/dynmatrix.hh>
31#include <dune/common/dynvector.hh>
32#include <dune/common/float_cmp.hh>
45template<
class ScalarType,
class Gr
idGeometry,
bool isNetwork>
54template<
class AdvectionType,
class Gr
idGeometry,
bool isNetwork>
63template<
class ScalarType,
class Gr
idGeometry>
66 int(GridGeometry::GridView::dimensionworld) ) >;
72template<
class AdvectionType,
class Gr
idGeometry>
75 using Scalar =
typename AdvectionType::Scalar;
76 using FVElementGeometry =
typename GridGeometry::LocalView;
77 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
78 using Element =
typename GridGeometry::GridView::template Codim<0>::Entity;
87 static constexpr int insideTijIdx = 0;
88 static constexpr int outsideTijIdx = 1;
89 static constexpr int facetTijIdx = 2;
98 template<
class Problem,
class ElementVolumeVariables >
100 const Element& element,
101 const FVElementGeometry& fvGeometry,
102 const ElementVolumeVariables& elemVolVars,
103 const SubControlVolumeFace &scvf)
105 tij_ = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, g_);
113 {
return tij_[insideTijIdx]; }
117 {
return tij_[insideTijIdx]; }
121 {
return tij_[outsideTijIdx]; }
125 {
return tij_[facetTijIdx]; }
132 std::array<Scalar, 3> tij_;
133 GravityCoefficients g_;
140template<
class ScalarType,
class Gr
idGeometry>
146 using FVElementGeometry =
typename GridGeometry::LocalView;
147 using SubControlVolume =
typename GridGeometry::SubControlVolume;
148 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;
155 template<
class VolumeVariables,
class FacetVolVars>
156 static ScalarType computeFacetTransmissibility_(
const VolumeVariables& insideVolVars,
157 const FacetVolVars& facetVolVars,
158 const SubControlVolumeFace& scvf)
160 return 2.0*scvf.area()*insideVolVars.extrusionFactor()
161 /facetVolVars.extrusionFactor()
162 *
vtmv(scvf.unitOuterNormal(), facetVolVars.permeability(), scvf.unitOuterNormal());
173 using TijContainer =
typename Cache::AdvectionTransmissibilityContainer;
177 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache >
179 const Element& element,
180 const FVElementGeometry& fvGeometry,
181 const ElementVolumeVariables& elemVolVars,
182 const SubControlVolumeFace& scvf,
184 const ElementFluxVarsCache& elemFluxVarsCache)
186 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
187 return TpfaDarcysLaw::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
190 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
191 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
192 const auto pInside = insideVolVars.pressure(phaseIdx);
193 const auto pFacet = facetVolVars.pressure(phaseIdx);
196 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
197 Scalar flux = fluxVarsCache.advectionTijInside()*pInside + fluxVarsCache.advectionTijFacet()*pFacet;
200 static const Scalar gravity = getParamFromGroup<bool>(problem.paramGroup(),
"Problem.EnableGravity");
204 const auto& g = problem.spatialParams().gravity(scvf.ipGlobal());
205 const auto rho = 0.5*(insideVolVars.density(phaseIdx) + facetVolVars.density(phaseIdx));
206 const auto rhoTimesArea = rho*scvf.area();
207 const auto alpha_inside = rhoTimesArea*insideVolVars.extrusionFactor()
208 *
vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g);
210 flux += alpha_inside;
211 if (!scvf.boundary())
213 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
216 if ( problem.interiorBoundaryTypes(element, scvf).hasOnlyNeumann() )
218 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(),
"FacetCoupling.Xi", 1.0);
219 const auto alpha_facet = rhoTimesArea*insideVolVars.extrusionFactor()
220 *
vtmv(scvf.unitOuterNormal(), facetVolVars.permeability(), g);
221 const auto alpha_outside = rhoTimesArea*outsideVolVars.extrusionFactor()
222 *
vtmv(scvf.unitOuterNormal(), outsideVolVars.permeability(), g);
224 flux -= fluxVarsCache.gravityCoefficients()[0]*(xi*alpha_inside - alpha_facet + (1.0 - xi)*alpha_outside);
225 flux += fluxVarsCache.gravityCoefficients()[1]*(xi*alpha_outside - alpha_facet + (1.0 - xi)*alpha_inside);
229 flux += fluxVarsCache.advectionTijOutside()*outsideVolVars.pressure(phaseIdx);
235 return scvf.boundary() ? flux
236 : flux + fluxVarsCache.advectionTijOutside()*elemVolVars[scvf.outsideScvIdx()].pressure(phaseIdx);
241 template<
class Problem,
class ElementVolumeVariables >
243 const Element& element,
244 const FVElementGeometry& fvGeometry,
245 const ElementVolumeVariables& elemVolVars,
246 const SubControlVolumeFace& scvf)
248 typename Cache::GravityCoefficients g;
249 return calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, g);
254 template<
class Problem,
class ElementVolumeVariables >
256 const Element& element,
257 const FVElementGeometry& fvGeometry,
258 const ElementVolumeVariables& elemVolVars,
259 const SubControlVolumeFace& scvf,
260 typename Cache::GravityCoefficients& g)
263 if (!problem.couplingManager().isCoupled(element, scvf))
266 tij[Cache::insideTijIdx] = TpfaDarcysLaw::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
271 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(),
"FacetCoupling.Xi", 1.0);
273 const auto insideScvIdx = scvf.insideScvIdx();
274 const auto& insideScv = fvGeometry.scv(insideScvIdx);
275 const auto& insideVolVars = elemVolVars[insideScvIdx];
278 insideVolVars.permeability(),
279 insideVolVars.extrusionFactor());
282 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
285 if (iBcTypes.hasOnlyNeumann())
287 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
288 const auto wFacet = computeFacetTransmissibility_(insideVolVars, facetVolVars, scvf);
297 if (!scvf.boundary())
299 const auto outsideScvIdx = scvf.outsideScvIdx();
300 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
302 fvGeometry.scv(outsideScvIdx),
303 outsideVolVars.permeability(),
304 outsideVolVars.extrusionFactor());
306 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
313 const Scalar xiMinusOne = (xi - 1.0);
314 const Scalar a01 = xiMinusOne*wOut;
315 const Scalar a11 = xi*wOut + wFacet;
316 const Scalar detA = (xi*wIn + wFacet)*a11 - xiMinusOne*wIn*a01;
317 g[0] = wIn*a11/detA; g[1] = -wIn*a01/detA;
320 const Scalar factor = wIn * wFacet / ( wIn * wOut * ( 2.0 * xi - 1.0 ) + wFacet * ( xi * ( wIn + wOut ) + wFacet ) );
321 tij[Cache::insideTijIdx] = factor * ( wOut * xi + wFacet );
322 tij[Cache::outsideTijIdx] = factor * ( wOut * ( 1.0 - xi ) );
323 tij[Cache::facetTijIdx] = factor * ( - wOut - wFacet );
327 g[0] = wIn/(wIn+wFacet); g[1] = 0.0;
328 tij[Cache::insideTijIdx] = wFacet*g[0];
329 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
330 tij[Cache::outsideTijIdx] = 0.0;
336 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
337 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
338 tij[Cache::outsideTijIdx] = 0.0;
341 else if (iBcTypes.hasOnlyDirichlet())
343 tij[Cache::insideTijIdx] = wIn;
344 tij[Cache::outsideTijIdx] = 0.0;
345 tij[Cache::facetTijIdx] = -wIn;
348 DUNE_THROW(Dune::NotImplemented,
"Interior boundary types other than pure Dirichlet or Neumann");
358template<
class AdvectionType,
class Gr
idGeometry>
361 using Scalar =
typename AdvectionType::Scalar;
362 using FVElementGeometry =
typename GridGeometry::LocalView;
363 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
364 using Element =
typename GridGeometry::GridView::template Codim<0>::Entity;
373 static constexpr int insideTijIdx = 0;
374 static constexpr int facetTijIdx = 1;
380 template<
class Problem,
class ElementVolumeVariables >
382 const Element& element,
383 const FVElementGeometry& fvGeometry,
384 const ElementVolumeVariables& elemVolVars,
385 const SubControlVolumeFace &scvf)
387 tij_ = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
395 {
return tij_[insideTijIdx]; }
399 {
return tij_[insideTijIdx]; }
403 {
return tij_[facetTijIdx]; }
406 AdvectionTransmissibilityContainer tij_;
413template<
class ScalarType,
class Gr
idGeometry>
419 using FVElementGeometry =
typename GridGeometry::LocalView;
420 using SubControlVolume =
typename GridGeometry::SubControlVolume;
421 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
423 using GridView =
typename GridGeometry::GridView;
424 using Element =
typename GridView::template Codim<0>::Entity;
425 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
435 using TijContainer =
typename Cache::AdvectionTransmissibilityContainer;
438 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache >
440 const Element& element,
441 const FVElementGeometry& fvGeometry,
442 const ElementVolumeVariables& elemVolVars,
443 const SubControlVolumeFace& scvf,
445 const ElementFluxVarsCache& elemFluxVarsCache)
447 static const Scalar gravity = getParamFromGroup<bool>(problem.paramGroup(),
"Problem.EnableGravity");
449 DUNE_THROW(Dune::NotImplemented,
"Gravity for darcys law with facet coupling on surface grids");
451 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
452 return TpfaDarcysLaw::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
455 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
456 const auto pInside = insideVolVars.pressure(phaseIdx);
457 const auto pFacet = problem.couplingManager().getLowDimVolVars(element, scvf).pressure(phaseIdx);
460 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
462 return fluxVarsCache.advectionTijInside()*pInside + fluxVarsCache.advectionTijFacet()*pFacet;
464 return fluxVarsCache.advectionTijInside()*pInside + fluxVarsCache.advectionTijFacet()*pFacet;
469 template<
class Problem,
class ElementVolumeVariables >
471 const Element& element,
472 const FVElementGeometry& fvGeometry,
473 const ElementVolumeVariables& elemVolVars,
474 const SubControlVolumeFace& scvf)
477 if (!problem.couplingManager().isCoupled(element, scvf))
480 tij[Cache::insideTijIdx] = TpfaDarcysLaw::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
485 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(),
"FacetCoupling.Xi", 1.0);
490 if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
491 DUNE_THROW(Dune::InvalidStateException,
"Xi != 1.0 cannot be used on surface grids");
493 const auto area = scvf.area();
494 const auto insideScvIdx = scvf.insideScvIdx();
495 const auto& insideScv = fvGeometry.scv(insideScvIdx);
496 const auto& insideVolVars = elemVolVars[insideScvIdx];
500 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
503 if (iBcTypes.hasOnlyNeumann())
505 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
510 const auto wFacet = 2.0*area*insideVolVars.extrusionFactor()
511 /sqrt(facetVolVars.extrusionFactor())
512 *
vtmv(scvf.unitOuterNormal(), facetVolVars.permeability(), scvf.unitOuterNormal());
515 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
516 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
518 else if (iBcTypes.hasOnlyDirichlet())
520 tij[Cache::insideTijIdx] = wIn;
521 tij[Cache::facetTijIdx] = -wIn;
524 DUNE_THROW(Dune::NotImplemented,
"Interior boundary types other than pure Dirichlet or Neumann");
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
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:47
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:840
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Class that fills the cache corresponding to tpfa Darcy's Law.
Definition: flux/cctpfa/darcyslaw.hh:69
Specialization of the CCTpfaDarcysLaw grids where dim=dimWorld.
Definition: flux/cctpfa/darcyslaw.hh:129
Specialization of the CCTpfaDarcysLaw grids where dim < dimWorld (network/surface grids)
Definition: flux/cctpfa/darcyslaw.hh:283
Forward declaration of the implementation.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:46
The cache corresponding to tpfa Darcy's Law with facet coupling.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:55
const GravityCoefficients & gravityCoefficients() const
return the coefficients for the computation of gravity at the scvf
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:128
std::array< Scalar, 3 > AdvectionTransmissibilityContainer
Export transmissibility storage type.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:92
Scalar advectionTijInside() const
returns the transmissibility associated with the inside cell
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:116
std::array< Scalar, 2 > GravityCoefficients
Export the type used for the gravity coefficients.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:95
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:99
Scalar advectionTijOutside() const
returns the transmissibility associated with the outside cell
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:120
Scalar advectionTijFacet() const
returns the transmissibility associated with the outside cell
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:124
Scalar advectionTij() const
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:112
Specialization of the CCTpfaDarcysLaw grids where dim=dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:142
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:178
ScalarType Scalar
state the scalar type of the law
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:167
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:242
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:255
typename Cache::AdvectionTransmissibilityContainer TijContainer
export the type used to store transmissibilities
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:173
std::array< Scalar, 2 > AdvectionTransmissibilityContainer
Export transmissibility storage type.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:377
Scalar advectionTij() const
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:394
Scalar advectionTijInside() const
returns the transmissibility associated with the inside cell
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:398
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:381
Scalar advectionTijFacet() const
returns the transmissibility associated with the outside cell
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:402
Specialization of the CCTpfaDarcysLaw grids where 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:470
typename Cache::AdvectionTransmissibilityContainer TijContainer
export the type used to store transmissibilities
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:435
ScalarType Scalar
state the scalar type of the law
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:429
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:439
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.