3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
flux/cctpfa/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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_DISCRETIZATION_CC_TPFA_FICKS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_TPFA_FICKS_LAW_HH
26
27#include <dune/common/fvector.hh>
28
32
36
38
39namespace Dumux {
40
41// forward declaration
42template<class TypeTag, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
43class FicksLawImplementation;
44
49template <class TypeTag, ReferenceSystemFormulation referenceSystem>
50class FicksLawImplementation<TypeTag, DiscretizationMethod::cctpfa, referenceSystem>
51{
55 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
56 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
58 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
59 using Element = typename GridView::template Codim<0>::Entity;
60 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
64
66 static const int dim = GridView::dimension;
67 static const int dimWorld = GridView::dimensionworld;
68 static const int numPhases = ModelTraits::numFluidPhases();
69 static const int numComponents = ModelTraits::numFluidComponents();
70
71 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
72
74 class TpfaFicksLawCacheFiller
75 {
76 public:
79 template<class FluxVariablesCacheFiller>
80 static void fill(FluxVariablesCache& scvfFluxVarsCache,
81 unsigned int phaseIdx, unsigned int compIdx,
82 const Problem& problem,
83 const Element& element,
84 const FVElementGeometry& fvGeometry,
85 const ElementVolumeVariables& elemVolVars,
86 const SubControlVolumeFace& scvf,
87 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
88 {
89 scvfFluxVarsCache.updateDiffusion(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
90 }
91 };
92
94 class TpfaFicksLawCache
95 {
96 public:
97 using Filler = TpfaFicksLawCacheFiller;
98
99 void updateDiffusion(const Problem& problem,
100 const Element& element,
101 const FVElementGeometry& fvGeometry,
102 const ElementVolumeVariables& elemVolVars,
103 const SubControlVolumeFace &scvf,
104 const unsigned int phaseIdx,
105 const unsigned int compIdx)
106 {
107 tij_[phaseIdx][compIdx] = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
108 }
109
110 const Scalar& diffusionTij(unsigned int phaseIdx, unsigned int compIdx) const
111 { return tij_[phaseIdx][compIdx]; }
112
113 private:
114 std::array< std::array<Scalar, numComponents>, numPhases> tij_;
115 };
116
117public:
122 { return referenceSystem; }
123
125 using Cache = TpfaFicksLawCache;
126
128
130 static ComponentFluxVector flux(const Problem& problem,
131 const Element& element,
132 const FVElementGeometry& fvGeometry,
133 const ElementVolumeVariables& elemVolVars,
134 const SubControlVolumeFace& scvf,
135 const int phaseIdx,
136 const ElementFluxVariablesCache& elemFluxVarsCache)
137 {
138 ComponentFluxVector componentFlux(0.0);
139 for (int compIdx = 0; compIdx < numComponents; compIdx++)
140 {
141 if constexpr (!FluidSystem::isTracerFluidSystem())
142 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
143 continue;
144
145 // diffusion tensors are always solution dependent
146 Scalar tij = elemFluxVarsCache[scvf].diffusionTij(phaseIdx, compIdx);
147
148 // get inside/outside volume variables
149 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
150 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
151
152 // the inside and outside mass/mole fractions fractions
153 const Scalar xInside = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
154 const Scalar massOrMoleFractionOutside = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
155 const Scalar xOutside = scvf.numOutsideScvs() == 1 ? massOrMoleFractionOutside
156 : branchingFacetX(problem, element, fvGeometry, elemVolVars,
157 elemFluxVarsCache, scvf, xInside, tij, phaseIdx, compIdx);
158
159 const Scalar rhoInside = massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
160 const Scalar rhoOutside = massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
161
162 const Scalar rho = scvf.numOutsideScvs() == 1 ? 0.5*(rhoInside + rhoOutside)
163 : branchingFacetDensity(elemVolVars, scvf, phaseIdx, rhoInside);
164
165 componentFlux[compIdx] = rho*tij*(xInside - xOutside);
166 if constexpr (!FluidSystem::isTracerFluidSystem())
167 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
168 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
169 }
170
171 return componentFlux;
172 }
173
175 static Scalar calculateTransmissibility(const Problem& problem,
176 const Element& element,
177 const FVElementGeometry& fvGeometry,
178 const ElementVolumeVariables& elemVolVars,
179 const SubControlVolumeFace& scvf,
180 const int phaseIdx, const int compIdx)
181 {
182
183
184 const auto insideScvIdx = scvf.insideScvIdx();
185 const auto& insideScv = fvGeometry.scv(insideScvIdx);
186 const auto& insideVolVars = elemVolVars[insideScvIdx];
187 const auto getDiffCoeff = [&](const auto& vv)
188 {
190 if constexpr (FluidSystem::isTracerFluidSystem())
191 return Deprecated::template effectiveDiffusionCoefficient<EffDiffModel>(vv, 0, 0, compIdx);
192 else
193 return Deprecated::template effectiveDiffusionCoefficient<EffDiffModel>(vv, phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
194 };
195
196 const auto insideD = getDiffCoeff(insideVolVars);
197
198 const Scalar ti = computeTpfaTransmissibility(scvf, insideScv, insideD, insideVolVars.extrusionFactor());
199
200 // for the boundary (dirichlet) or at branching points we only need ti
201 Scalar tij;
202 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
203 tij = scvf.area()*ti;
204
205 // otherwise we compute a tpfa harmonic mean
206 else
207 {
208 const auto outsideScvIdx = scvf.outsideScvIdx();
209 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
210 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
211 const auto outsideD = getDiffCoeff(outsideVolVars);
212
213 Scalar tj;
214 if (dim == dimWorld)
215 // assume the normal vector from outside is anti parallel so we save flipping a vector
216 tj = -1.0*computeTpfaTransmissibility(scvf, outsideScv, outsideD, outsideVolVars.extrusionFactor());
217 else
218 tj = computeTpfaTransmissibility(fvGeometry.flipScvf(scvf.index()),
219 outsideScv,
220 outsideD,
221 outsideVolVars.extrusionFactor());
222
223 // check if we are dividing by zero!
224 if (ti*tj <= 0.0)
225 tij = 0;
226 else
227 tij = scvf.area()*(ti * tj)/(ti + tj);
228 }
229
230 return tij;
231 }
232
233private:
235 static Scalar branchingFacetX(const Problem& problem,
236 const Element& element,
237 const FVElementGeometry& fvGeometry,
238 const ElementVolumeVariables& elemVolVars,
239 const ElementFluxVariablesCache& elemFluxVarsCache,
240 const SubControlVolumeFace& scvf,
241 const Scalar insideX, const Scalar insideTi,
242 const int phaseIdx, const int compIdx)
243 {
244 Scalar sumTi(insideTi);
245 Scalar sumXTi(insideTi*insideX);
246
247 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
248 {
249 const auto outsideScvIdx = scvf.outsideScvIdx(i);
250 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
251 const Scalar massOrMoleFractionOutside = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
252 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
253
254 const Scalar outsideTi = elemFluxVarsCache[flippedScvf].diffusionTij(phaseIdx, compIdx);
255 sumTi += outsideTi;
256 sumXTi += outsideTi*massOrMoleFractionOutside;
257 }
258
259 return sumTi > 0 ? sumXTi/sumTi : 0;
260 }
261
263 static Scalar branchingFacetDensity(const ElementVolumeVariables& elemVolVars,
264 const SubControlVolumeFace& scvf,
265 const int phaseIdx,
266 const Scalar insideRho)
267 {
268 Scalar rho(insideRho);
269 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
270 {
271 const auto outsideScvIdx = scvf.outsideScvIdx(i);
272 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
273 const Scalar rhoOutside = massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
274 rho += rhoOutside;
275 }
276 return rho/(scvf.numOutsideScvs()+1);
277 }
278};
279
280} // end namespace Dumux
281
282#endif
Helpers for deprecation.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
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
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition: referencesystemformulation.hh:45
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
forward declaration of the method-specific implemetation
Definition: flux/box/fickslaw.hh:43
Fick's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: flux/cctpfa/fickslaw.hh:51
static Scalar calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const int compIdx)
compute diffusive transmissibilities
Definition: flux/cctpfa/fickslaw.hh:175
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Return the reference system.
Definition: flux/cctpfa/fickslaw.hh:121
TpfaFicksLawCache Cache
state the type for the corresponding cache and its filler
Definition: flux/cctpfa/fickslaw.hh:125
static ComponentFluxVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
return diffusive fluxes for all components in a phase
Definition: flux/cctpfa/fickslaw.hh:130
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
Definition: fickiandiffusioncoefficients.hh:44
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...