25#ifndef DUMUX_DISCRETIZATION_BOX_FICKS_LAW_HH
26#define DUMUX_DISCRETIZATION_BOX_FICKS_LAW_HH
37template<
class TypeTag, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
44template <
class TypeTag, ReferenceSystemFormulation referenceSystem>
52 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
53 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
60 using Element =
typename GridView::template Codim<0>::Entity;
64 enum { dim = GridView::dimension} ;
65 enum { dimWorld = GridView::dimensionworld} ;
68 numPhases = ModelTraits::numFluidPhases(),
69 numComponents = ModelTraits::numFluidComponents()
71 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
72 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
77 {
return referenceSystem; }
79 static ComponentFluxVector
flux(
const Problem& problem,
80 const Element& element,
81 const FVElementGeometry& fvGeometry,
82 const ElementVolumeVariables& elemVolVars,
83 const SubControlVolumeFace& scvf,
85 const ElementFluxVariablesCache& elemFluxVarsCache)
87 ComponentFluxVector componentFlux(0.0);
90 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
91 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
94 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
95 const auto& shapeValues = fluxVarsCache.shapeValues();
99 for (
auto&& scv : scvs(fvGeometry))
100 rho +=
massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx)*shapeValues[scv.indexInElement()][0];
102 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
104 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
109 auto insideD = EffDiffModel::effectiveDiffusivity(insideVolVars.porosity(),
110 insideVolVars.saturation(phaseIdx),
111 insideVolVars.diffusionCoefficient(phaseIdx, compIdx));
112 auto outsideD = EffDiffModel::effectiveDiffusivity(outsideVolVars.porosity(),
113 outsideVolVars.saturation(phaseIdx),
114 outsideVolVars.diffusionCoefficient(phaseIdx, compIdx));
117 insideD *= insideVolVars.extrusionFactor();
118 outsideD *= outsideVolVars.extrusionFactor();
121 const auto D = problem.spatialParams().harmonicMean(insideD, outsideD, scvf.unitOuterNormal());
124 Dune::FieldVector<Scalar, dimWorld> gradX(0.0);
125 for (
auto&& scv : scvs(fvGeometry))
126 gradX.axpy(
massOrMoleFraction(elemVolVars[scv], referenceSystem, phaseIdx, compIdx), fluxVarsCache.gradN(scv.indexInElement()));
129 componentFlux[compIdx] = -1.0*rho*
vtmv(scvf.unitOuterNormal(), D, gradX)*scvf.area();
130 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
131 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
133 return componentFlux;
137 static std::array<std::vector<Scalar>, numComponents>
139 const Element& element,
140 const FVElementGeometry& fvGeometry,
141 const ElementVolumeVariables& elemVolVars,
142 const SubControlVolumeFace& scvf,
143 const FluxVarCache& fluxVarCache,
147 const auto& shapeValues = fluxVarCache.shapeValues();
148 for (
auto&& scv : scvs(fvGeometry))
149 rho +=
massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx)*shapeValues[scv.indexInElement()][0];
151 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
152 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
154 std::array<std::vector<Scalar>, numComponents> ti;
155 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
157 if(compIdx == FluidSystem::getMainComponent(phaseIdx))
162 auto insideD = EffDiffModel::effectiveDiffusivity(insideVolVars.porosity(),
163 insideVolVars.saturation(phaseIdx),
164 insideVolVars.diffusionCoefficient(phaseIdx, compIdx));
165 auto outsideD = EffDiffModel::effectiveDiffusivity(outsideVolVars.porosity(),
166 outsideVolVars.saturation(phaseIdx),
167 outsideVolVars.diffusionCoefficient(phaseIdx, compIdx));
170 insideD *= insideVolVars.extrusionFactor();
171 outsideD *= outsideVolVars.extrusionFactor();
174 const auto D = problem.spatialParams().harmonicMean(insideD, outsideD, scvf.unitOuterNormal());
176 ti[compIdx].resize(fvGeometry.numScv());
177 for (
auto&& scv : scvs(fvGeometry))
178 ti[compIdx][scv.indexInElement()] = -rho*
vtmv(scvf.unitOuterNormal(), D, fluxVarCache.gradN(scv.indexInElement()))*scvf.area();
Define some often used mathematical functions.
The available discretization methods in Dumux.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
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
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition: referencesystemformulation.hh:45
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
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
Property tag Indices
Definition: porousmediumflow/sequential/properties.hh:59
forward declaration of the method-specific implemetation
Definition: box/fickslaw.hh:38
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: box/fickslaw.hh:76
static ComponentFluxVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Definition: box/fickslaw.hh:79
static std::array< std::vector< Scalar >, numComponents > calculateTransmissibilities(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const FluxVarCache &fluxVarCache, const int phaseIdx)
Definition: box/fickslaw.hh:138
Declares all properties used in Dumux.