diffusionmodule.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 This file is part of the Open Porous Media project (OPM).
5
6 OPM 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 2 of the License, or
9 (at your option) any later version.
10
11 OPM 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 OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef OPM_DIFFUSION_MODULE_HH
29#define OPM_DIFFUSION_MODULE_HH
30
33
34#include <opm/material/common/Valgrind.hpp>
35
36#include <dune/common/fvector.hh>
37
38namespace Opm {
39
46template <class TypeTag, bool enableDiffusion>
48
52template <class TypeTag>
53class DiffusionModule<TypeTag, /*enableDiffusion=*/false>
54{
58
59public:
63 static void registerParameters()
64 {}
65
70 template <class Context>
71 static void addDiffusiveFlux(RateVector&,
72 const Context&,
73 unsigned,
74 unsigned)
75 {}
76};
77
81template <class TypeTag>
82class DiffusionModule<TypeTag, /*enableDiffusion=*/true>
83{
89
90 enum { numPhases = FluidSystem::numPhases };
91 enum { numComponents = FluidSystem::numComponents };
92 enum { conti0EqIdx = Indices::conti0EqIdx };
93
94 using Toolbox = Opm::MathToolbox<Evaluation>;
95
96public:
100 static void registerParameters()
101 {}
102
107 template <class Context>
108 static void addDiffusiveFlux(RateVector& flux, const Context& context,
109 unsigned spaceIdx, unsigned timeIdx)
110 {
111 const auto& extQuants = context.extensiveQuantities(spaceIdx, timeIdx);
112
113 const auto& fluidStateI = context.intensiveQuantities(extQuants.interiorIndex(), timeIdx).fluidState();
114 const auto& fluidStateJ = context.intensiveQuantities(extQuants.exteriorIndex(), timeIdx).fluidState();
115
116 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
117 // arithmetic mean of the phase's molar density
118 Evaluation rhoMolar = fluidStateI.molarDensity(phaseIdx);
119 rhoMolar += Toolbox::value(fluidStateJ.molarDensity(phaseIdx));
120 rhoMolar /= 2;
121
122 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
123 // mass flux due to molecular diffusion
124 flux[conti0EqIdx + compIdx] +=
125 -rhoMolar
126 * extQuants.moleFractionGradientNormal(phaseIdx, compIdx)
127 * extQuants.effectiveDiffusionCoefficient(phaseIdx, compIdx);
128 }
129 }
130};
131
139template <class TypeTag, bool enableDiffusion>
141
145template <class TypeTag>
146class DiffusionIntensiveQuantities<TypeTag, /*enableDiffusion=*/false>
147{
151
152public:
157 Scalar tortuosity(unsigned) const
158 {
159 throw std::logic_error("Method tortuosity() does not make sense "
160 "if diffusion is disabled");
161 }
162
167 Scalar diffusionCoefficient(unsigned, unsigned) const
168 {
169 throw std::logic_error("Method diffusionCoefficient() does not "
170 "make sense if diffusion is disabled");
171 }
172
177 Scalar effectiveDiffusionCoefficient(unsigned, unsigned) const
178 {
179 throw std::logic_error("Method effectiveDiffusionCoefficient() "
180 "does not make sense if diffusion is disabled");
181 }
182
183protected:
188 template <class FluidState>
189 void update_(FluidState&,
190 typename FluidSystem::template ParameterCache<typename FluidState::Scalar>&,
191 const ElementContext&,
192 unsigned,
193 unsigned)
194 { }
195};
196
200template <class TypeTag>
201class DiffusionIntensiveQuantities<TypeTag, /*enableDiffusion=*/true>
202{
207
208 enum { numPhases = FluidSystem::numPhases };
209 enum { numComponents = FluidSystem::numComponents };
210
211public:
216 Evaluation diffusionCoefficient(unsigned phaseIdx, unsigned compIdx) const
217 { return diffusionCoefficient_[phaseIdx][compIdx]; }
218
223 Evaluation tortuosity(unsigned phaseIdx) const
224 { return tortuosity_[phaseIdx]; }
225
230 Evaluation effectiveDiffusionCoefficient(unsigned phaseIdx, unsigned compIdx) const
231 { return tortuosity_[phaseIdx] * diffusionCoefficient_[phaseIdx][compIdx]; }
232
233protected:
238 template <class FluidState>
239 void update_(FluidState& fluidState,
240 typename FluidSystem::template ParameterCache<typename FluidState::Scalar>& paramCache,
241 const ElementContext& elemCtx,
242 unsigned dofIdx,
243 unsigned timeIdx)
244 {
245 using Toolbox = Opm::MathToolbox<Evaluation>;
246
247 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
248 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
249 if (!elemCtx.model().phaseIsConsidered(phaseIdx))
250 continue;
251
252 // TODO: let the problem do this (this is a constitutive
253 // relation of which the model should be free of from the
254 // abstraction POV!)
255 // Based on Millington, R. J., & Quirk, J. P. (1961).
256 const Evaluation& base =
257 Toolbox::max(0.0001,
258 intQuants.porosity()
259 * intQuants.fluidState().saturation(phaseIdx));
260 tortuosity_[phaseIdx] =
261 1.0 / (intQuants.porosity() * intQuants.porosity())
262 * Toolbox::pow(base, 10.0/3.0);
263
264 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
265 diffusionCoefficient_[phaseIdx][compIdx] =
266 FluidSystem::diffusionCoefficient(fluidState,
267 paramCache,
268 phaseIdx,
269 compIdx);
270 }
271 }
272 }
273
274private:
275 Evaluation tortuosity_[numPhases];
276 Evaluation diffusionCoefficient_[numPhases][numComponents];
277};
278
285template <class TypeTag, bool enableDiffusion>
287
291template <class TypeTag>
292class DiffusionExtensiveQuantities<TypeTag, /*enableDiffusion=*/false>
293{
297
298protected:
303 void update_(const ElementContext&,
304 unsigned,
305 unsigned)
306 {}
307
308 template <class Context, class FluidState>
309 void updateBoundary_(const Context&,
310 unsigned,
311 unsigned,
312 const FluidState&)
313 {}
314
315public:
322 const Evaluation& moleFractionGradientNormal(unsigned,
323 unsigned) const
324 {
325 throw std::logic_error("The method moleFractionGradient() does not "
326 "make sense if diffusion is disabled.");
327 }
328
336 const Evaluation& effectiveDiffusionCoefficient(unsigned,
337 unsigned) const
338 {
339 throw std::logic_error("The method effectiveDiffusionCoefficient() "
340 "does not make sense if diffusion is disabled.");
341 }
342};
343
347template <class TypeTag>
348class DiffusionExtensiveQuantities<TypeTag, /*enableDiffusion=*/true>
349{
354
355 enum { dimWorld = GridView::dimensionworld };
356 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
357 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
358
359 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
360 using DimEvalVector = Dune::FieldVector<Evaluation, dimWorld>;
361
362protected:
367 void update_(const ElementContext& elemCtx, unsigned faceIdx, unsigned timeIdx)
368 {
369 const auto& gradCalc = elemCtx.gradientCalculator();
370 Opm::MoleFractionCallback<TypeTag> moleFractionCallback(elemCtx);
371
372 const auto& face = elemCtx.stencil(timeIdx).interiorFace(faceIdx);
373 const auto& normal = face.normal();
374 const auto& extQuants = elemCtx.extensiveQuantities(faceIdx, timeIdx);
375
376 const auto& intQuantsInside = elemCtx.intensiveQuantities(extQuants.interiorIndex(), timeIdx);
377 const auto& intQuantsOutside = elemCtx.intensiveQuantities(extQuants.exteriorIndex(), timeIdx);
378
379 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
380 if (!elemCtx.model().phaseIsConsidered(phaseIdx))
381 continue;
382
383 moleFractionCallback.setPhaseIndex(phaseIdx);
384 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
385 moleFractionCallback.setComponentIndex(compIdx);
386
387 DimEvalVector moleFractionGradient(0.0);
388 gradCalc.calculateGradient(moleFractionGradient,
389 elemCtx,
390 faceIdx,
391 moleFractionCallback);
392
393 moleFractionGradientNormal_[phaseIdx][compIdx] = 0.0;
394 for (unsigned i = 0; i < normal.size(); ++i)
395 moleFractionGradientNormal_[phaseIdx][compIdx] +=
396 normal[i]*moleFractionGradient[i];
397 Opm::Valgrind::CheckDefined(moleFractionGradientNormal_[phaseIdx][compIdx]);
398
399 // use the arithmetic average for the effective
400 // diffusion coefficients.
401 effectiveDiffusionCoefficient_[phaseIdx][compIdx] =
402 (intQuantsInside.effectiveDiffusionCoefficient(phaseIdx, compIdx)
403 +
404 intQuantsOutside.effectiveDiffusionCoefficient(phaseIdx, compIdx))
405 / 2;
406
407 Opm::Valgrind::CheckDefined(effectiveDiffusionCoefficient_[phaseIdx][compIdx]);
408 }
409 }
410 }
411
412 template <class Context, class FluidState>
413 void updateBoundary_(const Context& context,
414 unsigned bfIdx,
415 unsigned timeIdx,
416 const FluidState& fluidState)
417 {
418 const auto& stencil = context.stencil(timeIdx);
419 const auto& face = stencil.boundaryFace(bfIdx);
420
421 const auto& elemCtx = context.elementContext();
422 unsigned insideScvIdx = face.interiorIndex();
423 const auto& insideScv = stencil.subControlVolume(insideScvIdx);
424
425 const auto& intQuantsInside = elemCtx.intensiveQuantities(insideScvIdx, timeIdx);
426 const auto& fluidStateInside = intQuantsInside.fluidState();
427
428 // distance between the center of the SCV and center of the boundary face
429 DimVector distVec = face.integrationPos();
430 distVec -= context.element().geometry().global(insideScv.localGeometry().center());
431
432 Scalar dist = distVec * face.normal();
433
434 // if the following assertation triggers, the center of the
435 // center of the interior SCV was not inside the element!
436 assert(dist > 0);
437
438 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
439 if (!elemCtx.model().phaseIsConsidered(phaseIdx))
440 continue;
441
442 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
443 // calculate mole fraction gradient using two-point
444 // gradients
445 moleFractionGradientNormal_[phaseIdx][compIdx] =
446 (fluidState.moleFraction(phaseIdx, compIdx)
447 -
448 fluidStateInside.moleFraction(phaseIdx, compIdx))
449 / dist;
450 Opm::Valgrind::CheckDefined(moleFractionGradientNormal_[phaseIdx][compIdx]);
451
452 // use effective diffusion coefficients of the interior finite
453 // volume.
454 effectiveDiffusionCoefficient_[phaseIdx][compIdx] =
455 intQuantsInside.effectiveDiffusionCoefficient(phaseIdx, compIdx);
456
457 Opm::Valgrind::CheckDefined(effectiveDiffusionCoefficient_[phaseIdx][compIdx]);
458 }
459 }
460 }
461
462public:
469 const Evaluation& moleFractionGradientNormal(unsigned phaseIdx, unsigned compIdx) const
470 { return moleFractionGradientNormal_[phaseIdx][compIdx]; }
471
479 const Evaluation& effectiveDiffusionCoefficient(unsigned phaseIdx, unsigned compIdx) const
480 { return effectiveDiffusionCoefficient_[phaseIdx][compIdx]; }
481
482private:
483 Evaluation moleFractionGradientNormal_[numPhases][numComponents];
484 Evaluation effectiveDiffusionCoefficient_[numPhases][numComponents];
485};
486
487} // namespace Opm
488
489#endif
void update_(const ElementContext &, unsigned, unsigned)
Update the quantities required to calculate the diffusive mass fluxes.
Definition: diffusionmodule.hh:303
const Evaluation & effectiveDiffusionCoefficient(unsigned, unsigned) const
The effective diffusion coeffcient of a component in a fluid phase at the face's integration point.
Definition: diffusionmodule.hh:336
const Evaluation & moleFractionGradientNormal(unsigned, unsigned) const
The the gradient of the mole fraction times the face normal.
Definition: diffusionmodule.hh:322
void updateBoundary_(const Context &, unsigned, unsigned, const FluidState &)
Definition: diffusionmodule.hh:309
void updateBoundary_(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Definition: diffusionmodule.hh:413
const Evaluation & moleFractionGradientNormal(unsigned phaseIdx, unsigned compIdx) const
The the gradient of the mole fraction times the face normal.
Definition: diffusionmodule.hh:469
const Evaluation & effectiveDiffusionCoefficient(unsigned phaseIdx, unsigned compIdx) const
The effective diffusion coeffcient of a component in a fluid phase at the face's integration point.
Definition: diffusionmodule.hh:479
void update_(const ElementContext &elemCtx, unsigned faceIdx, unsigned timeIdx)
Update the quantities required to calculate the diffusive mass fluxes.
Definition: diffusionmodule.hh:367
Provides the quantities required to calculate diffusive mass fluxes.
Definition: diffusionmodule.hh:286
Scalar diffusionCoefficient(unsigned, unsigned) const
Returns the molecular diffusion coefficient for a component in a phase.
Definition: diffusionmodule.hh:167
Scalar effectiveDiffusionCoefficient(unsigned, unsigned) const
Returns the effective molecular diffusion coefficient of the porous medium for a component in a phase...
Definition: diffusionmodule.hh:177
void update_(FluidState &, typename FluidSystem::template ParameterCache< typename FluidState::Scalar > &, const ElementContext &, unsigned, unsigned)
Update the quantities required to calculate diffusive mass fluxes.
Definition: diffusionmodule.hh:189
Scalar tortuosity(unsigned) const
Returns the tortuousity of the sub-domain of a fluid phase in the porous medium.
Definition: diffusionmodule.hh:157
Evaluation tortuosity(unsigned phaseIdx) const
Returns the tortuousity of the sub-domain of a fluid phase in the porous medium.
Definition: diffusionmodule.hh:223
void update_(FluidState &fluidState, typename FluidSystem::template ParameterCache< typename FluidState::Scalar > &paramCache, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Update the quantities required to calculate diffusive mass fluxes.
Definition: diffusionmodule.hh:239
Evaluation effectiveDiffusionCoefficient(unsigned phaseIdx, unsigned compIdx) const
Returns the effective molecular diffusion coefficient of the porous medium for a component in a phase...
Definition: diffusionmodule.hh:230
Evaluation diffusionCoefficient(unsigned phaseIdx, unsigned compIdx) const
Returns the molecular diffusion coefficient for a component in a phase.
Definition: diffusionmodule.hh:216
Provides the volumetric quantities required for the calculation of molecular diffusive fluxes.
Definition: diffusionmodule.hh:140
static void addDiffusiveFlux(RateVector &, const Context &, unsigned, unsigned)
Adds the diffusive mass flux flux to the flux vector over a flux integration point.
Definition: diffusionmodule.hh:71
static void registerParameters()
Register all run-time parameters for the diffusion module.
Definition: diffusionmodule.hh:63
static void addDiffusiveFlux(RateVector &flux, const Context &context, unsigned spaceIdx, unsigned timeIdx)
Adds the mass flux due to molecular diffusion to the flux vector over the flux integration point.
Definition: diffusionmodule.hh:108
static void registerParameters()
Register all run-time parameters for the diffusion module.
Definition: diffusionmodule.hh:100
Provides the auxiliary methods required for consideration of the diffusion equation.
Definition: diffusionmodule.hh:47
Callback class for a mole fraction of a component in a phase.
Definition: quantitycallbacks.hh:425
void setComponentIndex(unsigned compIdx)
Set the index of the component for which the mole fraction should be returned.
Definition: quantitycallbacks.hh:460
void setPhaseIndex(unsigned phaseIdx)
Set the index of the fluid phase for which a mole fraction should be returned.
Definition: quantitycallbacks.hh:453
Declare the properties used by the infrastructure code of the finite volume discretizations.
Definition: blackoilboundaryratevector.hh:37
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:242
This method contains all callback classes for quantities that are required by some extensive quantiti...