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