24#ifndef DUMUX_DISCRETIZATION_BOX_FACET_COUPLING_FICKS_LAW_HH
25#define DUMUX_DISCRETIZATION_BOX_FACET_COUPLING_FICKS_LAW_HH
30#include <dune/common/exceptions.hh>
31#include <dune/common/fvector.hh>
32#include <dune/common/float_cmp.hh>
50template<
class TypeTag, ReferenceSystemFormulation referenceSystem = ReferenceSystemFormulation::massAveraged>
57 using FVElementGeometry =
typename GridGeometry::LocalView;
58 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
60 using GridView =
typename GridGeometry::GridView;
61 using Element =
typename GridView::template Codim<0>::Entity;
63 static constexpr int dim = GridView::dimension;
64 static constexpr int dimWorld = GridView::dimensionworld;
69 static constexpr int numComponents = ModelTraits::numFluidComponents();
72 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
76 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache>
77 static ComponentFluxVector
flux(
const Problem& problem,
78 const Element& element,
79 const FVElementGeometry& fvGeometry,
80 const ElementVolumeVariables& elemVolVars,
81 const SubControlVolumeFace& scvf,
83 const ElementFluxVarsCache& elemFluxVarCache)
86 if (!scvf.interiorBoundary())
87 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarCache);
89 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(),
"FacetCoupling.Xi", 1.0);
90 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
91 DUNE_THROW(Dune::NotImplemented,
"Xi != 1.0 cannot be used with the Box-Facet-Coupling scheme");
94 const auto& fluxVarCache = elemFluxVarCache[scvf];
95 const auto& shapeValues = fluxVarCache.shapeValues();
96 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
97 const auto& insideVolVars = elemVolVars[insideScv];
101 for (
const auto& scv : scvs(fvGeometry))
102 rho +=
massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx)*shapeValues[scv.indexInElement()][0];
105 ComponentFluxVector componentFlux(0.0);
106 const auto bcTypes = problem.interiorBoundaryTypes(element, scvf);
107 if (bcTypes.hasOnlyNeumann())
110 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
116 const auto a = facetVolVars.extrusionFactor();
117 auto preGradX = scvf.unitOuterNormal();
118 preGradX *= dim == dimWorld ? 0.5*a : 0.5*sqrt(a);
119 preGradX /= preGradX.two_norm2();
121 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
123 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
128 for (
const auto& scv : scvs(fvGeometry))
129 x +=
massOrMoleFraction(elemVolVars[scv], referenceSystem, phaseIdx, compIdx)*shapeValues[scv.indexInElement()][0];
132 auto gradX = preGradX;
135 componentFlux[compIdx] = -1.0*rho*Extrusion::area(scvf)
136 *insideVolVars.extrusionFactor()
137 *
vtmv(scvf.unitOuterNormal(),
138 facetVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
141 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
142 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
145 return componentFlux;
149 else if (bcTypes.hasOnlyDirichlet())
151 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
153 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
157 std::vector<Scalar> xFractions(element.subEntities(dim));
158 for (
const auto& scv : scvs(fvGeometry))
159 xFractions[scv.localDofIndex()] =
massOrMoleFraction(elemVolVars[scv], referenceSystem, phaseIdx, compIdx);
162 for (
const auto& scvfJ : scvfs(fvGeometry))
163 if (scvfJ.interiorBoundary() && scvfJ.facetIndexInElement() == scvf.facetIndexInElement())
164 xFractions[ fvGeometry.scv(scvfJ.insideScvIdx()).localDofIndex() ]
165 =
massOrMoleFraction(problem.couplingManager().getLowDimVolVars(element, scvfJ), referenceSystem, phaseIdx, compIdx);
168 Dune::FieldVector<Scalar, dimWorld> gradX(0.0);
169 for (
const auto& scv : scvs(fvGeometry))
170 gradX.axpy(xFractions[scv.localDofIndex()], fluxVarCache.gradN(scv.indexInElement()));
173 componentFlux[compIdx] = -1.0*rho*Extrusion::area(scvf)
174 *insideVolVars.extrusionFactor()
175 *
vtmv(scvf.unitOuterNormal(),
176 insideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
179 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
180 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
183 return componentFlux;
188 DUNE_THROW(Dune::NotImplemented,
"Mixed boundary types are not supported");
192 template<
class Problem,
class ElementVolumeVariables,
class FluxVarCache>
194 const Element& element,
195 const FVElementGeometry& fvGeometry,
196 const ElementVolumeVariables& elemVolVars,
197 const SubControlVolumeFace& scvf,
198 const FluxVarCache& fluxVarCache,
199 unsigned int phaseIdx)
201 DUNE_THROW(Dune::NotImplemented,
"transmissibilty computation for BoxFacetCouplingFicksLaw");
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
Helper classes to compute the integration elements.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
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
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:849
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
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
forward declaration of the method-specific implemetation
Definition: flux/box/fickslaw.hh:43
Specialization of Fick's Law for the box method.
Definition: flux/box/fickslaw.hh:51
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:92
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:53
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:193
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:77
Declares all properties used in Dumux.
This file contains the data which is required to calculate diffusive fluxes due to molecular diffusio...