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
39#include <array>
40#include <stdexcept>
41
42namespace Opm {
43
50template <class TypeTag, bool enableDiffusion>
52
56template <class TypeTag>
57class DiffusionModule<TypeTag, /*enableDiffusion=*/false>
58{
60
61public:
65 static void registerParameters()
66 {}
67
72 template <class Context>
73 static void addDiffusiveFlux(RateVector&,
74 const Context&,
75 unsigned,
76 unsigned)
77 {}
78};
79
83template <class TypeTag>
84class DiffusionModule<TypeTag, /*enableDiffusion=*/true>
85{
91
92 enum { numPhases = FluidSystem::numPhases };
93 enum { numComponents = FluidSystem::numComponents };
94 enum { conti0EqIdx = Indices::conti0EqIdx };
95
96 using Toolbox = Opm::MathToolbox<Evaluation>;
97
98public:
102 static void registerParameters()
103 {}
104
109 template <class Context>
110 static void addDiffusiveFlux(RateVector& flux, const Context& context,
111 unsigned spaceIdx, unsigned timeIdx)
112 {
113 const auto& extQuants = context.extensiveQuantities(spaceIdx, timeIdx);
114
115 const auto& fluidStateI = context.intensiveQuantities(extQuants.interiorIndex(), timeIdx).fluidState();
116 const auto& fluidStateJ = context.intensiveQuantities(extQuants.exteriorIndex(), timeIdx).fluidState();
117
118 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
119 // arithmetic mean of the phase's molar density
120 Evaluation rhoMolar = fluidStateI.molarDensity(phaseIdx);
121 rhoMolar += Toolbox::value(fluidStateJ.molarDensity(phaseIdx));
122 rhoMolar /= 2;
123
124 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
125 // mass flux due to molecular diffusion
126 flux[conti0EqIdx + compIdx] +=
127 -rhoMolar
128 * extQuants.moleFractionGradientNormal(phaseIdx, compIdx)
129 * extQuants.effectiveDiffusionCoefficient(phaseIdx, compIdx);
130 }
131 }
132 }
133};
134
142template <class TypeTag, bool enableDiffusion>
144
148template <class TypeTag>
149class DiffusionIntensiveQuantities<TypeTag, /*enableDiffusion=*/false>
150{
154
155public:
160 Scalar tortuosity(unsigned) const
161 {
162 throw std::logic_error("Method tortuosity() does not make sense "
163 "if diffusion is disabled");
164 }
165
170 Scalar diffusionCoefficient(unsigned, unsigned) const
171 {
172 throw std::logic_error("Method diffusionCoefficient() does not "
173 "make sense if diffusion is disabled");
174 }
175
180 Scalar effectiveDiffusionCoefficient(unsigned, unsigned) const
181 {
182 throw std::logic_error("Method effectiveDiffusionCoefficient() "
183 "does not make sense if diffusion is disabled");
184 }
185
186protected:
191 template <class FluidState>
192 void update_(FluidState&,
193 typename FluidSystem::template ParameterCache<typename FluidState::Scalar>&,
194 const ElementContext&,
195 unsigned,
196 unsigned)
197 { }
198};
199
203template <class TypeTag>
204class DiffusionIntensiveQuantities<TypeTag, /*enableDiffusion=*/true>
205{
210
211 enum { numPhases = FluidSystem::numPhases };
212 enum { numComponents = FluidSystem::numComponents };
213
214public:
219 Evaluation diffusionCoefficient(unsigned phaseIdx, unsigned compIdx) const
220 { return diffusionCoefficient_[phaseIdx][compIdx]; }
221
226 Evaluation tortuosity(unsigned phaseIdx) const
227 { return tortuosity_[phaseIdx]; }
228
233 Evaluation effectiveDiffusionCoefficient(unsigned phaseIdx, unsigned compIdx) const
234 { return tortuosity_[phaseIdx] * diffusionCoefficient_[phaseIdx][compIdx]; }
235
236protected:
241 template <class FluidState>
242 void update_(FluidState& fluidState,
243 typename FluidSystem::template ParameterCache<typename FluidState::Scalar>& paramCache,
244 const ElementContext& elemCtx,
245 unsigned dofIdx,
246 unsigned timeIdx)
247 {
248 using Toolbox = Opm::MathToolbox<Evaluation>;
249
250 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
251 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
252 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
253 continue;
254 }
255
256 // TODO: let the problem do this (this is a constitutive
257 // relation of which the model should be free of from the
258 // abstraction POV!)
259 // Based on Millington, R. J., & Quirk, J. P. (1961).
260 const Evaluation& base =
261 Toolbox::max(0.0001,
262 intQuants.porosity()
263 * intQuants.fluidState().saturation(phaseIdx));
264 tortuosity_[phaseIdx] =
265 1.0 / (intQuants.porosity() * intQuants.porosity())
266 * Toolbox::pow(base, 10.0/3.0);
267
268 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
269 diffusionCoefficient_[phaseIdx][compIdx] =
270 FluidSystem::diffusionCoefficient(fluidState,
271 paramCache,
272 phaseIdx,
273 compIdx);
274 }
275 }
276 }
277
278private:
279 std::array<Evaluation, numPhases> tortuosity_;
280 std::array<std::array<Evaluation, numComponents>, numPhases> diffusionCoefficient_;
281};
282
289template <class TypeTag, bool enableDiffusion>
291
295template <class TypeTag>
296class DiffusionExtensiveQuantities<TypeTag, /*enableDiffusion=*/false>
297{
300
301protected:
306 void update_(const ElementContext&,
307 unsigned,
308 unsigned)
309 {}
310
311 template <class Context, class FluidState>
312 void updateBoundary_(const Context&,
313 unsigned,
314 unsigned,
315 const FluidState&)
316 {}
317
318public:
325 const Evaluation& moleFractionGradientNormal(unsigned,
326 unsigned) const
327 {
328 throw std::logic_error("The method moleFractionGradient() does not "
329 "make sense if diffusion is disabled.");
330 }
331
339 const Evaluation& effectiveDiffusionCoefficient(unsigned,
340 unsigned) const
341 {
342 throw std::logic_error("The method effectiveDiffusionCoefficient() "
343 "does not make sense if diffusion is disabled.");
344 }
345};
346
350template <class TypeTag>
351class DiffusionExtensiveQuantities<TypeTag, /*enableDiffusion=*/true>
352{
357
358 enum { dimWorld = GridView::dimensionworld };
359 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
360 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
361
362 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
363 using DimEvalVector = Dune::FieldVector<Evaluation, dimWorld>;
364
365protected:
370 void update_(const ElementContext& elemCtx, unsigned faceIdx, unsigned timeIdx)
371 {
372 const auto& gradCalc = elemCtx.gradientCalculator();
373 Opm::MoleFractionCallback<TypeTag> moleFractionCallback(elemCtx);
374
375 const auto& face = elemCtx.stencil(timeIdx).interiorFace(faceIdx);
376 const auto& normal = face.normal();
377 const auto& extQuants = elemCtx.extensiveQuantities(faceIdx, timeIdx);
378
379 const auto& intQuantsInside = elemCtx.intensiveQuantities(extQuants.interiorIndex(), timeIdx);
380 const auto& intQuantsOutside = elemCtx.intensiveQuantities(extQuants.exteriorIndex(), timeIdx);
381
382 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
383 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
384 continue;
385 }
386
387 moleFractionCallback.setPhaseIndex(phaseIdx);
388 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
389 moleFractionCallback.setComponentIndex(compIdx);
390
391 DimEvalVector moleFractionGradient(0.0);
392 gradCalc.calculateGradient(moleFractionGradient,
393 elemCtx,
394 faceIdx,
395 moleFractionCallback);
396
397 moleFractionGradientNormal_[phaseIdx][compIdx] = 0.0;
398 for (unsigned i = 0; i < normal.size(); ++i) {
399 moleFractionGradientNormal_[phaseIdx][compIdx] +=
400 normal[i]*moleFractionGradient[i];
401 }
402 Opm::Valgrind::CheckDefined(moleFractionGradientNormal_[phaseIdx][compIdx]);
403
404 // use the arithmetic average for the effective
405 // diffusion coefficients.
406 effectiveDiffusionCoefficient_[phaseIdx][compIdx] =
407 (intQuantsInside.effectiveDiffusionCoefficient(phaseIdx, compIdx)
408 +
409 intQuantsOutside.effectiveDiffusionCoefficient(phaseIdx, compIdx))
410 / 2;
411
412 Opm::Valgrind::CheckDefined(effectiveDiffusionCoefficient_[phaseIdx][compIdx]);
413 }
414 }
415 }
416
417 template <class Context, class FluidState>
418 void updateBoundary_(const Context& context,
419 unsigned bfIdx,
420 unsigned timeIdx,
421 const FluidState& fluidState)
422 {
423 const auto& stencil = context.stencil(timeIdx);
424 const auto& face = stencil.boundaryFace(bfIdx);
425
426 const auto& elemCtx = context.elementContext();
427 const unsigned insideScvIdx = face.interiorIndex();
428 const auto& insideScv = stencil.subControlVolume(insideScvIdx);
429
430 const auto& intQuantsInside = elemCtx.intensiveQuantities(insideScvIdx, timeIdx);
431 const auto& fluidStateInside = intQuantsInside.fluidState();
432
433 // distance between the center of the SCV and center of the boundary face
434 DimVector distVec = face.integrationPos();
435 distVec -= context.element().geometry().global(insideScv.localGeometry().center());
436
437 const Scalar dist = distVec * face.normal();
438
439 // if the following assertation triggers, the center of the
440 // center of the interior SCV was not inside the element!
441 assert(dist > 0);
442
443 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
444 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
445 continue;
446 }
447
448 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
449 // calculate mole fraction gradient using two-point
450 // gradients
451 moleFractionGradientNormal_[phaseIdx][compIdx] =
452 (fluidState.moleFraction(phaseIdx, compIdx)
453 -
454 fluidStateInside.moleFraction(phaseIdx, compIdx))
455 / dist;
456 Opm::Valgrind::CheckDefined(moleFractionGradientNormal_[phaseIdx][compIdx]);
457
458 // use effective diffusion coefficients of the interior finite
459 // volume.
460 effectiveDiffusionCoefficient_[phaseIdx][compIdx] =
461 intQuantsInside.effectiveDiffusionCoefficient(phaseIdx, compIdx);
462
463 Opm::Valgrind::CheckDefined(effectiveDiffusionCoefficient_[phaseIdx][compIdx]);
464 }
465 }
466 }
467
468public:
475 const Evaluation& moleFractionGradientNormal(unsigned phaseIdx, unsigned compIdx) const
476 { return moleFractionGradientNormal_[phaseIdx][compIdx]; }
477
485 const Evaluation& effectiveDiffusionCoefficient(unsigned phaseIdx, unsigned compIdx) const
486 { return effectiveDiffusionCoefficient_[phaseIdx][compIdx]; }
487
488private:
489 std::array<std::array<Evaluation, numComponents>, numPhases> moleFractionGradientNormal_;
490 std::array<std::array<Evaluation, numComponents>, numPhases> effectiveDiffusionCoefficient_;
491};
492
493} // namespace Opm
494
495#endif
void update_(const ElementContext &, unsigned, unsigned)
Update the quantities required to calculate the diffusive mass fluxes.
Definition: diffusionmodule.hh:306
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:339
const Evaluation & moleFractionGradientNormal(unsigned, unsigned) const
The the gradient of the mole fraction times the face normal.
Definition: diffusionmodule.hh:325
void updateBoundary_(const Context &, unsigned, unsigned, const FluidState &)
Definition: diffusionmodule.hh:312
void updateBoundary_(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Definition: diffusionmodule.hh:418
const Evaluation & moleFractionGradientNormal(unsigned phaseIdx, unsigned compIdx) const
The the gradient of the mole fraction times the face normal.
Definition: diffusionmodule.hh:475
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:485
void update_(const ElementContext &elemCtx, unsigned faceIdx, unsigned timeIdx)
Update the quantities required to calculate the diffusive mass fluxes.
Definition: diffusionmodule.hh:370
Provides the quantities required to calculate diffusive mass fluxes.
Definition: diffusionmodule.hh:290
Scalar diffusionCoefficient(unsigned, unsigned) const
Returns the molecular diffusion coefficient for a component in a phase.
Definition: diffusionmodule.hh:170
Scalar effectiveDiffusionCoefficient(unsigned, unsigned) const
Returns the effective molecular diffusion coefficient of the porous medium for a component in a phase...
Definition: diffusionmodule.hh:180
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:192
Scalar tortuosity(unsigned) const
Returns the tortuousity of the sub-domain of a fluid phase in the porous medium.
Definition: diffusionmodule.hh:160
Evaluation tortuosity(unsigned phaseIdx) const
Returns the tortuousity of the sub-domain of a fluid phase in the porous medium.
Definition: diffusionmodule.hh:226
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:242
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:233
Evaluation diffusionCoefficient(unsigned phaseIdx, unsigned compIdx) const
Returns the molecular diffusion coefficient for a component in a phase.
Definition: diffusionmodule.hh:219
Provides the volumetric quantities required for the calculation of molecular diffusive fluxes.
Definition: diffusionmodule.hh:143
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:73
static void registerParameters()
Register all run-time parameters for the diffusion module.
Definition: diffusionmodule.hh:65
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:110
static void registerParameters()
Register all run-time parameters for the diffusion module.
Definition: diffusionmodule.hh:102
Provides the auxiliary methods required for consideration of the diffusion equation.
Definition: diffusionmodule.hh:51
Callback class for a mole fraction of a component in a phase.
Definition: quantitycallbacks.hh:422
void setComponentIndex(unsigned compIdx)
Set the index of the component for which the mole fraction should be returned.
Definition: quantitycallbacks.hh:457
void setPhaseIndex(unsigned phaseIdx)
Set the index of the fluid phase for which a mole fraction should be returned.
Definition: quantitycallbacks.hh:450
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: blackoilbioeffectsmodules.hh:43
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:233
This method contains all callback classes for quantities that are required by some extensive quantiti...