version 3.8
flux/box/fickslaw.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_DISCRETIZATION_BOX_FICKS_LAW_HH
14#define DUMUX_DISCRETIZATION_BOX_FICKS_LAW_HH
15
16#include <dune/common/fvector.hh>
17#include <dune/common/fmatrix.hh>
18
19#include <dumux/common/math.hh>
23
27
28namespace Dumux {
29
30// forward declaration
31template<class TypeTag, class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
33
38template <class TypeTag, ReferenceSystemFormulation referenceSystem>
39class FicksLawImplementation<TypeTag, DiscretizationMethods::Box, referenceSystem>
40{
46 using FVElementGeometry = typename GridGeometry::LocalView;
47 using SubControlVolume = typename GridGeometry::SubControlVolume;
48 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
49 using Extrusion = Extrusion_t<GridGeometry>;
50 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
52 using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView;
53 using FluxVarCache = typename GridFluxVariablesCache::FluxVariablesCache;
56 using Element = typename GridView::template Codim<0>::Entity;
58
59 enum { dim = GridView::dimension} ;
60 enum { dimWorld = GridView::dimensionworld} ;
61 enum
62 {
63 numPhases = ModelTraits::numFluidPhases(),
64 numComponents = ModelTraits::numFluidComponents()
65 };
66 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
67 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
68
69public:
71
72 //return the reference system
74 { return referenceSystem; }
75
82 static ComponentFluxVector flux(const Problem& problem,
83 const Element& element,
84 const FVElementGeometry& fvGeometry,
85 const ElementVolumeVariables& elemVolVars,
86 const SubControlVolumeFace& scvf,
87 const int phaseIdx,
88 const ElementFluxVariablesCache& elemFluxVarsCache)
89 {
90 // get inside and outside diffusion tensors and calculate the harmonic mean
91 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
92 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
93
94 // evaluate gradX at integration point and interpolate density
95 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
96 const auto& shapeValues = fluxVarsCache.shapeValues();
97
98 // density interpolation
99 Scalar rho(0.0);
100 for (auto&& scv : scvs(fvGeometry))
101 rho += massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx)*shapeValues[scv.indexInElement()][0];
102
103 ComponentFluxVector componentFlux(0.0);
104 for (int compIdx = 0; compIdx < numComponents; compIdx++)
105 {
106 if constexpr (!FluidSystem::isTracerFluidSystem())
107 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
108 continue;
109
110 const auto diffCoeff = averageDiffusionCoefficient_(phaseIdx, compIdx, insideVolVars, outsideVolVars, problem, scvf);
111
112 // compute the diffusive flux
113 const auto massOrMoleFrac = [&](const SubControlVolume& scv){ return massOrMoleFraction(elemVolVars[scv], referenceSystem, phaseIdx, compIdx); };
114 componentFlux[compIdx] = discreteFlux_(fvGeometry, scvf, fluxVarsCache, massOrMoleFrac, diffCoeff, rho);
115
116 // if the main component is balanced subtract the same flux from there (conservation)
117 if constexpr (!FluidSystem::isTracerFluidSystem())
118 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
119 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
120 }
121 return componentFlux;
122 }
123
124 // compute transmissibilities ti for analytical jacobians
125 static std::array<std::vector<Scalar>, numComponents>
126 calculateTransmissibilities(const Problem& problem,
127 const Element& element,
128 const FVElementGeometry& fvGeometry,
129 const ElementVolumeVariables& elemVolVars,
130 const SubControlVolumeFace& scvf,
131 const FluxVarCache& fluxVarCache,
132 const int phaseIdx)
133 {
134 Scalar rho(0.0);
135 const auto& shapeValues = fluxVarCache.shapeValues();
136 for (auto&& scv : scvs(fvGeometry))
137 rho += massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx)*shapeValues[scv.indexInElement()][0];
138
139 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
140 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
141
142 std::array<std::vector<Scalar>, numComponents> ti;
143 for (int compIdx = 0; compIdx < numComponents; compIdx++)
144 {
145 if constexpr (!FluidSystem::isTracerFluidSystem())
146 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
147 continue;
148
149 const auto diffCoeff = averageDiffusionCoefficient_(phaseIdx, compIdx, insideVolVars, outsideVolVars, problem, scvf);
150
151 ti[compIdx].resize(fvGeometry.numScv());
152 for (auto&& scv : scvs(fvGeometry))
153 ti[compIdx][scv.indexInElement()] = -rho*vtmv(scvf.unitOuterNormal(), diffCoeff, fluxVarCache.gradN(scv.indexInElement()))*Extrusion::area(fvGeometry, scvf);
154 }
155
156 return ti;
157 }
158
159private:
160 static Scalar averageDiffusionCoefficient_(const int phaseIdx, const int compIdx,
161 const VolumeVariables& insideVV, const VolumeVariables& outsideVV,
162 const Problem& problem,
163 const SubControlVolumeFace& scvf)
164 {
165 // effective diffusion tensors
166 auto [insideDiffCoeff, outsideDiffCoeff] = diffusionCoefficientsAtInterface_(phaseIdx, compIdx, insideVV, outsideVV);
167
168 // scale by extrusion factor
169 insideDiffCoeff *= insideVV.extrusionFactor();
170 outsideDiffCoeff *= outsideVV.extrusionFactor();
171
172 // the resulting averaged diffusion tensor
173 return faceTensorAverage(insideDiffCoeff, outsideDiffCoeff, scvf.unitOuterNormal());
174 }
175
176 static std::pair<Scalar, Scalar>
177 diffusionCoefficientsAtInterface_([[maybe_unused]] const int phaseIdx, const int compIdx,
178 const VolumeVariables& insideVV, const VolumeVariables& outsideVV)
179 {
180 if constexpr (!FluidSystem::isTracerFluidSystem())
181 {
182 const auto mainCompIdx = FluidSystem::getMainComponent(phaseIdx);
183 const auto insideDiffCoeff = insideVV.effectiveDiffusionCoefficient(phaseIdx, mainCompIdx, compIdx);
184 const auto outsideDiffCoeff = outsideVV.effectiveDiffusionCoefficient(phaseIdx, mainCompIdx, compIdx);
185 return { std::move(insideDiffCoeff), std::move(outsideDiffCoeff) };
186 }
187 else
188 {
189 const auto insideDiffCoeff = insideVV.effectiveDiffusionCoefficient(0, 0, compIdx);
190 const auto outsideDiffCoeff = outsideVV.effectiveDiffusionCoefficient(0, 0, compIdx);
191 return { std::move(insideDiffCoeff), std::move(outsideDiffCoeff) };
192 }
193 }
194
195 template<class EvaluateVariable, class Tensor>
196 static Scalar discreteFlux_(const FVElementGeometry& fvGeometry,
197 const SubControlVolumeFace& scvf,
198 const FluxVarCache& fluxVarsCache,
199 const EvaluateVariable& massOrMoleFraction,
200 const Tensor& D, const Scalar preFactor)
201 {
202 Dune::FieldVector<Scalar, dimWorld> gradX(0.0);
203 for (auto&& scv : scvs(fvGeometry))
204 gradX.axpy(massOrMoleFraction(scv), fluxVarsCache.gradN(scv.indexInElement()));
205 return -1.0*preFactor*vtmv(scvf.unitOuterNormal(), D, gradX)*Extrusion::area(fvGeometry, scvf);
206 }
207};
208
209} // end namespace Dumux
210
211#endif
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
Definition: fickiandiffusioncoefficients.hh:32
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: flux/box/fickslaw.hh:126
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: flux/box/fickslaw.hh:73
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.
A free function to average a Tensor at an interface.
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
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:851
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
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition: referencesystemformulation.hh:33
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Define some often used mathematical functions.
The available discretization methods in Dumux.
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:94
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
Scalar faceTensorAverage(const Scalar T1, const Scalar T2, const Dune::FieldVector< Scalar, dim > &normal)
Average of a discontinuous scalar field at discontinuity interface (for compatibility reasons with th...
Definition: facetensoraverage.hh:29
The reference frameworks and formulations available for splitting total fluxes into a advective and d...