transfluxmodule.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*/
31#ifndef EWOMS_TRANS_FLUX_MODULE_HH
32#define EWOMS_TRANS_FLUX_MODULE_HH
33
34#include <dune/common/fmatrix.hh>
35#include <dune/common/fvector.hh>
36
37#include <opm/material/common/MathToolbox.hpp>
38#include <opm/material/common/Valgrind.hpp>
39
44
45#include <array>
46#include <cassert>
47#include <cmath>
48#include <stdexcept>
49#include <type_traits>
50
51namespace Opm {
52
53template <class TypeTag>
54class TransIntensiveQuantities;
55
56template <class TypeTag>
57class TransExtensiveQuantities;
58
59template <class TypeTag>
60class TransBaseProblem;
61
65template <class TypeTag>
67{
71
75 static void registerParameters()
76 {}
77};
78
83template <class TypeTag>
85{};
86
90template <class TypeTag>
92{
94
95protected:
96 void update_(const ElementContext&, unsigned, unsigned)
97 {}
98};
99
103template <class TypeTag>
105{
107
115
116 enum { dimWorld = GridView::dimensionworld };
117 enum { numPhases = FluidSystem::numPhases };
118
119 typedef MathToolbox<Evaluation> Toolbox;
120 typedef Dune::FieldVector<Evaluation, dimWorld> EvalDimVector;
121 typedef Dune::FieldMatrix<Scalar, dimWorld, dimWorld> DimMatrix;
122
123public:
127 const DimMatrix& intrinsicPermeability() const
128 {
129 throw std::logic_error("The ECL transmissibility module does not "
130 "provide an explicit intrinsic permeability");
131 }
132
139 const EvalDimVector& potentialGrad(unsigned) const
140 {
141 throw std::logic_error("The ECL transmissibility module does not "
142 "provide explicit potential gradients");
143 }
144
151 const Evaluation& pressureDifference(unsigned phaseIdx) const
152 { return pressureDifference_[phaseIdx]; }
153
160 const EvalDimVector& filterVelocity(unsigned) const
161 {
162 throw std::logic_error("The ECL transmissibility module does not "
163 "provide explicit filter velocities");
164 }
165
175 const Evaluation& volumeFlux(unsigned phaseIdx) const
176 { return volumeFlux_[phaseIdx]; }
177
178protected:
186 unsigned upstreamIndex_(unsigned phaseIdx) const
187 {
188 assert(phaseIdx < numPhases);
189
190 return upIdx_[phaseIdx];
191 }
192
200 unsigned downstreamIndex_(unsigned phaseIdx) const
201 {
202 assert(phaseIdx < numPhases);
203
204 return dnIdx_[phaseIdx];
205 }
206
207 void updateSolvent(const ElementContext& elemCtx, unsigned scvfIdx, unsigned timeIdx)
208 { asImp_().updateVolumeFluxTrans(elemCtx, scvfIdx, timeIdx); }
209
210 void updatePolymer(const ElementContext& elemCtx, unsigned scvfIdx, unsigned timeIdx)
211 { asImp_().updateShearMultipliers(elemCtx, scvfIdx, timeIdx); }
212
216 void calculateGradients_(const ElementContext& elemCtx, unsigned scvfIdx, unsigned timeIdx)
217 {
218 Valgrind::SetUndefined(*this);
219
220 // only valid for element center finite volume discretization
221 static_assert(std::is_same_v<Discretization, EcfvDiscretization<TypeTag>>);
222
223 const auto& stencil = elemCtx.stencil(timeIdx);
224 const auto& scvf = stencil.interiorFace(scvfIdx);
225
226 interiorDofIdx_ = scvf.interiorIndex();
227 exteriorDofIdx_ = scvf.exteriorIndex();
228 assert(interiorDofIdx_ != exteriorDofIdx_);
229
230 const unsigned I = stencil.globalSpaceIndex(interiorDofIdx_);
231 const unsigned J = stencil.globalSpaceIndex(exteriorDofIdx_);
232
233 const Scalar trans = transmissibility_(elemCtx, scvfIdx, timeIdx);
234
235 // estimate the gravity correction: for performance reasons we use a simplified
236 // approach for this flux module that assumes that gravity is constant and always
237 // acts into the downwards direction. (i.e., no centrifuge experiments, sorry.)
238 const Scalar g = elemCtx.problem().gravity()[dimWorld - 1];
239
240 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx_, timeIdx);
241 const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx_, timeIdx);
242
243 const Scalar zIn = dofCenterDepth_(elemCtx, interiorDofIdx_, timeIdx);
244 const Scalar zEx = dofCenterDepth_(elemCtx, exteriorDofIdx_, timeIdx);
245 // the distances from the DOF's depths. (i.e., the additional depth of the
246 // exterior DOF)
247 const Scalar distZ = zIn - zEx;
248
249 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
250 if (!FluidSystem::phaseIsActive(phaseIdx)) {
251 continue;
252 }
253
254 // check shortcut: if the mobility of the phase is zero in the interior as
255 // well as the exterior DOF, we can skip looking at the phase.
256 if (intQuantsIn.mobility(phaseIdx) <= 0.0 &&
257 intQuantsEx.mobility(phaseIdx) <= 0.0)
258 {
259 upIdx_[phaseIdx] = interiorDofIdx_;
260 dnIdx_[phaseIdx] = exteriorDofIdx_;
261 pressureDifference_[phaseIdx] = 0.0;
262 volumeFlux_[phaseIdx] = 0.0;
263 continue;
264 }
265
266 // do the gravity correction: compute the hydrostatic pressure for the
267 // external at the depth of the internal one
268 const Evaluation& rhoIn = intQuantsIn.fluidState().density(phaseIdx);
269 const Scalar rhoEx = Toolbox::value(intQuantsEx.fluidState().density(phaseIdx));
270 const Evaluation rhoAvg = (rhoIn + rhoEx) / 2;
271
272 const Evaluation& pressureInterior = intQuantsIn.fluidState().pressure(phaseIdx);
273 Evaluation pressureExterior = Toolbox::value(intQuantsEx.fluidState().pressure(phaseIdx));
274
275 pressureExterior += rhoAvg * (distZ * g);
276
277 pressureDifference_[phaseIdx] = pressureExterior - pressureInterior;
278
279 // decide the upstream index for the phase. for this we make sure that the
280 // degree of freedom which is regarded upstream if both pressures are equal
281 // is always the same: if the pressure is equal, the DOF with the lower
282 // global index is regarded to be the upstream one.
283 if (pressureDifference_[phaseIdx] > 0.0) {
284 upIdx_[phaseIdx] = exteriorDofIdx_;
285 dnIdx_[phaseIdx] = interiorDofIdx_;
286 }
287 else if (pressureDifference_[phaseIdx] < 0.0) {
288 upIdx_[phaseIdx] = interiorDofIdx_;
289 dnIdx_[phaseIdx] = exteriorDofIdx_;
290 }
291 else {
292 // if the pressure difference is zero, we chose the DOF which has the
293 // larger volume associated to it as upstream DOF
294 const Scalar Vin = elemCtx.dofVolume(interiorDofIdx_, /*timeIdx=*/0);
295 const Scalar Vex = elemCtx.dofVolume(exteriorDofIdx_, /*timeIdx=*/0);
296 if (Vin > Vex) {
297 upIdx_[phaseIdx] = interiorDofIdx_;
298 dnIdx_[phaseIdx] = exteriorDofIdx_;
299 }
300 else if (Vin < Vex) {
301 upIdx_[phaseIdx] = exteriorDofIdx_;
302 dnIdx_[phaseIdx] = interiorDofIdx_;
303 }
304 else {
305 assert(Vin == Vex);
306 // if the volumes are also equal, we pick the DOF which exhibits the
307 // smaller global index
308 if (I < J) {
309 upIdx_[phaseIdx] = interiorDofIdx_;
310 dnIdx_[phaseIdx] = exteriorDofIdx_;
311 }
312 else {
313 upIdx_[phaseIdx] = exteriorDofIdx_;
314 dnIdx_[phaseIdx] = interiorDofIdx_;
315 }
316 }
317 }
318
319 // this is slightly hacky because in the automatic differentiation case, it
320 // only works for the element centered finite volume method. for ebos this
321 // does not matter, though.
322 const unsigned upstreamIdx = upstreamIndex_(phaseIdx);
323 const auto& up = elemCtx.intensiveQuantities(upstreamIdx, timeIdx);
324
325 if (upstreamIdx == interiorDofIdx_) {
326 volumeFlux_[phaseIdx] =
327 pressureDifference_[phaseIdx] * up.mobility(phaseIdx) * (-trans);
328 }
329 else {
330 volumeFlux_[phaseIdx] =
331 pressureDifference_[phaseIdx] * (Toolbox::value(up.mobility(phaseIdx)) * (-trans));
332 }
333 }
334 }
335
339 template <class FluidState>
340 void calculateBoundaryGradients_(const ElementContext& elemCtx,
341 unsigned scvfIdx,
342 unsigned timeIdx,
343 const FluidState& exFluidState)
344 {
345 const auto& stencil = elemCtx.stencil(timeIdx);
346 const auto& scvf = stencil.boundaryFace(scvfIdx);
347
348 interiorDofIdx_ = scvf.interiorIndex();
349
350 const Scalar trans = transmissibilityBoundary_(elemCtx, scvfIdx, timeIdx);
351
352 // estimate the gravity correction: for performance reasons we use a simplified
353 // approach for this flux module that assumes that gravity is constant and always
354 // acts into the downwards direction. (i.e., no centrifuge experiments, sorry.)
355 const Scalar g = elemCtx.problem().gravity()[dimWorld - 1];
356
357 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx_, timeIdx);
358
359 // this is quite hacky because the dune grid interface does not provide a
360 // cellCenterDepth() method (so we ask the problem to provide it). The "good"
361 // solution would be to take the Z coordinate of the element centroids, but since
362 // ECL seems to like to be inconsistent on that front, it needs to be done like
363 // here...
364 const Scalar zIn = dofCenterDepth_(elemCtx, interiorDofIdx_, timeIdx);
365 const Scalar zEx = scvf.integrationPos()[dimWorld - 1];
366
367 // the distances from the DOF's depths. (i.e., the additional depth of the
368 // exterior DOF)
369 const Scalar distZ = zIn - zEx;
370
371 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
372 if (!FluidSystem::phaseIsActive(phaseIdx)) {
373 continue;
374 }
375
376 // do the gravity correction: compute the hydrostatic pressure for the
377 // integration position
378 const Evaluation& rhoIn = intQuantsIn.fluidState().density(phaseIdx);
379 const auto& rhoEx = exFluidState.density(phaseIdx);
380 const Evaluation rhoAvg = (rhoIn + rhoEx) / 2;
381
382 const Evaluation& pressureInterior = intQuantsIn.fluidState().pressure(phaseIdx);
383 Evaluation pressureExterior = exFluidState.pressure(phaseIdx);
384 pressureExterior += rhoAvg * (distZ * g);
385
386 pressureDifference_[phaseIdx] = pressureExterior - pressureInterior;
387
388 // decide the upstream index for the phase. for this we make sure that the
389 // degree of freedom which is regarded upstream if both pressures are equal
390 // is always the same: if the pressure is equal, the DOF with the lower
391 // global index is regarded to be the upstream one.
392 if (pressureDifference_[phaseIdx] > 0.0) {
393 upIdx_[phaseIdx] = -1;
394 dnIdx_[phaseIdx] = interiorDofIdx_;
395 }
396 else {
397 upIdx_[phaseIdx] = interiorDofIdx_;
398 dnIdx_[phaseIdx] = -1;
399 }
400
401 const short upstreamIdx = upstreamIndex_(phaseIdx);
402 if (upstreamIdx == interiorDofIdx_) {
403 // this is slightly hacky because in the automatic differentiation case, it
404 // only works for the element centered finite volume method. for ebos this
405 // does not matter, though.
406 const auto& up = elemCtx.intensiveQuantities(upstreamIdx, timeIdx);
407
408 volumeFlux_[phaseIdx] =
409 pressureDifference_[phaseIdx] * up.mobility(phaseIdx) * (-trans);
410 }
411 else {
412 // compute the phase mobility using the material law parameters of the
413 // interior element. \todo {this could probably be done more efficiently}
414 const auto& matParams =
415 elemCtx.problem().materialLawParams(elemCtx,
416 interiorDofIdx_,
417 /*timeIdx=*/0);
418 std::array<typename FluidState::Scalar,numPhases> kr;
419 MaterialLaw::relativePermeabilities(kr, matParams, exFluidState);
420
421 const auto& mob = kr[phaseIdx] / exFluidState.viscosity(phaseIdx);
422 volumeFlux_[phaseIdx] =
423 pressureDifference_[phaseIdx] * mob * (-trans);
424 }
425 }
426 }
427
431 void calculateFluxes_(const ElementContext&, unsigned, unsigned)
432 {}
433
434 void calculateBoundaryFluxes_(const ElementContext&, unsigned, unsigned)
435 {}
436
437private:
438 Scalar transmissibility_(const ElementContext& elemCtx, unsigned scvfIdx, unsigned timeIdx) const
439 {
440 const auto& stencil = elemCtx.stencil(timeIdx);
441 const auto& face = stencil.interiorFace(scvfIdx);
442 const auto& interiorPos = stencil.subControlVolume(face.interiorIndex()).globalPos();
443 const auto& exteriorPos = stencil.subControlVolume(face.exteriorIndex()).globalPos();
444 const auto distVec0 = face.integrationPos() - interiorPos;
445 const auto distVec1 = face.integrationPos() - exteriorPos;
446 const Scalar ndotDistIn = std::abs(face.normal() * distVec0);
447 const Scalar ndotDistExt = std::abs(face.normal() * distVec1);
448
449 const Scalar distSquaredIn = distVec0 * distVec0;
450 const Scalar distSquaredExt = distVec1 * distVec1;
451 const auto& K0mat = elemCtx.problem().intrinsicPermeability(elemCtx, face.interiorIndex(), timeIdx);
452 const auto& K1mat = elemCtx.problem().intrinsicPermeability(elemCtx, face.exteriorIndex(), timeIdx);
453
454 // the permeability per definition aligns with the grid
455 // we only support diagonal permeability tensor
456 // and can therefore neglect off-diagonal values
457 int idx = 0;
458 Scalar val = 0.0;
459 for (unsigned i = 0; i < dimWorld; ++i) {
460 if (std::abs(face.normal()[i]) > val) {
461 val = std::abs(face.normal()[i]);
462 idx = i;
463 }
464 }
465 const Scalar& K0 = K0mat[idx][idx];
466 const Scalar& K1 = K1mat[idx][idx];
467 const Scalar T0 = K0 * ndotDistIn / distSquaredIn;
468 const Scalar T1 = K1 * ndotDistExt / distSquaredExt;
469 return T0 * T1 / (T0 + T1);
470 }
471 Scalar transmissibilityBoundary_(const ElementContext& elemCtx, unsigned scvfIdx, unsigned timeIdx) const
472 {
473 const auto& stencil = elemCtx.stencil(timeIdx);
474 const auto& face = stencil.interiorFace(scvfIdx);
475 const auto& interiorPos = stencil.subControlVolume(face.interiorIndex()).globalPos();
476 const auto distVec0 = face.integrationPos() - interiorPos;
477 const Scalar ndotDistIn = face.normal() * distVec0;
478 const Scalar distSquaredIn = distVec0 * distVec0;
479 const auto& K0mat = elemCtx.problem().intrinsicPermeability(elemCtx, face.interiorIndex(), timeIdx);
480
481 // the permeability per definition aligns with the grid
482 // we only support diagonal permeability tensor
483 // and can therefore neglect off-diagonal values
484 int idx = 0;
485 Scalar val = 0.0;
486 for (unsigned i = 0; i < dimWorld; ++i) {
487 if (std::abs(face.normal()[i]) > val) {
488 val = std::abs(face.normal()[i]);
489 idx = i;
490 }
491 }
492 const Scalar& K0 = K0mat[idx][idx];
493 const Scalar T0 = K0 * ndotDistIn / distSquaredIn;
494 return T0;
495 }
496
497 template <class Context>
498 Scalar dofCenterDepth_(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
499 {
500 const auto& pos = context.pos(spaceIdx, timeIdx);
501 return pos[dimWorld-1];
502 }
503
504 Implementation& asImp_()
505 { return *static_cast<Implementation*>(this); }
506
507 const Implementation& asImp_() const
508 { return *static_cast<const Implementation*>(this); }
509
510 // the volumetric flux of all phases [m^3/s]
511 std::array<Evaluation, numPhases> volumeFlux_;
512
513 // the difference in effective pressure between the exterior and the interior degree
514 // of freedom [Pa]
515 std::array<Evaluation, numPhases> pressureDifference_;
516
517 // the local indices of the interior and exterior degrees of freedom
518 unsigned short interiorDofIdx_{};
519 unsigned short exteriorDofIdx_{};
520 std::array<short, numPhases> upIdx_{};
521 std::array<short, numPhases> dnIdx_{};
522};
523
524} // namespace Opm
525
526#endif
Provides the defaults for the parameters required by the transmissibility based volume flux calculati...
Definition: transfluxmodule.hh:85
Provides the transmissibility based flux module.
Definition: transfluxmodule.hh:105
const Evaluation & volumeFlux(unsigned phaseIdx) const
Return the volume flux of a fluid phase at the face's integration point .
Definition: transfluxmodule.hh:175
unsigned downstreamIndex_(unsigned phaseIdx) const
Returns the local index of the degree of freedom in which is in downstream direction.
Definition: transfluxmodule.hh:200
const EvalDimVector & filterVelocity(unsigned) const
Return the filter velocity of a fluid phase at the face's integration point [m/s].
Definition: transfluxmodule.hh:160
void calculateGradients_(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Update the required gradients for interior faces.
Definition: transfluxmodule.hh:216
void updateSolvent(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: transfluxmodule.hh:207
const Evaluation & pressureDifference(unsigned phaseIdx) const
Return the gravity corrected pressure difference between the interior and the exterior of a face.
Definition: transfluxmodule.hh:151
void calculateFluxes_(const ElementContext &, unsigned, unsigned)
Update the volumetric fluxes for all fluid phases on the interior faces of the context.
Definition: transfluxmodule.hh:431
void calculateBoundaryFluxes_(const ElementContext &, unsigned, unsigned)
Definition: transfluxmodule.hh:434
const EvalDimVector & potentialGrad(unsigned) const
Return the pressure potential gradient of a fluid phase at the face's integration point [Pa/m].
Definition: transfluxmodule.hh:139
void calculateBoundaryGradients_(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx, const FluidState &exFluidState)
Update the required gradients for boundary faces.
Definition: transfluxmodule.hh:340
const DimMatrix & intrinsicPermeability() const
Return the intrinsic permeability tensor at a face [m^2].
Definition: transfluxmodule.hh:127
unsigned upstreamIndex_(unsigned phaseIdx) const
Returns the local index of the degree of freedom in which is in upstream direction.
Definition: transfluxmodule.hh:186
void updatePolymer(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: transfluxmodule.hh:210
Provides the intensive quantities for the transmissibility based flux module.
Definition: transfluxmodule.hh:92
void update_(const ElementContext &, unsigned, unsigned)
Definition: transfluxmodule.hh:96
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
Specifies a flux module which uses transmissibilities.
Definition: transfluxmodule.hh:67
static void registerParameters()
Register all run-time parameters for the flux module.
Definition: transfluxmodule.hh:75