3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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
29
32
34
35namespace Dumux {
36
37// forward declaration
38template<class TypeTag, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
39class FicksLawImplementation;
40
45template <class TypeTag, ReferenceSystemFormulation referenceSystem>
46class FicksLawImplementation<TypeTag, DiscretizationMethod::cctpfa, referenceSystem>
47{
51 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
52 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
53 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
56 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
57 using Element = typename GridView::template Codim<0>::Entity;
58 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
62
64 using Indices = typename ModelTraits::Indices;
65
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 DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
72 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
73
75 class TpfaFicksLawCacheFiller
76 {
77 public:
80 template<class FluxVariablesCacheFiller>
81 static void fill(FluxVariablesCache& scvfFluxVarsCache,
82 unsigned int phaseIdx, unsigned int compIdx,
83 const Problem& problem,
84 const Element& element,
85 const FVElementGeometry& fvGeometry,
86 const ElementVolumeVariables& elemVolVars,
87 const SubControlVolumeFace& scvf,
88 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
89 {
90 scvfFluxVarsCache.updateDiffusion(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
91 }
92 };
93
95 class TpfaFicksLawCache
96 {
97 public:
98 using Filler = TpfaFicksLawCacheFiller;
99
100 void updateDiffusion(const Problem& problem,
101 const Element& element,
102 const FVElementGeometry& fvGeometry,
103 const ElementVolumeVariables& elemVolVars,
104 const SubControlVolumeFace &scvf,
105 const unsigned int phaseIdx,
106 const unsigned int compIdx)
107 {
108 tij_[phaseIdx][compIdx] = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
109 }
110
111 const Scalar& diffusionTij(unsigned int phaseIdx, unsigned int compIdx) const
112 { return tij_[phaseIdx][compIdx]; }
113
114 private:
115 std::array< std::array<Scalar, numComponents>, numPhases> tij_;
116 };
117
118public:
123 { return referenceSystem; }
124
126 using Cache = TpfaFicksLawCache;
127
129 static ComponentFluxVector flux(const Problem& problem,
130 const Element& element,
131 const FVElementGeometry& fvGeometry,
132 const ElementVolumeVariables& elemVolVars,
133 const SubControlVolumeFace& scvf,
134 const int phaseIdx,
135 const ElementFluxVariablesCache& elemFluxVarsCache)
136 {
137 ComponentFluxVector componentFlux(0.0);
138 for (int compIdx = 0; compIdx < numComponents; compIdx++)
139 {
140 if(compIdx == FluidSystem::getMainComponent(phaseIdx))
141 continue;
142
143 // diffusion tensors are always solution dependent
144 Scalar tij = elemFluxVarsCache[scvf].diffusionTij(phaseIdx, compIdx);
145
146 // get inside/outside volume variables
147 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
148 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
149
150 // the inside and outside mass/mole fractions fractions
151 const Scalar xInside = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
152 const Scalar massOrMoleFractionOutside = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
153 const Scalar xOutside = scvf.numOutsideScvs() == 1 ? massOrMoleFractionOutside
154 : branchingFacetX(problem, element, fvGeometry, elemVolVars,
155 elemFluxVarsCache, scvf, xInside, tij, phaseIdx, compIdx);
156
157 const Scalar rhoInside = massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
158 const Scalar rhoOutside = massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
159
160 const Scalar rho = scvf.numOutsideScvs() == 1 ? 0.5*(rhoInside + rhoOutside)
161 : branchingFacetDensity(elemVolVars, scvf, phaseIdx, rhoInside);
162
163 componentFlux[compIdx] = rho*tij*(xInside - xOutside);
164 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
165 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
166 }
167
168 return componentFlux;
169 }
170
172 static Scalar calculateTransmissibility(const Problem& problem,
173 const Element& element,
174 const FVElementGeometry& fvGeometry,
175 const ElementVolumeVariables& elemVolVars,
176 const SubControlVolumeFace& scvf,
177 const int phaseIdx, const int compIdx)
178 {
179 Scalar tij;
180
181 const auto insideScvIdx = scvf.insideScvIdx();
182 const auto& insideScv = fvGeometry.scv(insideScvIdx);
183 const auto& insideVolVars = elemVolVars[insideScvIdx];
184
186 const auto insideD = EffDiffModel::effectiveDiffusivity(insideVolVars.porosity(),
187 insideVolVars.saturation(phaseIdx),
188 insideVolVars.diffusionCoefficient(phaseIdx, compIdx));
189 const Scalar ti = computeTpfaTransmissibility(scvf, insideScv, insideD, insideVolVars.extrusionFactor());
190
191 // for the boundary (dirichlet) or at branching points we only need ti
192 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
193 {
194 tij = scvf.area()*ti;
195 }
196 // otherwise we compute a tpfa harmonic mean
197 else
198 {
199 const auto outsideScvIdx = scvf.outsideScvIdx();
200 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
201 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
202
203 const auto outsideD = EffDiffModel::effectiveDiffusivity(outsideVolVars.porosity(),
204 outsideVolVars.saturation(phaseIdx),
205 outsideVolVars.diffusionCoefficient(phaseIdx, compIdx));
206
207 Scalar tj;
208 if (dim == dimWorld)
209 // assume the normal vector from outside is anti parallel so we save flipping a vector
210 tj = -1.0*computeTpfaTransmissibility(scvf, outsideScv, outsideD, outsideVolVars.extrusionFactor());
211 else
212 tj = computeTpfaTransmissibility(fvGeometry.flipScvf(scvf.index()),
213 outsideScv,
214 outsideD,
215 outsideVolVars.extrusionFactor());
216
217 // check if we are dividing by zero!
218 if (ti*tj <= 0.0)
219 tij = 0;
220 else
221 tij = scvf.area()*(ti * tj)/(ti + tj);
222 }
223
224 return tij;
225 }
226
227private:
229 static Scalar branchingFacetX(const Problem& problem,
230 const Element& element,
231 const FVElementGeometry& fvGeometry,
232 const ElementVolumeVariables& elemVolVars,
233 const ElementFluxVariablesCache& elemFluxVarsCache,
234 const SubControlVolumeFace& scvf,
235 const Scalar insideX, const Scalar insideTi,
236 const int phaseIdx, const int compIdx)
237 {
238 Scalar sumTi(insideTi);
239 Scalar sumXTi(insideTi*insideX);
240
241 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
242 {
243 const auto outsideScvIdx = scvf.outsideScvIdx(i);
244 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
245 const Scalar massOrMoleFractionOutside = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
246 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
247
248 const Scalar outsideTi = elemFluxVarsCache[flippedScvf].diffusionTij(phaseIdx, compIdx);
249 sumTi += outsideTi;
250 sumXTi += outsideTi*massOrMoleFractionOutside;
251 }
252
253 return sumTi > 0 ? sumXTi/sumTi : 0;
254 }
255
257 static Scalar branchingFacetDensity(const ElementVolumeVariables& elemVolVars,
258 const SubControlVolumeFace& scvf,
259 const int phaseIdx,
260 const Scalar insideRho)
261 {
262 Scalar rho(insideRho);
263 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
264 {
265 const auto outsideScvIdx = scvf.outsideScvIdx(i);
266 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
267 const Scalar rhoOutside = massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
268 rho += rhoOutside;
269 }
270 return rho/(scvf.numOutsideScvs()+1);
271 }
272};
273
274} // end namespace Dumux
275
276#endif
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...
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
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
Fick's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: cctpfa/fickslaw.hh:47
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: cctpfa/fickslaw.hh:172
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Return the reference system.
Definition: cctpfa/fickslaw.hh:122
TpfaFicksLawCache Cache
state the type for the corresponding cache and its filler
Definition: cctpfa/fickslaw.hh:126
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: cctpfa/fickslaw.hh:129
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...