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
31#include <dune/common/fvector.hh>
32
33#include <opm/material/common/Valgrind.hpp>
34
38
39namespace Opm {
40
47template <class TypeTag, bool enableDiffusion>
49
53template <class TypeTag>
54class DiffusionModule<TypeTag, /*enableDiffusion=*/false>
55{
59
60public:
64 static void registerParameters()
65 {}
66
71 template <class Context>
72 static void addDiffusiveFlux(RateVector&,
73 const Context&,
74 unsigned,
75 unsigned)
76 {}
77};
78
82template <class TypeTag>
83class DiffusionModule<TypeTag, /*enableDiffusion=*/true>
84{
90
91 enum { numPhases = FluidSystem::numPhases };
92 enum { numComponents = FluidSystem::numComponents };
93 enum { conti0EqIdx = Indices::conti0EqIdx };
94
95 using Toolbox = Opm::MathToolbox<Evaluation>;
96
97public:
101 static void registerParameters()
102 {}
103
108 template <class Context>
109 static void addDiffusiveFlux(RateVector& flux, const Context& context,
110 unsigned spaceIdx, unsigned timeIdx)
111 {
112 const auto& extQuants = context.extensiveQuantities(spaceIdx, timeIdx);
113
114 const auto& fluidStateI = context.intensiveQuantities(extQuants.interiorIndex(), timeIdx).fluidState();
115 const auto& fluidStateJ = context.intensiveQuantities(extQuants.exteriorIndex(), timeIdx).fluidState();
116
117 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
118 // arithmetic mean of the phase's molar density
119 Evaluation rhoMolar = fluidStateI.molarDensity(phaseIdx);
120 rhoMolar += Toolbox::value(fluidStateJ.molarDensity(phaseIdx));
121 rhoMolar /= 2;
122
123 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
124 // mass flux due to molecular diffusion
125 flux[conti0EqIdx + compIdx] +=
126 -rhoMolar
127 * extQuants.moleFractionGradientNormal(phaseIdx, compIdx)
128 * extQuants.effectiveDiffusionCoefficient(phaseIdx, compIdx);
129 }
130 }
131};
132
140template <class TypeTag, bool enableDiffusion>
142
146template <class TypeTag>
147class DiffusionIntensiveQuantities<TypeTag, /*enableDiffusion=*/false>
148{
152
153public:
158 Scalar tortuosity(unsigned) const
159 {
160 throw std::logic_error("Method tortuosity() does not make sense "
161 "if diffusion is disabled");
162 }
163
168 Scalar diffusionCoefficient(unsigned, unsigned) const
169 {
170 throw std::logic_error("Method diffusionCoefficient() does not "
171 "make sense if diffusion is disabled");
172 }
173
178 Scalar effectiveDiffusionCoefficient(unsigned, unsigned) const
179 {
180 throw std::logic_error("Method effectiveDiffusionCoefficient() "
181 "does not make sense if diffusion is disabled");
182 }
183
184protected:
189 template <class FluidState>
190 void update_(FluidState&,
191 typename FluidSystem::template ParameterCache<typename FluidState::Scalar>&,
192 const ElementContext&,
193 unsigned,
194 unsigned)
195 { }
196};
197
201template <class TypeTag>
202class DiffusionIntensiveQuantities<TypeTag, /*enableDiffusion=*/true>
203{
208
209 enum { numPhases = FluidSystem::numPhases };
210 enum { numComponents = FluidSystem::numComponents };
211
212public:
217 Evaluation diffusionCoefficient(unsigned phaseIdx, unsigned compIdx) const
218 { return diffusionCoefficient_[phaseIdx][compIdx]; }
219
224 Evaluation tortuosity(unsigned phaseIdx) const
225 { return tortuosity_[phaseIdx]; }
226
231 Evaluation effectiveDiffusionCoefficient(unsigned phaseIdx, unsigned compIdx) const
232 { return tortuosity_[phaseIdx] * diffusionCoefficient_[phaseIdx][compIdx]; }
233
234protected:
239 template <class FluidState>
240 void update_(FluidState& fluidState,
241 typename FluidSystem::template ParameterCache<typename FluidState::Scalar>& paramCache,
242 const ElementContext& elemCtx,
243 unsigned dofIdx,
244 unsigned timeIdx)
245 {
246 using Toolbox = Opm::MathToolbox<Evaluation>;
247
248 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
249 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
250 if (!elemCtx.model().phaseIsConsidered(phaseIdx))
251 continue;
252
253 // TODO: let the problem do this (this is a constitutive
254 // relation of which the model should be free of from the
255 // abstraction POV!)
256 // Based on Millington, R. J., & Quirk, J. P. (1961).
257 const Evaluation& base =
258 Toolbox::max(0.0001,
259 intQuants.porosity()
260 * intQuants.fluidState().saturation(phaseIdx));
261 tortuosity_[phaseIdx] =
262 1.0 / (intQuants.porosity() * intQuants.porosity())
263 * Toolbox::pow(base, 10.0/3.0);
264
265 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
266 diffusionCoefficient_[phaseIdx][compIdx] =
267 FluidSystem::diffusionCoefficient(fluidState,
268 paramCache,
269 phaseIdx,
270 compIdx);
271 }
272 }
273 }
274
275private:
276 Evaluation tortuosity_[numPhases];
277 Evaluation diffusionCoefficient_[numPhases][numComponents];
278};
279
286template <class TypeTag, bool enableDiffusion>
288
292template <class TypeTag>
293class DiffusionExtensiveQuantities<TypeTag, /*enableDiffusion=*/false>
294{
298
299protected:
304 void update_(const ElementContext&,
305 unsigned,
306 unsigned)
307 {}
308
309 template <class Context, class FluidState>
310 void updateBoundary_(const Context&,
311 unsigned,
312 unsigned,
313 const FluidState&)
314 {}
315
316public:
323 const Evaluation& moleFractionGradientNormal(unsigned,
324 unsigned) const
325 {
326 throw std::logic_error("The method moleFractionGradient() does not "
327 "make sense if diffusion is disabled.");
328 }
329
337 const Evaluation& effectiveDiffusionCoefficient(unsigned,
338 unsigned) const
339 {
340 throw std::logic_error("The method effectiveDiffusionCoefficient() "
341 "does not make sense if diffusion is disabled.");
342 }
343};
344
348template <class TypeTag>
349class DiffusionExtensiveQuantities<TypeTag, /*enableDiffusion=*/true>
350{
355
356 enum { dimWorld = GridView::dimensionworld };
357 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
358 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
359
360 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
361 using DimEvalVector = Dune::FieldVector<Evaluation, dimWorld>;
362
363protected:
368 void update_(const ElementContext& elemCtx, unsigned faceIdx, unsigned timeIdx)
369 {
370 const auto& gradCalc = elemCtx.gradientCalculator();
371 Opm::MoleFractionCallback<TypeTag> moleFractionCallback(elemCtx);
372
373 const auto& face = elemCtx.stencil(timeIdx).interiorFace(faceIdx);
374 const auto& normal = face.normal();
375 const auto& extQuants = elemCtx.extensiveQuantities(faceIdx, timeIdx);
376
377 const auto& intQuantsInside = elemCtx.intensiveQuantities(extQuants.interiorIndex(), timeIdx);
378 const auto& intQuantsOutside = elemCtx.intensiveQuantities(extQuants.exteriorIndex(), timeIdx);
379
380 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
381 if (!elemCtx.model().phaseIsConsidered(phaseIdx))
382 continue;
383
384 moleFractionCallback.setPhaseIndex(phaseIdx);
385 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
386 moleFractionCallback.setComponentIndex(compIdx);
387
388 DimEvalVector moleFractionGradient(0.0);
389 gradCalc.calculateGradient(moleFractionGradient,
390 elemCtx,
391 faceIdx,
392 moleFractionCallback);
393
394 moleFractionGradientNormal_[phaseIdx][compIdx] = 0.0;
395 for (unsigned i = 0; i < normal.size(); ++i)
396 moleFractionGradientNormal_[phaseIdx][compIdx] +=
397 normal[i]*moleFractionGradient[i];
398 Opm::Valgrind::CheckDefined(moleFractionGradientNormal_[phaseIdx][compIdx]);
399
400 // use the arithmetic average for the effective
401 // diffusion coefficients.
402 effectiveDiffusionCoefficient_[phaseIdx][compIdx] =
403 (intQuantsInside.effectiveDiffusionCoefficient(phaseIdx, compIdx)
404 +
405 intQuantsOutside.effectiveDiffusionCoefficient(phaseIdx, compIdx))
406 / 2;
407
408 Opm::Valgrind::CheckDefined(effectiveDiffusionCoefficient_[phaseIdx][compIdx]);
409 }
410 }
411 }
412
413 template <class Context, class FluidState>
414 void updateBoundary_(const Context& context,
415 unsigned bfIdx,
416 unsigned timeIdx,
417 const FluidState& fluidState)
418 {
419 const auto& stencil = context.stencil(timeIdx);
420 const auto& face = stencil.boundaryFace(bfIdx);
421
422 const auto& elemCtx = context.elementContext();
423 unsigned insideScvIdx = face.interiorIndex();
424 const auto& insideScv = stencil.subControlVolume(insideScvIdx);
425
426 const auto& intQuantsInside = elemCtx.intensiveQuantities(insideScvIdx, timeIdx);
427 const auto& fluidStateInside = intQuantsInside.fluidState();
428
429 // distance between the center of the SCV and center of the boundary face
430 DimVector distVec = face.integrationPos();
431 distVec -= context.element().geometry().global(insideScv.localGeometry().center());
432
433 Scalar dist = distVec * face.normal();
434
435 // if the following assertation triggers, the center of the
436 // center of the interior SCV was not inside the element!
437 assert(dist > 0);
438
439 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
440 if (!elemCtx.model().phaseIsConsidered(phaseIdx))
441 continue;
442
443 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
444 // calculate mole fraction gradient using two-point
445 // gradients
446 moleFractionGradientNormal_[phaseIdx][compIdx] =
447 (fluidState.moleFraction(phaseIdx, compIdx)
448 -
449 fluidStateInside.moleFraction(phaseIdx, compIdx))
450 / dist;
451 Opm::Valgrind::CheckDefined(moleFractionGradientNormal_[phaseIdx][compIdx]);
452
453 // use effective diffusion coefficients of the interior finite
454 // volume.
455 effectiveDiffusionCoefficient_[phaseIdx][compIdx] =
456 intQuantsInside.effectiveDiffusionCoefficient(phaseIdx, compIdx);
457
458 Opm::Valgrind::CheckDefined(effectiveDiffusionCoefficient_[phaseIdx][compIdx]);
459 }
460 }
461 }
462
463public:
470 const Evaluation& moleFractionGradientNormal(unsigned phaseIdx, unsigned compIdx) const
471 { return moleFractionGradientNormal_[phaseIdx][compIdx]; }
472
480 const Evaluation& effectiveDiffusionCoefficient(unsigned phaseIdx, unsigned compIdx) const
481 { return effectiveDiffusionCoefficient_[phaseIdx][compIdx]; }
482
483private:
484 Evaluation moleFractionGradientNormal_[numPhases][numComponents];
485 Evaluation effectiveDiffusionCoefficient_[numPhases][numComponents];
486};
487
488} // namespace Opm
489
490#endif
void update_(const ElementContext &, unsigned, unsigned)
Update the quantities required to calculate the diffusive mass fluxes.
Definition: diffusionmodule.hh:304
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:337
const Evaluation & moleFractionGradientNormal(unsigned, unsigned) const
The the gradient of the mole fraction times the face normal.
Definition: diffusionmodule.hh:323
void updateBoundary_(const Context &, unsigned, unsigned, const FluidState &)
Definition: diffusionmodule.hh:310
void updateBoundary_(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Definition: diffusionmodule.hh:414
const Evaluation & moleFractionGradientNormal(unsigned phaseIdx, unsigned compIdx) const
The the gradient of the mole fraction times the face normal.
Definition: diffusionmodule.hh:470
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:480
void update_(const ElementContext &elemCtx, unsigned faceIdx, unsigned timeIdx)
Update the quantities required to calculate the diffusive mass fluxes.
Definition: diffusionmodule.hh:368
Provides the quantities required to calculate diffusive mass fluxes.
Definition: diffusionmodule.hh:287
Scalar diffusionCoefficient(unsigned, unsigned) const
Returns the molecular diffusion coefficient for a component in a phase.
Definition: diffusionmodule.hh:168
Scalar effectiveDiffusionCoefficient(unsigned, unsigned) const
Returns the effective molecular diffusion coefficient of the porous medium for a component in a phase...
Definition: diffusionmodule.hh:178
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:190
Scalar tortuosity(unsigned) const
Returns the tortuousity of the sub-domain of a fluid phase in the porous medium.
Definition: diffusionmodule.hh:158
Evaluation tortuosity(unsigned phaseIdx) const
Returns the tortuousity of the sub-domain of a fluid phase in the porous medium.
Definition: diffusionmodule.hh:224
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:240
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:231
Evaluation diffusionCoefficient(unsigned phaseIdx, unsigned compIdx) const
Returns the molecular diffusion coefficient for a component in a phase.
Definition: diffusionmodule.hh:217
Provides the volumetric quantities required for the calculation of molecular diffusive fluxes.
Definition: diffusionmodule.hh:141
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:72
static void registerParameters()
Register all run-time parameters for the diffusion module.
Definition: diffusionmodule.hh:64
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:109
static void registerParameters()
Register all run-time parameters for the diffusion module.
Definition: diffusionmodule.hh:101
Provides the auxiliary methods required for consideration of the diffusion equation.
Definition: diffusionmodule.hh:48
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.
Defines the common properties required by the porous medium multi-phase models.
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:235
This method contains all callback classes for quantities that are required by some extensive quantiti...