25#ifndef DUMUX_STOKES_DARCY_COUPLINGDATA_HH
26#define DUMUX_STOKES_DARCY_COUPLINGDATA_HH
54 harmonic, arithmethic, ffOnly, pmOnly
63 if (diffusionCoefficientAveragingType ==
"Harmonic")
65 else if (diffusionCoefficientAveragingType ==
"Arithmethic")
67 else if (diffusionCoefficientAveragingType ==
"FreeFlowOnly")
69 else if (diffusionCoefficientAveragingType ==
"PorousMediumOnly")
72 DUNE_THROW(Dune::IOError,
"Unknown DiffusionCoefficientAveragingType");
84template<
class FFFS,
class PMFS>
87 static_assert(FFFS::numPhases == 1,
"Only single-phase fluidsystems may be used for free flow.");
88 static constexpr bool value = std::is_same<typename FFFS::MultiPhaseFluidSystem, PMFS>::value;
99 static_assert(FS::numPhases == 1,
"Only single-phase fluidsystems may be used for free flow.");
100 static constexpr bool value = std::is_same<FS, FS>::value;
104template <
class TypeTag, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
112template<
class DiffLaw>
120template<
class T, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
132template<std::
size_t stokesIdx, std::
size_t darcyIdx,
class FFFS,
bool hasAdapter>
143template<std::
size_t stokesIdx, std::
size_t darcyIdx,
class FFFS>
149 template<std::
size_t i>
151 {
return coupledPhaseIdx; }
156 template<std::
size_t i>
158 {
return coupledCompdIdx; }
169template<std::
size_t stokesIdx, std::
size_t darcyIdx,
class FFFS>
175 static constexpr auto couplingPhaseIdx(Dune::index_constant<stokesIdx>,
int coupledPhaseIdx = 0)
181 static constexpr auto couplingPhaseIdx(Dune::index_constant<darcyIdx>,
int coupledPhaseIdx = 0)
182 {
return FFFS::multiphaseFluidsystemPhaseIdx; }
187 static constexpr auto couplingCompIdx(Dune::index_constant<stokesIdx>,
int coupledCompdIdx)
188 {
return coupledCompdIdx; }
193 static constexpr auto couplingCompIdx(Dune::index_constant<darcyIdx>,
int coupledCompdIdx)
194 {
return FFFS::compIdx(coupledCompdIdx); }
198template <
class TypeTag, DiscretizationMethod discMethod>
199class DarcysLawImplementation;
202template <
class TypeTag, DiscretizationMethod discMethod>
203class ForchheimersLawImplementation;
206template<
class MDTraits,
class CouplingManager,
bool enableEnergyBalance,
bool isCompositional>
214template<
class MDTraits,
class CouplingManager>
223template<
class MDTraits,
class CouplingManager>
226 using Scalar =
typename MDTraits::Scalar;
228 template<std::
size_t id>
using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
230 template<std::
size_t id>
using Element =
typename GridGeometry<id>::GridView::template Codim<0>::Entity;
231 template<std::
size_t id>
using FVElementGeometry =
typename GridGeometry<id>::LocalView;
232 template<std::
size_t id>
using SubControlVolumeFace =
typename GridGeometry<id>::LocalView::SubControlVolumeFace;
233 template<std::
size_t id>
using SubControlVolume =
typename GridGeometry<id>::LocalView::SubControlVolume;
240 template<std::
size_t id>
using GlobalPosition =
typename Element<id>::Geometry::GlobalCoordinate;
241 static constexpr auto stokesIdx = CouplingManager::stokesIdx;
242 static constexpr auto darcyIdx = CouplingManager::darcyIdx;
248 static constexpr bool adapterUsed = ModelTraits<darcyIdx>::numFluidPhases() > 1;
253 "All submodels must both be either isothermal or non-isothermal");
256 FluidSystem<darcyIdx>>::value,
257 "All submodels must use the same fluid system");
267 template<std::
size_t i>
268 static constexpr auto couplingPhaseIdx(Dune::index_constant<i>
id,
int coupledPhaseIdx = 0)
269 {
return IndexHelper::couplingPhaseIdx(
id, coupledPhaseIdx); }
274 template<std::
size_t i>
275 static constexpr auto couplingCompIdx(Dune::index_constant<i>
id,
int coupledCompdIdx)
276 {
return IndexHelper::couplingCompIdx(
id, coupledCompdIdx); }
282 {
return couplingManager_; }
287 Scalar
darcyPermeability(
const Element<stokesIdx>& element,
const SubControlVolumeFace<stokesIdx>& scvf)
const
289 const auto& stokesContext =
couplingManager().stokesCouplingContext(element, scvf);
290 return stokesContext.volVars.permeability();
300 template<
class ElementFaceVariables>
302 const FVElementGeometry<stokesIdx>& fvGeometry,
303 const ElementVolumeVariables<stokesIdx>& stokesElemVolVars,
304 const ElementFaceVariables& stokesElemFaceVars,
305 const SubControlVolumeFace<stokesIdx>& scvf)
const
309 Scalar momentumFlux(0.0);
310 const auto& stokesContext = couplingManager_.stokesCouplingContext(element, scvf);
314 const Scalar darcyPressure = stokesContext.volVars.pressure(darcyPhaseIdx);
316 if(numPhasesDarcy > 1)
317 momentumFlux = darcyPressure;
324 momentumFlux -= couplingManager_.
problem(stokesIdx).initial(scvf)[Indices<stokesIdx>::pressureIdx];
326 momentumFlux *= scvf.directionSign();
334 Scalar
advectiveFlux(
const Scalar insideQuantity,
const Scalar outsideQuantity,
const Scalar volumeFlow,
bool insideIsUpstream)
const
336 const Scalar upwindWeight = 1.0;
339 return (upwindWeight * insideQuantity + (1.0 - upwindWeight) * outsideQuantity) * volumeFlow;
341 return (upwindWeight * outsideQuantity + (1.0 - upwindWeight) * insideQuantity) * volumeFlow;
349 template<std::
size_t i, std::
size_t j>
351 Dune::index_constant<j> domainJ,
352 const Scalar insideDistance,
353 const Scalar outsideDistance,
354 const Scalar avgQuantityI,
355 const Scalar avgQuantityJ,
356 const DiffusionCoefficientAveragingType diffCoeffAvgType)
const
358 const Scalar totalDistance = insideDistance + outsideDistance;
359 if(diffCoeffAvgType == DiffusionCoefficientAveragingType::harmonic)
361 return harmonicMean(avgQuantityI, avgQuantityJ, insideDistance, outsideDistance)
364 else if(diffCoeffAvgType == DiffusionCoefficientAveragingType::arithmethic)
366 return arithmeticMean(avgQuantityI, avgQuantityJ, insideDistance, outsideDistance)
369 else if(diffCoeffAvgType == DiffusionCoefficientAveragingType::ffOnly)
370 return domainI == stokesIdx
371 ? avgQuantityI / totalDistance
372 : avgQuantityJ / totalDistance;
375 return domainI == darcyIdx
376 ? avgQuantityI / totalDistance
377 : avgQuantityJ / totalDistance;
383 template<
class Scv,
class Scvf>
386 return (scv.dofPosition() - scvf.ipGlobal()).two_norm();
392 template<std::
size_t i, std::
size_t j,
bool isNI = enableEnergyBalance,
typename std::enable_if_t<isNI,
int> = 0>
394 Dune::index_constant<j> domainJ,
395 const FVElementGeometry<i>& fvGeometryI,
396 const FVElementGeometry<j>& fvGeometryJ,
397 const SubControlVolumeFace<i>& scvfI,
398 const SubControlVolume<i>& scvI,
399 const SubControlVolume<j>& scvJ,
400 const VolumeVariables<i>& volVarsI,
401 const VolumeVariables<j>& volVarsJ,
402 const DiffusionCoefficientAveragingType diffCoeffAvgType)
const
404 const Scalar insideDistance =
getDistance_(scvI, scvfI);
405 const Scalar outsideDistance =
getDistance_(scvJ, scvfI);
407 const Scalar deltaT = volVarsJ.temperature() - volVarsI.temperature();
412 volVarsI.effectiveThermalConductivity(),
413 volVarsJ.effectiveThermalConductivity(),
416 return -tij * deltaT;
422 template<
class ElementFaceVariables,
class CouplingContext>
424 const SubControlVolumeFace<stokesIdx>& scvf,
425 const ElementFaceVariables& elemFaceVars,
426 const CouplingContext& context)
const
428 GlobalPosition<stokesIdx> velocity(0.0);
429 velocity[scvf.directionIndex()] = elemFaceVars[scvf].velocitySelf();
430 const auto& darcyScvf = context.fvGeometry.scvf(context.darcyScvfIdx);
438 const FVElementGeometry<darcyIdx>& fvGeometry,
439 const SubControlVolumeFace<darcyIdx>& scvf,
440 const VolumeVariables<darcyIdx>& volVars,
441 const typename Element<stokesIdx>::Geometry::GlobalCoordinate& couplingPhaseVelocity,
445 const Scalar cellCenterPressure = volVars.pressure(darcyPhaseIdx);
452 const auto velocity = couplingPhaseVelocity;
453 const Scalar mu = volVars.viscosity(darcyPhaseIdx);
454 const Scalar rho = volVars.density(darcyPhaseIdx);
455 const auto K = volVars.permeability();
456 const auto alpha =
vtmv(scvf.unitOuterNormal(), K, couplingManager_.
problem(darcyIdx).spatialParams().gravity(scvf.center()));
458 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
462 Scalar cF = couplingManager_.
problem(darcyIdx).spatialParams().forchCoeff(scvf);
464 const Scalar interfacePressure = ((-mu*(scvf.unitOuterNormal() * velocity))
465 + (-(scvf.unitOuterNormal() * velocity) * velocity.two_norm() * rho * sqrt(
darcyPermeability(element, scvf)) * cF)
466 + rho * alpha)/ti + cellCenterPressure;
467 return interfacePressure;
474 const FVElementGeometry<darcyIdx>& fvGeometry,
475 const SubControlVolumeFace<darcyIdx>& scvf,
476 const VolumeVariables<darcyIdx>& volVars,
477 const typename Element<stokesIdx>::Geometry::GlobalCoordinate& couplingPhaseVelocity,
481 const Scalar couplingPhaseCellCenterPressure = volVars.pressure(darcyPhaseIdx);
482 const Scalar couplingPhaseMobility = volVars.mobility(darcyPhaseIdx);
483 const Scalar couplingPhaseDensity = volVars.density(darcyPhaseIdx);
484 const auto K = volVars.permeability();
490 const auto alpha =
vtmv(scvf.unitOuterNormal(), K, couplingManager_.
problem(darcyIdx).spatialParams().gravity(scvf.center()));
492 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
495 return (-1/couplingPhaseMobility * (scvf.unitOuterNormal() * couplingPhaseVelocity) + couplingPhaseDensity * alpha)/ti
496 + couplingPhaseCellCenterPressure;
509template<
class MDTraits,
class CouplingManager,
bool enableEnergyBalance>
514 using Scalar =
typename MDTraits::Scalar;
515 static constexpr auto stokesIdx = CouplingManager::stokesIdx;
516 static constexpr auto darcyIdx = CouplingManager::darcyIdx;
517 static constexpr auto stokesFaceIdx = CouplingManager::stokesFaceIdx;
518 static constexpr auto stokesCellCenterIdx = CouplingManager::stokesCellCenterIdx;
521 template<std::
size_t id>
522 using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
525 template<std::
size_t id>
using Element =
typename GridGeometry<id>::GridView::template Codim<0>::Entity;
526 template<std::
size_t id>
using FVElementGeometry =
typename GridGeometry<id>::LocalView;
527 template<std::
size_t id>
using SubControlVolumeFace =
typename GridGeometry<id>::LocalView::SubControlVolumeFace;
528 template<std::
size_t id>
using SubControlVolume =
typename GridGeometry<id>::LocalView::SubControlVolume;
535 "Darcy Model must not be compositional");
540 using ParentType::ParentType;
541 using ParentType::couplingPhaseIdx;
547 const FVElementGeometry<darcyIdx>& fvGeometry,
548 const ElementVolumeVariables<darcyIdx>& darcyElemVolVars,
549 const SubControlVolumeFace<darcyIdx>& scvf)
const
551 const auto& darcyContext = this->couplingManager().darcyCouplingContext(element, scvf);
552 const Scalar velocity = darcyContext.velocity * scvf.unitOuterNormal();
553 const Scalar darcyDensity = darcyElemVolVars[scvf.insideScvIdx()].density(couplingPhaseIdx(darcyIdx));
554 const Scalar stokesDensity = darcyContext.volVars.density();
555 const bool insideIsUpstream = velocity > 0.0;
557 return massFlux_(velocity, darcyDensity, stokesDensity, insideIsUpstream);
564 const FVElementGeometry<stokesIdx>& fvGeometry,
565 const ElementVolumeVariables<stokesIdx>& stokesElemVolVars,
566 const ElementFaceVariables<stokesIdx>& stokesElemFaceVars,
567 const SubControlVolumeFace<stokesIdx>& scvf)
const
569 const auto& stokesContext = this->couplingManager().stokesCouplingContext(element, scvf);
570 const Scalar velocity = stokesElemFaceVars[scvf].velocitySelf();
571 const Scalar stokesDensity = stokesElemVolVars[scvf.insideScvIdx()].density();
572 const Scalar darcyDensity = stokesContext.volVars.density(couplingPhaseIdx(darcyIdx));
573 const bool insideIsUpstream =
sign(velocity) == scvf.directionSign();
575 return massFlux_(velocity * scvf.directionSign(), stokesDensity, darcyDensity, insideIsUpstream);
581 template<
bool isNI = enableEnergyBalance,
typename std::enable_if_t<isNI,
int> = 0>
583 const FVElementGeometry<darcyIdx>& fvGeometry,
584 const ElementVolumeVariables<darcyIdx>& darcyElemVolVars,
585 const SubControlVolumeFace<darcyIdx>& scvf,
586 const DiffusionCoefficientAveragingType diffCoeffAvgType = DiffusionCoefficientAveragingType::ffOnly)
const
588 const auto& darcyContext = this->couplingManager().darcyCouplingContext(element, scvf);
589 const auto& darcyVolVars = darcyElemVolVars[scvf.insideScvIdx()];
590 const auto& stokesVolVars = darcyContext.volVars;
592 const Scalar velocity = darcyContext.velocity * scvf.unitOuterNormal();
593 const bool insideIsUpstream = velocity > 0.0;
595 return energyFlux_(darcyIdx, stokesIdx, fvGeometry, darcyContext.fvGeometry, scvf,
596 darcyVolVars, stokesVolVars, velocity, insideIsUpstream, diffCoeffAvgType);
602 template<
bool isNI = enableEnergyBalance,
typename std::enable_if_t<isNI,
int> = 0>
604 const FVElementGeometry<stokesIdx>& fvGeometry,
605 const ElementVolumeVariables<stokesIdx>& stokesElemVolVars,
606 const ElementFaceVariables<stokesIdx>& stokesElemFaceVars,
607 const SubControlVolumeFace<stokesIdx>& scvf,
608 const DiffusionCoefficientAveragingType diffCoeffAvgType = DiffusionCoefficientAveragingType::ffOnly)
const
610 const auto& stokesContext = this->couplingManager().stokesCouplingContext(element, scvf);
611 const auto& stokesVolVars = stokesElemVolVars[scvf.insideScvIdx()];
612 const auto& darcyVolVars = stokesContext.volVars;
614 const Scalar velocity = stokesElemFaceVars[scvf].velocitySelf();
615 const bool insideIsUpstream =
sign(velocity) == scvf.directionSign();
617 return energyFlux_(stokesIdx, darcyIdx, fvGeometry, stokesContext.fvGeometry, scvf,
618 stokesVolVars, darcyVolVars, velocity * scvf.directionSign(), insideIsUpstream, diffCoeffAvgType);
626 Scalar massFlux_(
const Scalar velocity,
627 const Scalar insideDensity,
628 const Scalar outSideDensity,
629 bool insideIsUpstream)
const
631 return this->advectiveFlux(insideDensity, outSideDensity, velocity, insideIsUpstream);
637 template<std::
size_t i, std::
size_t j,
bool isNI = enableEnergyBalance,
typename std::enable_if_t<isNI,
int> = 0>
638 Scalar energyFlux_(Dune::index_constant<i> domainI,
639 Dune::index_constant<j> domainJ,
640 const FVElementGeometry<i>& insideFvGeometry,
641 const FVElementGeometry<j>& outsideFvGeometry,
642 const SubControlVolumeFace<i>& scvf,
643 const VolumeVariables<i>& insideVolVars,
644 const VolumeVariables<j>& outsideVolVars,
645 const Scalar velocity,
646 const bool insideIsUpstream,
647 const DiffusionCoefficientAveragingType diffCoeffAvgType)
const
651 const auto& insideScv = (*scvs(insideFvGeometry).begin());
652 const auto& outsideScv = (*scvs(outsideFvGeometry).begin());
655 const Scalar insideTerm = insideVolVars.density(couplingPhaseIdx(domainI)) * insideVolVars.enthalpy(couplingPhaseIdx(domainI));
656 const Scalar outsideTerm = outsideVolVars.density(couplingPhaseIdx(domainJ)) * outsideVolVars.enthalpy(couplingPhaseIdx(domainJ));
658 flux += this->advectiveFlux(insideTerm, outsideTerm, velocity, insideIsUpstream);
660 flux += this->conductiveEnergyFlux_(domainI, domainJ, insideFvGeometry, outsideFvGeometry, scvf, insideScv, outsideScv, insideVolVars, outsideVolVars, diffCoeffAvgType);
671template<
class MDTraits,
class CouplingManager,
bool enableEnergyBalance>
676 using Scalar =
typename MDTraits::Scalar;
677 static constexpr auto stokesIdx = CouplingManager::stokesIdx;
678 static constexpr auto darcyIdx = CouplingManager::darcyIdx;
679 static constexpr auto stokesFaceIdx = CouplingManager::stokesFaceIdx;
680 static constexpr auto stokesCellCenterIdx = CouplingManager::stokesCellCenterIdx;
683 template<std::
size_t id>
684 using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
687 template<std::
size_t id>
using Element =
typename GridGeometry<id>::GridView::template Codim<0>::Entity;
688 template<std::
size_t id>
using FVElementGeometry =
typename GridGeometry<id>::LocalView;
689 template<std::
size_t id>
using SubControlVolumeFace =
typename FVElementGeometry<id>::SubControlVolumeFace;
690 template<std::
size_t id>
using SubControlVolume =
typename GridGeometry<id>::LocalView::SubControlVolume;
703 static_assert(getPropValue<SubDomainTypeTag<darcyIdx>,
Properties::UseMoles>() == useMoles,
"Both submodels must either use moles or not");
704 static_assert(getPropValue<SubDomainTypeTag<darcyIdx>,
Properties::ReplaceCompEqIdx>() == replaceCompEqIdx,
"Both submodels must use the same replaceCompEqIdx");
706 "Both submodels must use the same reference system formulation for diffusion");
708 using NumEqVector = Dune::FieldVector<Scalar, numComponents>;
714 "Both submodels must use the same diffusion law.");
716 using ReducedComponentVector = Dune::FieldVector<Scalar, numComponents-1>;
717 using ReducedComponentMatrix = Dune::FieldMatrix<Scalar, numComponents-1, numComponents-1>;
721 using ParentType::ParentType;
722 using ParentType::couplingPhaseIdx;
723 using ParentType::couplingCompIdx;
729 const FVElementGeometry<darcyIdx>& fvGeometry,
730 const ElementVolumeVariables<darcyIdx>& darcyElemVolVars,
731 const SubControlVolumeFace<darcyIdx>& scvf,
732 const DiffusionCoefficientAveragingType diffCoeffAvgType = DiffusionCoefficientAveragingType::ffOnly)
const
734 NumEqVector flux(0.0);
735 const auto& darcyContext = this->couplingManager().darcyCouplingContext(element, scvf);
736 const auto& darcyVolVars = darcyElemVolVars[scvf.insideScvIdx()];
737 const auto& stokesVolVars = darcyContext.volVars;
738 const auto& outsideScv = (*scvs(darcyContext.fvGeometry).begin());
740 const Scalar velocity = darcyContext.velocity * scvf.unitOuterNormal();
741 const bool insideIsUpstream = velocity > 0.0;
743 return massFlux_(darcyIdx, stokesIdx, fvGeometry,
744 scvf, darcyVolVars, stokesVolVars,
745 outsideScv, velocity, insideIsUpstream,
753 const FVElementGeometry<stokesIdx>& fvGeometry,
754 const ElementVolumeVariables<stokesIdx>& stokesElemVolVars,
755 const ElementFaceVariables<stokesIdx>& stokesElemFaceVars,
756 const SubControlVolumeFace<stokesIdx>& scvf,
757 const DiffusionCoefficientAveragingType diffCoeffAvgType = DiffusionCoefficientAveragingType::ffOnly)
const
759 NumEqVector flux(0.0);
760 const auto& stokesContext = this->couplingManager().stokesCouplingContext(element, scvf);
761 const auto& stokesVolVars = stokesElemVolVars[scvf.insideScvIdx()];
762 const auto& darcyVolVars = stokesContext.volVars;
763 const auto& outsideScv = (*scvs(stokesContext.fvGeometry).begin());
765 const Scalar velocity = stokesElemFaceVars[scvf].velocitySelf();
766 const bool insideIsUpstream =
sign(velocity) == scvf.directionSign();
768 return massFlux_(stokesIdx, darcyIdx, fvGeometry,
769 scvf, stokesVolVars, darcyVolVars,
770 outsideScv, velocity * scvf.directionSign(),
771 insideIsUpstream, diffCoeffAvgType);
777 template<
bool isNI = enableEnergyBalance,
typename std::enable_if_t<isNI,
int> = 0>
779 const FVElementGeometry<darcyIdx>& fvGeometry,
780 const ElementVolumeVariables<darcyIdx>& darcyElemVolVars,
781 const SubControlVolumeFace<darcyIdx>& scvf,
782 const DiffusionCoefficientAveragingType diffCoeffAvgType = DiffusionCoefficientAveragingType::ffOnly)
const
784 const auto& darcyContext = this->couplingManager().darcyCouplingContext(element, scvf);
785 const auto& darcyVolVars = darcyElemVolVars[scvf.insideScvIdx()];
786 const auto& stokesVolVars = darcyContext.volVars;
788 const Scalar velocity = darcyContext.velocity * scvf.unitOuterNormal();
789 const bool insideIsUpstream = velocity > 0.0;
791 return energyFlux_(darcyIdx, stokesIdx, fvGeometry, darcyContext.fvGeometry, scvf,
792 darcyVolVars, stokesVolVars, velocity, insideIsUpstream, diffCoeffAvgType);
798 template<
bool isNI = enableEnergyBalance,
typename std::enable_if_t<isNI,
int> = 0>
800 const FVElementGeometry<stokesIdx>& fvGeometry,
801 const ElementVolumeVariables<stokesIdx>& stokesElemVolVars,
802 const ElementFaceVariables<stokesIdx>& stokesElemFaceVars,
803 const SubControlVolumeFace<stokesIdx>& scvf,
804 const DiffusionCoefficientAveragingType diffCoeffAvgType = DiffusionCoefficientAveragingType::ffOnly)
const
806 const auto& stokesContext = this->couplingManager().stokesCouplingContext(element, scvf);
807 const auto& stokesVolVars = stokesElemVolVars[scvf.insideScvIdx()];
808 const auto& darcyVolVars = stokesContext.volVars;
810 const Scalar velocity = stokesElemFaceVars[scvf].velocitySelf();
811 const bool insideIsUpstream =
sign(velocity) == scvf.directionSign();
813 return energyFlux_(stokesIdx, darcyIdx, fvGeometry, stokesContext.fvGeometry, scvf,
814 stokesVolVars, darcyVolVars, velocity * scvf.directionSign(), insideIsUpstream, diffCoeffAvgType);
822 template<std::
size_t i, std::
size_t j>
824 Dune::index_constant<j> domainJ,
825 const FVElementGeometry<i>& insideFvGeometry,
826 const SubControlVolumeFace<i>& scvf,
827 const VolumeVariables<i>& insideVolVars,
828 const VolumeVariables<j>& outsideVolVars,
829 const SubControlVolume<j>& outsideScv,
830 const Scalar velocity,
831 const bool insideIsUpstream,
832 const DiffusionCoefficientAveragingType diffCoeffAvgType)
const
834 NumEqVector flux(0.0);
835 NumEqVector diffusiveFlux(0.0);
837 auto moleOrMassFraction = [&](
const auto& volVars,
int phaseIdx,
int compIdx)
838 {
return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); };
840 auto moleOrMassDensity = [&](
const auto& volVars,
int phaseIdx)
841 {
return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); };
844 auto insideTerm = [&](
int compIdx)
845 {
return moleOrMassFraction(insideVolVars, couplingPhaseIdx(domainI), compIdx) * moleOrMassDensity(insideVolVars, couplingPhaseIdx(domainI)); };
847 auto outsideTerm = [&](
int compIdx)
848 {
return moleOrMassFraction(outsideVolVars, couplingPhaseIdx(domainJ), compIdx) * moleOrMassDensity(outsideVolVars, couplingPhaseIdx(domainJ)); };
851 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
853 const int domainICompIdx = couplingCompIdx(domainI, compIdx);
854 const int domainJCompIdx = couplingCompIdx(domainJ, compIdx);
855 flux[domainICompIdx] += this->advectiveFlux(insideTerm(domainICompIdx), outsideTerm(domainJCompIdx), velocity, insideIsUpstream);
859 const auto& insideScv = insideFvGeometry.scv(scvf.insideScvIdx());
862 diffusiveFlux += diffusiveMolecularFluxFicksLaw_(domainI, domainJ, scvf, insideScv, outsideScv, insideVolVars, outsideVolVars, diffCoeffAvgType);
864 diffusiveFlux += diffusiveMolecularFluxMaxwellStefan_(domainI, domainJ, scvf, insideScv, outsideScv, insideVolVars, outsideVolVars);
869 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
871 const int domainICompIdx = couplingCompIdx(domainI, compIdx);
872 diffusiveFlux[domainICompIdx] *= 1/FluidSystem<i>::molarMass(domainICompIdx);
877 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
879 const int domainICompIdx = couplingCompIdx(domainI, compIdx);
880 diffusiveFlux[domainICompIdx] *= FluidSystem<i>::molarMass(domainICompIdx);
884 flux += diffusiveFlux;
886 if (replaceCompEqIdx < numComponents)
887 flux[replaceCompEqIdx] = std::accumulate(flux.begin(), flux.end(), 0.0);
894 return FluidSystem<stokesIdx>::componentEnthalpy(volVars.fluidState(), 0, compIdx);
899 return FluidSystem<darcyIdx>::componentEnthalpy(volVars.fluidState(), phaseIdx, compIdx);
905 template<std::
size_t i, std::
size_t j>
907 Dune::index_constant<j> domainJ,
908 const SubControlVolumeFace<i>& scvfI,
909 const SubControlVolume<i>& scvI,
910 const SubControlVolume<j>& scvJ,
911 const VolumeVariables<i>& volVarsI,
912 const VolumeVariables<j>& volVarsJ)
const
915 NumEqVector diffusiveFlux(0.0);
917 const Scalar insideDistance = this->getDistance_(scvI, scvfI);
918 const Scalar outsideDistance = this->getDistance_(scvJ, scvfI);
920 ReducedComponentVector moleFracInside(0.0);
921 ReducedComponentVector moleFracOutside(0.0);
922 ReducedComponentVector reducedFlux(0.0);
923 ReducedComponentMatrix reducedDiffusionMatrixInside(0.0);
924 ReducedComponentMatrix reducedDiffusionMatrixOutside(0.0);
927 for (
int compIdx = 0; compIdx < numComponents-1; compIdx++)
929 const int domainICompIdx = couplingCompIdx(domainI, compIdx);
930 const int domainJCompIdx = couplingCompIdx(domainJ, compIdx);
932 assert(FluidSystem<i>::componentName(domainICompIdx) == FluidSystem<j>::componentName(domainJCompIdx));
935 const Scalar xInside = volVarsI.moleFraction(couplingPhaseIdx(domainI), domainICompIdx);
937 const Scalar xOutside = volVarsJ.moleFraction(couplingPhaseIdx(domainJ), domainJCompIdx);
938 moleFracInside[domainICompIdx] = xInside;
939 moleFracOutside[domainICompIdx] = xOutside;
947 for (
int compKIdx = 0; compKIdx < numComponents-1; compKIdx++)
949 const int domainICompKIdx = couplingCompIdx(domainI, compKIdx);
950 const Scalar xk = volVarsI.moleFraction(couplingPhaseIdx(domainI), domainICompKIdx);
951 const Scalar avgMolarMass = volVarsI.averageMolarMass(couplingPhaseIdx(domainI));
952 const Scalar Mn = FluidSystem<i>::molarMass(numComponents-1);
953 const Scalar tkn = Deprecated::template effectiveDiffusionCoefficient<EffDiffModel>(volVarsI, couplingPhaseIdx(domainI), domainICompKIdx, couplingCompIdx(domainI, numComponents-1));
956 reducedDiffusionMatrixInside[domainICompKIdx][domainICompKIdx] += xk*avgMolarMass/(tkn*Mn);
958 for (
int compLIdx = 0; compLIdx < numComponents; compLIdx++)
960 const int domainICompLIdx = couplingCompIdx(domainI, compLIdx);
963 if (domainICompKIdx == domainICompLIdx)
966 const Scalar xl = volVarsI.moleFraction(couplingPhaseIdx(domainI), domainICompLIdx);
967 const Scalar Mk = FluidSystem<i>::molarMass(domainICompKIdx);
968 const Scalar Ml = FluidSystem<i>::molarMass(domainICompLIdx);
969 const Scalar tkl = Deprecated::template effectiveDiffusionCoefficient<EffDiffModel>(volVarsI, couplingPhaseIdx(domainI), domainICompKIdx, domainICompLIdx);
970 reducedDiffusionMatrixInside[domainICompKIdx][domainICompKIdx] += xl*avgMolarMass/(tkl*Mk);
971 reducedDiffusionMatrixInside[domainICompKIdx][domainICompLIdx] += xk*(avgMolarMass/(tkn*Mn) - avgMolarMass/(tkl*Ml));
975 for (
int compKIdx = 0; compKIdx < numComponents-1; compKIdx++)
977 const int domainJCompKIdx = couplingCompIdx(domainJ, compKIdx);
978 const int domainICompKIdx = couplingCompIdx(domainI, compKIdx);
980 const Scalar xk = volVarsJ.moleFraction(couplingPhaseIdx(domainJ), domainJCompKIdx);
981 const Scalar avgMolarMass = volVarsJ.averageMolarMass(couplingPhaseIdx(domainJ));
982 const Scalar Mn = FluidSystem<j>::molarMass(numComponents-1);
983 const Scalar tkn = Deprecated::template effectiveDiffusionCoefficient<EffDiffModel>(volVarsJ, couplingPhaseIdx(domainJ), domainJCompKIdx, couplingCompIdx(domainJ, numComponents-1));
986 reducedDiffusionMatrixOutside[domainICompKIdx][domainICompKIdx] += xk*avgMolarMass/(tkn*Mn);
988 for (
int compLIdx = 0; compLIdx < numComponents; compLIdx++)
990 const int domainJCompLIdx = couplingCompIdx(domainJ, compLIdx);
991 const int domainICompLIdx = couplingCompIdx(domainI, compLIdx);
994 if (domainJCompLIdx == domainJCompKIdx)
997 const Scalar xl = volVarsJ.moleFraction(couplingPhaseIdx(domainJ), domainJCompLIdx);
998 const Scalar Mk = FluidSystem<j>::molarMass(domainJCompKIdx);
999 const Scalar Ml = FluidSystem<j>::molarMass(domainJCompLIdx);
1000 const Scalar tkl = Deprecated::template effectiveDiffusionCoefficient<EffDiffModel>(volVarsJ, couplingPhaseIdx(domainJ), domainJCompKIdx, domainJCompLIdx);
1001 reducedDiffusionMatrixOutside[domainICompKIdx][domainICompKIdx] += xl*avgMolarMass/(tkl*Mk);
1002 reducedDiffusionMatrixOutside[domainICompKIdx][domainICompLIdx] += xk*(avgMolarMass/(tkn*Mn) - avgMolarMass/(tkl*Ml));
1006 const Scalar omegai = 1/insideDistance;
1007 const Scalar omegaj = 1/outsideDistance;
1009 reducedDiffusionMatrixInside.invert();
1010 reducedDiffusionMatrixInside *= omegai*volVarsI.density(couplingPhaseIdx(domainI));
1011 reducedDiffusionMatrixOutside.invert();
1012 reducedDiffusionMatrixOutside *= omegaj*volVarsJ.density(couplingPhaseIdx(domainJ));
1015 ReducedComponentVector helperVector(0.0);
1016 ReducedComponentVector gradientVectori(0.0);
1017 ReducedComponentVector gradientVectorj(0.0);
1019 reducedDiffusionMatrixInside.mv(moleFracInside, gradientVectori);
1020 reducedDiffusionMatrixOutside.mv(moleFracOutside, gradientVectorj);
1022 auto gradientVectorij = (gradientVectori + gradientVectorj);
1025 reducedDiffusionMatrixOutside += reducedDiffusionMatrixInside;
1027 reducedDiffusionMatrixOutside.solve(helperVector, gradientVectorij);
1030 helperVector -=moleFracInside;
1031 reducedDiffusionMatrixInside.mv(helperVector, reducedFlux);
1035 for (
int compIdx = 0; compIdx < numComponents-1; compIdx++)
1037 const int domainICompIdx = couplingCompIdx(domainI, compIdx);
1038 diffusiveFlux[domainICompIdx] = reducedFlux[domainICompIdx];
1039 diffusiveFlux[couplingCompIdx(domainI, numComponents-1)] -= reducedFlux[domainICompIdx];
1041 return diffusiveFlux;
1044 template<std::
size_t i, std::
size_t j>
1046 Dune::index_constant<j> domainJ,
1047 const SubControlVolumeFace<i>& scvfI,
1048 const SubControlVolume<i>& scvI,
1049 const SubControlVolume<j>& scvJ,
1050 const VolumeVariables<i>& volVarsI,
1051 const VolumeVariables<j>& volVarsJ,
1052 const DiffusionCoefficientAveragingType diffCoeffAvgType)
const
1055 NumEqVector diffusiveFlux(0.0);
1057 const Scalar rhoInside =
massOrMolarDensity(volVarsI, referenceSystemFormulation, couplingPhaseIdx(domainI));
1058 const Scalar rhoOutside =
massOrMolarDensity(volVarsJ, referenceSystemFormulation, couplingPhaseIdx(domainJ));
1059 const Scalar avgDensity = 0.5 * rhoInside + 0.5 * rhoOutside;
1061 const Scalar insideDistance = this->getDistance_(scvI, scvfI);
1062 const Scalar outsideDistance = this->getDistance_(scvJ, scvfI);
1064 for (
int compIdx = 1; compIdx < numComponents; ++compIdx)
1066 const int domainIMainCompIdx = couplingPhaseIdx(domainI);
1067 const int domainJMainCompIdx = couplingPhaseIdx(domainJ);
1068 const int domainICompIdx = couplingCompIdx(domainI, compIdx);
1069 const int domainJCompIdx = couplingCompIdx(domainJ, compIdx);
1071 assert(FluidSystem<i>::componentName(domainICompIdx) == FluidSystem<j>::componentName(domainJCompIdx));
1073 const Scalar massOrMoleFractionInside =
massOrMoleFraction(volVarsI, referenceSystemFormulation, couplingPhaseIdx(domainI), domainICompIdx);
1074 const Scalar massOrMoleFractionOutside =
massOrMoleFraction(volVarsJ, referenceSystemFormulation, couplingPhaseIdx(domainJ), domainJCompIdx);
1076 const Scalar deltaMassOrMoleFrac = massOrMoleFractionOutside - massOrMoleFractionInside;
1077 const Scalar tij = this->transmissibility_(domainI,
1081 Deprecated::template effectiveDiffusionCoefficient<EffDiffModel>(volVarsI, couplingPhaseIdx(domainI), domainIMainCompIdx, domainICompIdx),
1082 Deprecated::template effectiveDiffusionCoefficient<EffDiffModel>(volVarsJ, couplingPhaseIdx(domainJ), domainJMainCompIdx, domainJCompIdx),
1084 diffusiveFlux[domainICompIdx] += -avgDensity * tij * deltaMassOrMoleFrac;
1087 const Scalar cumulativeFlux = std::accumulate(diffusiveFlux.begin(), diffusiveFlux.end(), 0.0);
1088 diffusiveFlux[couplingCompIdx(domainI, 0)] = -cumulativeFlux;
1090 return diffusiveFlux;
1096 template<std::
size_t i, std::
size_t j,
bool isNI = enableEnergyBalance,
typename std::enable_if_t<isNI,
int> = 0>
1098 Dune::index_constant<j> domainJ,
1099 const FVElementGeometry<i>& insideFvGeometry,
1100 const FVElementGeometry<j>& outsideFvGeometry,
1101 const SubControlVolumeFace<i>& scvf,
1102 const VolumeVariables<i>& insideVolVars,
1103 const VolumeVariables<j>& outsideVolVars,
1104 const Scalar velocity,
1105 const bool insideIsUpstream,
1106 const DiffusionCoefficientAveragingType diffCoeffAvgType)
const
1110 const auto& insideScv = (*scvs(insideFvGeometry).begin());
1111 const auto& outsideScv = (*scvs(outsideFvGeometry).begin());
1114 const Scalar insideTerm = insideVolVars.density(couplingPhaseIdx(domainI)) * insideVolVars.enthalpy(couplingPhaseIdx(domainI));
1115 const Scalar outsideTerm = outsideVolVars.density(couplingPhaseIdx(domainJ)) * outsideVolVars.enthalpy(couplingPhaseIdx(domainJ));
1117 flux += this->advectiveFlux(insideTerm, outsideTerm, velocity, insideIsUpstream);
1119 flux += this->conductiveEnergyFlux_(domainI, domainJ, insideFvGeometry, outsideFvGeometry, scvf, insideScv, outsideScv, insideVolVars, outsideVolVars, diffCoeffAvgType);
1121 auto diffusiveFlux = isFicksLaw ? diffusiveMolecularFluxFicksLaw_(domainI, domainJ, scvf, insideScv, outsideScv, insideVolVars, outsideVolVars, diffCoeffAvgType)
1122 : diffusiveMolecularFluxMaxwellStefan_(domainI, domainJ, scvf, insideScv, outsideScv, insideVolVars, outsideVolVars);
1125 for (
int compIdx = 0; compIdx < diffusiveFlux.size(); ++compIdx)
1127 const int domainICompIdx = couplingCompIdx(domainI, compIdx);
1128 const int domainJCompIdx = couplingCompIdx(domainJ, compIdx);
1130 const bool insideDiffFluxIsUpstream = diffusiveFlux[domainICompIdx] > 0;
1131 const Scalar componentEnthalpy = insideDiffFluxIsUpstream ?
1132 getComponentEnthalpy(insideVolVars, couplingPhaseIdx(domainI), domainICompIdx)
1133 : getComponentEnthalpy(outsideVolVars, couplingPhaseIdx(domainJ), domainJCompIdx);
1136 flux += diffusiveFlux[domainICompIdx] * componentEnthalpy;
1138 flux += diffusiveFlux[domainICompIdx] * FluidSystem<i>::molarMass(domainICompIdx) * componentEnthalpy;
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
Define some often used mathematical functions.
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:47
VolumeVariables::PrimaryVariables::value_type massOrMoleFraction(const VolumeVariables &volVars, ReferenceSystemFormulation referenceSys, const int phaseIdx, const int compIdx)
returns the mass or mole fraction to be used in Fick's law based on the reference system
Definition: referencesystemformulation.hh:66
VolumeVariables::PrimaryVariables::value_type massOrMolarDensity(const VolumeVariables &volVars, ReferenceSystemFormulation referenceSys, const int phaseIdx)
evaluates the density to be used in Fick's law based on the reference system
Definition: referencesystemformulation.hh:55
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:68
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
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:618
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:49
constexpr auto getPropValue()
get the value data member of a property
Definition: propertysystem.hh:153
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
Traits class encapsulating model specifications.
Definition: common/properties.hh:66
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:70
Property whether to use moles or kg as amount unit for balance equations.
Definition: common/properties.hh:100
The component balance index that should be replaced by the total mass/mole balance.
Definition: common/properties.hh:102
Definition: common/properties.hh:113
The type for a global container for the volume variables.
Definition: common/properties.hh:120
The type for the calculation the advective fluxes.
Definition: common/properties.hh:152
The type for the calculation of the molecular diffusion fluxes.
Definition: common/properties.hh:156
The type of the fluid system to use.
Definition: common/properties.hh:167
The employed model for the computation of the effective diffusivity.
Definition: common/properties.hh:175
Global vector containing face-related data.
Definition: common/properties.hh:229
Returns whether to normalize the pressure term in the momentum balance or not.
Definition: common/properties.hh:290
forward declaration of the method-specific implementation
Definition: flux/darcyslaw.hh:38
forward declaration of the method-specific implemetation
Definition: flux/box/fickslaw.hh:43
forward declare
Definition: forchheimerslaw.hh:38
This structs holds a set of options which allow to modify the Stokes-Darcy coupling mechanism during ...
Definition: couplingdata.hh:46
DiffusionCoefficientAveragingType
Defines which kind of averanging of diffusion coefficiencients (moleculat diffusion or thermal conduc...
Definition: couplingdata.hh:53
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: couplingdata.hh:61
This structs helps to check if the two sub models use the same fluidsystem. Specialization for the ca...
Definition: couplingdata.hh:86
static constexpr bool value
Definition: couplingdata.hh:88
This structs indicates that Fick's law is not used for diffusion.
Definition: couplingdata.hh:113
Helper struct to choose the correct index for phases and components. This is need if the porous-mediu...
Definition: couplingdata.hh:133
static constexpr auto couplingCompIdx(Dune::index_constant< i >, int coupledCompdIdx)
No adapter is used, just return the input index.
Definition: couplingdata.hh:157
static constexpr auto couplingPhaseIdx(Dune::index_constant< i >, int coupledPhaseIdx=0)
No adapter is used, just return the input index.
Definition: couplingdata.hh:150
static constexpr auto couplingPhaseIdx(Dune::index_constant< darcyIdx >, int coupledPhaseIdx=0)
The phase index of the porous-medium-flow model is given by the adapter fluidsytem (i....
Definition: couplingdata.hh:181
static constexpr auto couplingCompIdx(Dune::index_constant< stokesIdx >, int coupledCompdIdx)
The free-flow model does not need any change of the component index.
Definition: couplingdata.hh:187
static constexpr auto couplingCompIdx(Dune::index_constant< darcyIdx >, int coupledCompdIdx)
The component index of the porous-medium-flow model is mapped by the adapter fluidsytem.
Definition: couplingdata.hh:193
static constexpr auto couplingPhaseIdx(Dune::index_constant< stokesIdx >, int coupledPhaseIdx=0)
The free-flow model always uses phase index 0.
Definition: couplingdata.hh:175
Definition: couplingdata.hh:207
A base class which provides some common methods used for Stokes-Darcy coupling.
Definition: couplingdata.hh:225
static constexpr auto couplingPhaseIdx(Dune::index_constant< i > id, int coupledPhaseIdx=0)
Returns the corresponding phase index needed for coupling.
Definition: couplingdata.hh:268
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 acorss the interface.
Definition: couplingdata.hh:393
static constexpr auto couplingCompIdx(Dune::index_constant< i > id, int coupledCompdIdx)
Returns the corresponding component index needed for coupling.
Definition: couplingdata.hh:275
Scalar pressureAtInterface_(const Element< stokesIdx > &element, const SubControlVolumeFace< stokesIdx > &scvf, const ElementFaceVariables &elemFaceVars, const CouplingContext &context) const
Returns the pressure at the interface.
Definition: couplingdata.hh:423
Scalar advectiveFlux(const Scalar insideQuantity, const Scalar outsideQuantity, const Scalar volumeFlow, bool insideIsUpstream) const
Evaluate an advective flux across the interface and consider upwinding.
Definition: couplingdata.hh:334
Scalar computeCouplingPhasePressureAtInterface_(const Element< darcyIdx > &element, const FVElementGeometry< darcyIdx > &fvGeometry, const SubControlVolumeFace< darcyIdx > &scvf, const VolumeVariables< darcyIdx > &volVars, const typename Element< stokesIdx >::Geometry::GlobalCoordinate &couplingPhaseVelocity, DarcysLaw) const
Returns the pressure at the interface using Darcy's law for reconstruction.
Definition: couplingdata.hh:473
Scalar getDistance_(const Scv &scv, const Scvf &scvf) const
Returns the distance between an scvf and the corresponding scv center.
Definition: couplingdata.hh:384
Scalar darcyPermeability(const Element< stokesIdx > &element, const SubControlVolumeFace< stokesIdx > &scvf) const
Returns the intrinsic permeability of the coupled Darcy element.
Definition: couplingdata.hh:287
Scalar momentumCouplingCondition(const Element< stokesIdx > &element, const FVElementGeometry< stokesIdx > &fvGeometry, const ElementVolumeVariables< stokesIdx > &stokesElemVolVars, const ElementFaceVariables &stokesElemFaceVars, const SubControlVolumeFace< stokesIdx > &scvf) const
Returns the momentum flux across the coupling boundary.
Definition: couplingdata.hh:301
Scalar computeCouplingPhasePressureAtInterface_(const Element< darcyIdx > &element, const FVElementGeometry< darcyIdx > &fvGeometry, const SubControlVolumeFace< darcyIdx > &scvf, const VolumeVariables< darcyIdx > &volVars, const typename Element< stokesIdx >::Geometry::GlobalCoordinate &couplingPhaseVelocity, ForchheimersLaw) const
Returns the pressure at the interface using Forchheimers's law for reconstruction.
Definition: couplingdata.hh:437
const CouplingManager & couplingManager() const
Returns a reference to the coupling manager.
Definition: couplingdata.hh:281
StokesDarcyCouplingDataImplementationBase(const CouplingManager &couplingmanager)
Definition: couplingdata.hh:262
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: couplingdata.hh:350
Scalar energyCouplingCondition(const Element< stokesIdx > &element, const FVElementGeometry< stokesIdx > &fvGeometry, const ElementVolumeVariables< stokesIdx > &stokesElemVolVars, const ElementFaceVariables< stokesIdx > &stokesElemFaceVars, const SubControlVolumeFace< stokesIdx > &scvf, const DiffusionCoefficientAveragingType diffCoeffAvgType=DiffusionCoefficientAveragingType::ffOnly) const
Returns the energy flux across the coupling boundary as seen from the free-flow domain.
Definition: couplingdata.hh:603
Scalar massCouplingCondition(const Element< stokesIdx > &element, const FVElementGeometry< stokesIdx > &fvGeometry, const ElementVolumeVariables< stokesIdx > &stokesElemVolVars, const ElementFaceVariables< stokesIdx > &stokesElemFaceVars, const SubControlVolumeFace< stokesIdx > &scvf) const
Returns the mass flux across the coupling boundary as seen from the free-flow domain.
Definition: couplingdata.hh:563
Scalar energyCouplingCondition(const Element< darcyIdx > &element, const FVElementGeometry< darcyIdx > &fvGeometry, const ElementVolumeVariables< darcyIdx > &darcyElemVolVars, const SubControlVolumeFace< darcyIdx > &scvf, const DiffusionCoefficientAveragingType diffCoeffAvgType=DiffusionCoefficientAveragingType::ffOnly) const
Returns the energy flux across the coupling boundary as seen from the Darcy domain.
Definition: couplingdata.hh:582
Scalar massCouplingCondition(const Element< darcyIdx > &element, const FVElementGeometry< darcyIdx > &fvGeometry, const ElementVolumeVariables< darcyIdx > &darcyElemVolVars, const SubControlVolumeFace< darcyIdx > &scvf) const
Returns the mass flux across the coupling boundary as seen from the Darcy domain.
Definition: couplingdata.hh:546
NumEqVector massCouplingCondition(const Element< stokesIdx > &element, const FVElementGeometry< stokesIdx > &fvGeometry, const ElementVolumeVariables< stokesIdx > &stokesElemVolVars, const ElementFaceVariables< stokesIdx > &stokesElemFaceVars, const SubControlVolumeFace< stokesIdx > &scvf, const DiffusionCoefficientAveragingType diffCoeffAvgType=DiffusionCoefficientAveragingType::ffOnly) const
Returns the mass flux across the coupling boundary as seen from the free-flow domain.
Definition: couplingdata.hh:752
Scalar energyCouplingCondition(const Element< darcyIdx > &element, const FVElementGeometry< darcyIdx > &fvGeometry, const ElementVolumeVariables< darcyIdx > &darcyElemVolVars, const SubControlVolumeFace< darcyIdx > &scvf, const DiffusionCoefficientAveragingType diffCoeffAvgType=DiffusionCoefficientAveragingType::ffOnly) const
Returns the energy flux across the coupling boundary as seen from the Darcy domain.
Definition: couplingdata.hh:778
Scalar getComponentEnthalpy(const VolumeVariables< darcyIdx > &volVars, int phaseIdx, int compIdx) const
Definition: couplingdata.hh:897
Scalar getComponentEnthalpy(const VolumeVariables< stokesIdx > &volVars, int phaseIdx, int compIdx) const
Definition: couplingdata.hh:892
NumEqVector massFlux_(Dune::index_constant< i > domainI, Dune::index_constant< j > domainJ, const FVElementGeometry< i > &insideFvGeometry, const SubControlVolumeFace< i > &scvf, const VolumeVariables< i > &insideVolVars, const VolumeVariables< j > &outsideVolVars, const SubControlVolume< j > &outsideScv, const Scalar velocity, const bool insideIsUpstream, const DiffusionCoefficientAveragingType diffCoeffAvgType) const
Evaluate the compositional mole/mass flux across the interface.
Definition: couplingdata.hh:823
Scalar energyFlux_(Dune::index_constant< i > domainI, Dune::index_constant< j > domainJ, const FVElementGeometry< i > &insideFvGeometry, const FVElementGeometry< j > &outsideFvGeometry, const SubControlVolumeFace< i > &scvf, const VolumeVariables< i > &insideVolVars, const VolumeVariables< j > &outsideVolVars, const Scalar velocity, const bool insideIsUpstream, const DiffusionCoefficientAveragingType diffCoeffAvgType) const
Evaluate the energy flux across the interface.
Definition: couplingdata.hh:1097
NumEqVector massCouplingCondition(const Element< darcyIdx > &element, const FVElementGeometry< darcyIdx > &fvGeometry, const ElementVolumeVariables< darcyIdx > &darcyElemVolVars, const SubControlVolumeFace< darcyIdx > &scvf, const DiffusionCoefficientAveragingType diffCoeffAvgType=DiffusionCoefficientAveragingType::ffOnly) const
Returns the mass flux across the coupling boundary as seen from the Darcy domain.
Definition: couplingdata.hh:728
Scalar energyCouplingCondition(const Element< stokesIdx > &element, const FVElementGeometry< stokesIdx > &fvGeometry, const ElementVolumeVariables< stokesIdx > &stokesElemVolVars, const ElementFaceVariables< stokesIdx > &stokesElemFaceVars, const SubControlVolumeFace< stokesIdx > &scvf, const DiffusionCoefficientAveragingType diffCoeffAvgType=DiffusionCoefficientAveragingType::ffOnly) const
Returns the energy flux across the coupling boundary as seen from the free-flow domain.
Definition: couplingdata.hh:799
NumEqVector diffusiveMolecularFluxMaxwellStefan_(Dune::index_constant< i > domainI, Dune::index_constant< j > domainJ, const SubControlVolumeFace< i > &scvfI, const SubControlVolume< i > &scvI, const SubControlVolume< j > &scvJ, const VolumeVariables< i > &volVarsI, const VolumeVariables< j > &volVarsJ) const
Evaluate the diffusive mole/mass flux across the interface.
Definition: couplingdata.hh:906
NumEqVector diffusiveMolecularFluxFicksLaw_(Dune::index_constant< i > domainI, Dune::index_constant< j > domainJ, const SubControlVolumeFace< i > &scvfI, const SubControlVolume< i > &scvI, const SubControlVolume< j > &scvJ, const VolumeVariables< i > &volVarsI, const VolumeVariables< j > &volVarsJ, const DiffusionCoefficientAveragingType diffCoeffAvgType) const
Definition: couplingdata.hh:1045
Definition: multidomain/couplingmanager.hh:46
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:264
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...
The interface of the coupling manager for multi domain problems.
Declares all properties used in Dumux.