12#ifndef DUMUX_DISCRETIZATION_BOX_FACET_COUPLING_FICKS_LAW_HH
13#define DUMUX_DISCRETIZATION_BOX_FACET_COUPLING_FICKS_LAW_HH
18#include <dune/common/exceptions.hh>
19#include <dune/common/fvector.hh>
20#include <dune/common/float_cmp.hh>
38template<
class TypeTag, ReferenceSystemFormulation referenceSystem = ReferenceSystemFormulation::massAveraged>
45 using FVElementGeometry =
typename GridGeometry::LocalView;
46 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
48 using GridView =
typename GridGeometry::GridView;
49 using Element =
typename GridView::template Codim<0>::Entity;
51 static constexpr int dim = GridView::dimension;
52 static constexpr int dimWorld = GridView::dimensionworld;
57 static constexpr int numComponents = ModelTraits::numFluidComponents();
60 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
64 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache>
65 static ComponentFluxVector
flux(
const Problem& problem,
66 const Element& element,
67 const FVElementGeometry& fvGeometry,
68 const ElementVolumeVariables& elemVolVars,
69 const SubControlVolumeFace& scvf,
71 const ElementFluxVarsCache& elemFluxVarCache)
74 if (!scvf.interiorBoundary())
75 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarCache);
77 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(),
"FacetCoupling.Xi", 1.0);
78 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
79 DUNE_THROW(Dune::NotImplemented,
"Xi != 1.0 cannot be used with the Box-Facet-Coupling scheme");
82 const auto& fluxVarCache = elemFluxVarCache[scvf];
83 const auto& shapeValues = fluxVarCache.shapeValues();
84 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
85 const auto& insideVolVars = elemVolVars[insideScv];
89 for (
const auto& scv : scvs(fvGeometry))
90 rho +=
massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx)*shapeValues[scv.indexInElement()][0];
93 ComponentFluxVector componentFlux(0.0);
94 const auto bcTypes = problem.interiorBoundaryTypes(element, scvf);
95 if (bcTypes.hasOnlyNeumann())
98 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
104 const auto a = facetVolVars.extrusionFactor();
105 auto preGradX = scvf.unitOuterNormal();
106 preGradX *= dim == dimWorld ? 0.5*a : 0.5*sqrt(a);
107 preGradX /= preGradX.two_norm2();
109 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
111 if constexpr (!FluidSystem::isTracerFluidSystem())
112 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
117 for (
const auto& scv : scvs(fvGeometry))
118 x +=
massOrMoleFraction(elemVolVars[scv], referenceSystem, phaseIdx, compIdx)*shapeValues[scv.indexInElement()][0];
121 auto gradX = preGradX;
124 componentFlux[compIdx] = -1.0*rho*Extrusion::area(fvGeometry, scvf)
125 *insideVolVars.extrusionFactor()
126 *
vtmv(scvf.unitOuterNormal(),
127 facetVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
130 if constexpr (!FluidSystem::isTracerFluidSystem())
131 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
132 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
135 return componentFlux;
139 else if (bcTypes.hasOnlyDirichlet())
141 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
143 if constexpr (!FluidSystem::isTracerFluidSystem())
144 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
148 std::vector<Scalar> xFractions(element.subEntities(dim));
149 for (
const auto& scv : scvs(fvGeometry))
150 xFractions[scv.localDofIndex()] =
massOrMoleFraction(elemVolVars[scv], referenceSystem, phaseIdx, compIdx);
153 for (
const auto& scvfJ : scvfs(fvGeometry))
154 if (scvfJ.interiorBoundary() && scvfJ.facetIndexInElement() == scvf.facetIndexInElement())
155 xFractions[ fvGeometry.scv(scvfJ.insideScvIdx()).localDofIndex() ]
156 =
massOrMoleFraction(problem.couplingManager().getLowDimVolVars(element, scvfJ), referenceSystem, phaseIdx, compIdx);
159 Dune::FieldVector<Scalar, dimWorld> gradX(0.0);
160 for (
const auto& scv : scvs(fvGeometry))
161 gradX.axpy(xFractions[scv.localDofIndex()], fluxVarCache.gradN(scv.indexInElement()));
164 componentFlux[compIdx] = -1.0*rho*Extrusion::area(fvGeometry, scvf)
165 *insideVolVars.extrusionFactor()
166 *
vtmv(scvf.unitOuterNormal(),
167 insideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
170 if constexpr (!FluidSystem::isTracerFluidSystem())
171 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
172 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
175 return componentFlux;
180 DUNE_THROW(Dune::NotImplemented,
"Mixed boundary types are not supported");
184 template<
class Problem,
class ElementVolumeVariables,
class FluxVarCache>
186 const Element& element,
187 const FVElementGeometry& fvGeometry,
188 const ElementVolumeVariables& elemVolVars,
189 const SubControlVolumeFace& scvf,
190 const FluxVarCache& fluxVarCache,
191 unsigned int phaseIdx)
193 DUNE_THROW(Dune::NotImplemented,
"transmissibilty computation for BoxFacetCouplingFicksLaw");
Ficks's law for the box scheme in the context of coupled models where coupling occurs across the face...
Definition: multidomain/facet/box/fickslaw.hh:41
static std::vector< Scalar > calculateTransmissibilities(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const FluxVarCache &fluxVarCache, unsigned int phaseIdx)
Definition: multidomain/facet/box/fickslaw.hh:185
static ComponentFluxVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVarsCache &elemFluxVarCache)
Definition: multidomain/facet/box/fickslaw.hh:65
Specialization of Fick's Law for the box method.
Definition: flux/box/fickslaw.hh:40
static ComponentFluxVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Returns the diffusive fluxes of all components within a fluid phase across the given sub-control volu...
Definition: flux/box/fickslaw.hh:82
forward declaration of the method-specific implementation
Definition: flux/box/fickslaw.hh:32
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
This file contains the data which is required to calculate diffusive fluxes due to molecular diffusio...
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:880
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:54
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:43
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.