25#ifndef DUMUX_MD_FREEFLOW_POROUSMEDIUM_COUPLINGCONDITIONS_HH
26#define DUMUX_MD_FREEFLOW_POROUSMEDIUM_COUPLINGCONDITIONS_HH
30#include <dune/common/exceptions.hh>
60 harmonic, arithmetic, ffOnly, pmOnly
69 if (diffusionCoefficientAveragingType ==
"Harmonic")
71 else if (diffusionCoefficientAveragingType ==
"Arithmetic")
73 else if (diffusionCoefficientAveragingType ==
"FreeFlowOnly")
75 else if (diffusionCoefficientAveragingType ==
"PorousMediumOnly")
78 DUNE_THROW(Dune::IOError,
"Unknown DiffusionCoefficientAveragingType");
83template<
class MDTraits,
class CouplingManager,
bool enableEnergyBalance,
bool isCompositional>
91template<
class MDTraits,
class CouplingManager>
103template<
class MDTraits,
class CouplingManager>
106 using Scalar =
typename MDTraits::Scalar;
108 template<std::
size_t id>
using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
110 template<std::
size_t id>
using Element =
typename GridGeometry<id>::GridView::template Codim<0>::Entity;
111 template<std::
size_t id>
using FVElementGeometry =
typename GridGeometry<id>::LocalView;
112 template<std::
size_t id>
using SubControlVolumeFace =
typename GridGeometry<id>::LocalView::SubControlVolumeFace;
113 template<std::
size_t id>
using SubControlVolume =
typename GridGeometry<id>::LocalView::SubControlVolume;
120 template<std::
size_t id>
using GlobalPosition =
typename Element<id>::Geometry::GlobalCoordinate;
121 template<std::
size_t id>
using NumEqVector =
typename Problem<id>::Traits::NumEqVector;
133 static constexpr bool adapterUsed = ModelTraits<porousMediumIndex>::numFluidPhases() > 1;
138 "All submodels must both be either isothermal or non-isothermal");
141 FluidSystem<porousMediumIndex>>::value,
142 "All submodels must use the same fluid system");
150 template<std::
size_t i>
151 static constexpr auto couplingPhaseIdx(Dune::index_constant<i>
id,
int coupledPhaseIdx = 0)
152 {
return IndexHelper::couplingPhaseIdx(
id, coupledPhaseIdx); }
157 template<std::
size_t i>
159 {
return IndexHelper::couplingCompIdx(
id, coupledCompIdx); }
164 template<
class Context>
166 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
167 const Context& context)
168 {
return context.volVars.permeability(); }
177 template<
class Context>
179 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
180 const ElementVolumeVariables<freeFlowMomentumIndex>& elemVolVars,
181 const Context& context)
184 NumEqVector<freeFlowMomentumIndex> momentumFlux(0.0);
186 if (!scvf.isFrontal())
191 if(numPhasesDarcy > 1)
192 momentumFlux[scvf.normalAxis()] = context.volVars.pressure(pmPhaseIdx);
199 momentumFlux[scvf.normalAxis()] -= elemVolVars.gridVolVars().problem().referencePressure(fvGeometry.element(), fvGeometry, scvf);
201 momentumFlux[scvf.normalAxis()] *= scvf.directionSign();
209 static Scalar
advectiveFlux(
const Scalar insideQuantity,
const Scalar outsideQuantity,
const Scalar volumeFlow,
bool insideIsUpstream)
211 const Scalar upwindWeight = 1.0;
214 return (upwindWeight * insideQuantity + (1.0 - upwindWeight) * outsideQuantity) * volumeFlow;
216 return (upwindWeight * outsideQuantity + (1.0 - upwindWeight) * insideQuantity) * volumeFlow;
224 template<std::
size_t i, std::
size_t j>
226 Dune::index_constant<j> domainJ,
227 const Scalar insideDistance,
228 const Scalar outsideDistance,
229 const Scalar avgQuantityI,
230 const Scalar avgQuantityJ,
231 const DiffusionCoefficientAveragingType diffCoeffAvgType)
const
233 const Scalar totalDistance = insideDistance + outsideDistance;
234 if(diffCoeffAvgType == DiffusionCoefficientAveragingType::harmonic)
236 return harmonicMean(avgQuantityI, avgQuantityJ, insideDistance, outsideDistance)
239 else if(diffCoeffAvgType == DiffusionCoefficientAveragingType::arithmetic)
241 return arithmeticMean(avgQuantityI, avgQuantityJ, insideDistance, outsideDistance)
244 else if(diffCoeffAvgType == DiffusionCoefficientAveragingType::ffOnly)
246 ? avgQuantityI / totalDistance
247 : avgQuantityJ / totalDistance;
251 ? avgQuantityI / totalDistance
252 : avgQuantityJ / totalDistance;
258 template<
class Scv,
class Scvf>
261 return (scv.dofPosition() - scvf.ipGlobal()).two_norm();
267 template<std::
size_t i, std::
size_t j,
bool isNI = enableEnergyBalance,
typename std::enable_if_t<isNI,
int> = 0>
269 Dune::index_constant<j> domainJ,
270 const FVElementGeometry<i>& fvGeometryI,
271 const FVElementGeometry<j>& fvGeometryJ,
272 const SubControlVolumeFace<i>& scvfI,
273 const SubControlVolume<i>& scvI,
274 const SubControlVolume<j>& scvJ,
275 const VolumeVariables<i>& volVarsI,
276 const VolumeVariables<j>& volVarsJ,
277 const DiffusionCoefficientAveragingType diffCoeffAvgType)
const
279 const Scalar insideDistance =
getDistance_(scvI, scvfI);
280 const Scalar outsideDistance =
getDistance_(scvJ, scvfI);
282 const Scalar deltaT = volVarsJ.temperature() - volVarsI.temperature();
285 insideDistance, outsideDistance,
286 volVarsI.effectiveThermalConductivity(), volVarsJ.effectiveThermalConductivity(),
290 return -tij * deltaT;
296 template<
class CouplingContext>
298 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
299 const ElementVolumeVariables<freeFlowMomentumIndex>& elemVolVars,
300 const CouplingContext& context)
302 GlobalPosition<freeFlowMomentumIndex> velocity(0.0);
303 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
304 velocity[scv.dofAxis()] = elemVolVars[scv].velocity();
305 const auto& pmScvf = context.fvGeometry.scvf(context.porousMediumScvfIdx);
313 const SubControlVolumeFace<porousMediumIndex>& scvf,
314 const VolumeVariables<porousMediumIndex>& volVars,
315 const typename Element<freeFlowMomentumIndex>::Geometry::GlobalCoordinate& couplingPhaseVelocity,
318 DUNE_THROW(Dune::NotImplemented,
"Forchheimer's law pressure reconstruction");
325 const SubControlVolumeFace<porousMediumIndex>& scvf,
326 const VolumeVariables<porousMediumIndex>& volVars,
327 const typename Element<freeFlowMomentumIndex>::Geometry::GlobalCoordinate& gravity,
328 const typename Element<freeFlowMomentumIndex>::Geometry::GlobalCoordinate& couplingPhaseVelocity,
332 const Scalar couplingPhaseCellCenterPressure = volVars.pressure(pmPhaseIdx);
333 const Scalar couplingPhaseMobility = volVars.mobility(pmPhaseIdx);
334 const Scalar couplingPhaseDensity = volVars.density(pmPhaseIdx);
335 const auto K = volVars.permeability();
341 const auto alpha =
vtmv(scvf.unitOuterNormal(), K, gravity);
343 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
346 return (-1/couplingPhaseMobility * (scvf.unitOuterNormal() * couplingPhaseVelocity) + couplingPhaseDensity * alpha)/ti
347 + couplingPhaseCellCenterPressure;
355template<
class MDTraits,
class CouplingManager,
bool enableEnergyBalance>
360 using Scalar =
typename MDTraits::Scalar;
363 template<std::
size_t id>
364 using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
367 template<std::
size_t id>
using Element =
typename GridGeometry<id>::GridView::template Codim<0>::Entity;
368 template<std::
size_t id>
using FVElementGeometry =
typename GridGeometry<id>::LocalView;
369 template<std::
size_t id>
using SubControlVolumeFace =
typename GridGeometry<id>::LocalView::SubControlVolumeFace;
370 template<std::
size_t id>
using SubControlVolume =
typename GridGeometry<id>::LocalView::SubControlVolume;
376 "Darcy Model must not be compositional");
381 using ParentType::ParentType;
382 using ParentType::couplingPhaseIdx;
387 template<std::
size_t i, std::
size_t j,
class CouplingContext>
389 Dune::index_constant<j> domainJ,
390 const FVElementGeometry<i>& fvGeometry,
391 const SubControlVolumeFace<i>& scvf,
392 const ElementVolumeVariables<i>& insideVolVars,
393 const CouplingContext& context)
395 const Scalar normalVelocity = context.velocity * scvf.unitOuterNormal();
397 const Scalar stokesDensity = context.volVars.density();
398 const bool insideIsUpstream = normalVelocity > 0.0;
399 return ParentType::advectiveFlux(darcyDensity, stokesDensity, normalVelocity, insideIsUpstream);
405 template<
bool isNI = enableEnergyBalance,
typename std::enable_if_t<isNI,
int> = 0>
407 const FVElementGeometry<ParentType::porousMediumIndex>& fvGeometry,
408 const ElementVolumeVariables<ParentType::porousMediumIndex>& darcyElemVolVars,
409 const SubControlVolumeFace<ParentType::porousMediumIndex>& scvf,
410 const DiffusionCoefficientAveragingType diffCoeffAvgType = DiffusionCoefficientAveragingType::ffOnly)
const
412 DUNE_THROW(Dune::NotImplemented,
"Energy coupling condition");
418 template<
bool isNI = enableEnergyBalance,
typename std::enable_if_t<isNI,
int> = 0>
420 const FVElementGeometry<ParentType::freeFlowMassIndex>& fvGeometry,
421 const ElementVolumeVariables<ParentType::freeFlowMassIndex>& stokesElemVolVars,
422 const SubControlVolumeFace<ParentType::freeFlowMassIndex>& scvf,
423 const DiffusionCoefficientAveragingType diffCoeffAvgType = DiffusionCoefficientAveragingType::ffOnly)
const
425 DUNE_THROW(Dune::NotImplemented,
"Energy coupling condition");
435 template<std::
size_t i, std::
size_t j,
bool isNI = enableEnergyBalance,
typename std::enable_if_t<isNI,
int> = 0>
436 Scalar energyFlux_(Dune::index_constant<i> domainI,
437 Dune::index_constant<j> domainJ,
438 const FVElementGeometry<i>& insideFvGeometry,
439 const FVElementGeometry<j>& outsideFvGeometry,
440 const SubControlVolumeFace<i>& scvf,
441 const VolumeVariables<i>& insideVolVars,
442 const VolumeVariables<j>& outsideVolVars,
443 const Scalar velocity,
444 const bool insideIsUpstream,
445 const DiffusionCoefficientAveragingType diffCoeffAvgType)
const
449 const auto& insideScv = (*scvs(insideFvGeometry).begin());
450 const auto& outsideScv = (*scvs(outsideFvGeometry).begin());
453 const Scalar insideTerm = insideVolVars.density(couplingPhaseIdx(domainI)) * insideVolVars.enthalpy(couplingPhaseIdx(domainI));
454 const Scalar outsideTerm = outsideVolVars.density(couplingPhaseIdx(domainJ)) * outsideVolVars.enthalpy(couplingPhaseIdx(domainJ));
456 flux += this->advectiveFlux(insideTerm, outsideTerm, velocity, insideIsUpstream);
458 flux += this->conductiveEnergyFlux_(domainI, domainJ, insideFvGeometry, outsideFvGeometry, scvf, insideScv, outsideScv, insideVolVars, outsideVolVars, diffCoeffAvgType);
Define some often used mathematical functions.
Darcy's law specialized for different discretization schemes This file contains the data which is req...
Forchheimer's law specialized for different discretization schemes This file contains the data which ...
Fick's law specialized for different discretization schemes. This file contains the data which is req...
Index helpers for the free-flow/porous-medium-flow coupling.
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
constexpr Scalar harmonicMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) harmonic mean of two scalar values.
Definition: math.hh:69
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
constexpr Scalar arithmeticMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) arithmetic mean of two scalar values.
Definition: math.hh:50
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
static constexpr auto freeFlowMassIndex
Definition: multidomain/boundary/freeflowporenetwork/couplingmanager.hh:47
static constexpr auto freeFlowMomentumIndex
Definition: multidomain/boundary/freeflowporenetwork/couplingmanager.hh:46
static constexpr auto porousMediumIndex
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:47
Traits class encapsulating model specifications.
Definition: common/properties.hh:51
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:55
Definition: common/properties.hh:100
The type for a global container for the volume variables.
Definition: common/properties.hh:107
The type for the calculation the advective fluxes.
Definition: common/properties.hh:139
The type of the fluid system to use.
Definition: common/properties.hh:160
forward declaration of the method-specific implementation
Definition: flux/box/darcyslaw.hh:44
forward declare
Definition: forchheimerslaw_fwd.hh:39
This structs holds a set of options which allow to modify the Stokes-Darcy coupling mechanism during ...
Definition: freeflowporousmedium/couplingconditions.hh:52
static DiffusionCoefficientAveragingType stringToEnum(DiffusionCoefficientAveragingType, const std::string &diffusionCoefficientAveragingType)
Convenience function to convert user input given as std::string to the corresponding enum class used ...
Definition: freeflowporousmedium/couplingconditions.hh:67
DiffusionCoefficientAveragingType
Defines which kind of averanging of diffusion coefficients (moleculat diffusion or thermal conductanc...
Definition: freeflowporousmedium/couplingconditions.hh:59
Definition: freeflowporousmedium/couplingconditions.hh:84
A base class which provides some common methods used for Stokes-Darcy coupling.
Definition: freeflowporousmedium/couplingconditions.hh:105
static constexpr auto freeFlowMomentumIndex
Definition: freeflowporousmedium/couplingconditions.hh:124
static Scalar advectiveFlux(const Scalar insideQuantity, const Scalar outsideQuantity, const Scalar volumeFlow, bool insideIsUpstream)
Evaluate an advective flux across the interface and consider upwinding.
Definition: freeflowporousmedium/couplingconditions.hh:209
Scalar conductiveEnergyFlux_(Dune::index_constant< i > domainI, Dune::index_constant< j > domainJ, const FVElementGeometry< i > &fvGeometryI, const FVElementGeometry< j > &fvGeometryJ, const SubControlVolumeFace< i > &scvfI, const SubControlVolume< i > &scvI, const SubControlVolume< j > &scvJ, const VolumeVariables< i > &volVarsI, const VolumeVariables< j > &volVarsJ, const DiffusionCoefficientAveragingType diffCoeffAvgType) const
Returns the conductive energy flux across the interface.
Definition: freeflowporousmedium/couplingconditions.hh:268
static NumEqVector< freeFlowMomentumIndex > momentumCouplingCondition(const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const ElementVolumeVariables< freeFlowMomentumIndex > &elemVolVars, const Context &context)
Returns the momentum flux across the coupling boundary.
Definition: freeflowporousmedium/couplingconditions.hh:178
static constexpr auto freeFlowMassIndex
Definition: freeflowporousmedium/couplingconditions.hh:125
static Scalar computeCouplingPhasePressureAtInterface_(const FVElementGeometry< porousMediumIndex > &fvGeometry, const SubControlVolumeFace< porousMediumIndex > &scvf, const VolumeVariables< porousMediumIndex > &volVars, const typename Element< freeFlowMomentumIndex >::Geometry::GlobalCoordinate &couplingPhaseVelocity, ForchheimersLaw)
Returns the pressure at the interface using Forchheimers's law for reconstruction.
Definition: freeflowporousmedium/couplingconditions.hh:312
static Scalar computeCouplingPhasePressureAtInterface_(const FVElementGeometry< porousMediumIndex > &fvGeometry, const SubControlVolumeFace< porousMediumIndex > &scvf, const VolumeVariables< porousMediumIndex > &volVars, const typename Element< freeFlowMomentumIndex >::Geometry::GlobalCoordinate &gravity, const typename Element< freeFlowMomentumIndex >::Geometry::GlobalCoordinate &couplingPhaseVelocity, DarcysLaw)
Returns the pressure at the interface using Darcy's law for reconstruction.
Definition: freeflowporousmedium/couplingconditions.hh:324
Scalar getDistance_(const Scv &scv, const Scvf &scvf) const
Returns the distance between an scvf and the corresponding scv center.
Definition: freeflowporousmedium/couplingconditions.hh:259
static auto darcyPermeability(const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const Context &context)
Returns the intrinsic permeability of the coupled Darcy element.
Definition: freeflowporousmedium/couplingconditions.hh:165
static constexpr auto couplingPhaseIdx(Dune::index_constant< i > id, int coupledPhaseIdx=0)
Returns the corresponding phase index needed for coupling.
Definition: freeflowporousmedium/couplingconditions.hh:151
static constexpr auto couplingCompIdx(Dune::index_constant< i > id, int coupledCompIdx)
Returns the corresponding component index needed for coupling.
Definition: freeflowporousmedium/couplingconditions.hh:158
static constexpr auto porousMediumIndex
Definition: freeflowporousmedium/couplingconditions.hh:126
static Scalar pressureAtInterface_(const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const ElementVolumeVariables< freeFlowMomentumIndex > &elemVolVars, const CouplingContext &context)
Returns the pressure at the interface.
Definition: freeflowporousmedium/couplingconditions.hh:297
Scalar transmissibility_(Dune::index_constant< i > domainI, Dune::index_constant< j > domainJ, const Scalar insideDistance, const Scalar outsideDistance, const Scalar avgQuantityI, const Scalar avgQuantityJ, const DiffusionCoefficientAveragingType diffCoeffAvgType) const
Returns the transmissibility used for either molecular diffusion or thermal conductivity.
Definition: freeflowporousmedium/couplingconditions.hh:225
static Scalar massCouplingCondition(Dune::index_constant< i > domainI, Dune::index_constant< j > domainJ, const FVElementGeometry< i > &fvGeometry, const SubControlVolumeFace< i > &scvf, const ElementVolumeVariables< i > &insideVolVars, const CouplingContext &context)
Returns the mass flux across the coupling boundary as seen from the Darcy domain.
Definition: freeflowporousmedium/couplingconditions.hh:388
Scalar energyCouplingCondition(const Element< ParentType::porousMediumIndex > &element, const FVElementGeometry< ParentType::porousMediumIndex > &fvGeometry, const ElementVolumeVariables< ParentType::porousMediumIndex > &darcyElemVolVars, const SubControlVolumeFace< ParentType::porousMediumIndex > &scvf, const DiffusionCoefficientAveragingType diffCoeffAvgType=DiffusionCoefficientAveragingType::ffOnly) const
Returns the energy flux across the coupling boundary as seen from the Darcy domain.
Definition: freeflowporousmedium/couplingconditions.hh:406
Scalar energyCouplingCondition(const Element< ParentType::freeFlowMassIndex > &element, const FVElementGeometry< ParentType::freeFlowMassIndex > &fvGeometry, const ElementVolumeVariables< ParentType::freeFlowMassIndex > &stokesElemVolVars, const SubControlVolumeFace< ParentType::freeFlowMassIndex > &scvf, const DiffusionCoefficientAveragingType diffCoeffAvgType=DiffusionCoefficientAveragingType::ffOnly) const
Returns the energy flux across the coupling boundary as seen from the free-flow domain.
Definition: freeflowporousmedium/couplingconditions.hh:419
Helper struct to choose the correct index for phases and components. This is need if the porous-mediu...
Definition: indexhelper.hh:42
This structs helps to check if the two sub models use the same fluidsystem. Specialization for the ca...
Definition: multidomain/boundary/freeflowporousmedium/traits.hh:42
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:60
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...
Traits for the free-flow/porous-medium-flow coupling.