darcyfluxmodule.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*/
30#ifndef EWOMS_DARCY_FLUX_MODULE_HH
31#define EWOMS_DARCY_FLUX_MODULE_HH
32
33#include <dune/common/fmatrix.hh>
34#include <dune/common/fvector.hh>
35
36#include <opm/common/Exceptions.hpp>
37
38#include <opm/material/common/Valgrind.hpp>
39
43
44#include <array>
45#include <cassert>
46#include <cmath>
47#include <string>
48#include <type_traits>
49
50namespace Opm {
51
52template <class TypeTag>
53class DarcyIntensiveQuantities;
54
55template <class TypeTag>
56class DarcyExtensiveQuantities;
57
58template <class TypeTag>
59class DarcyBaseProblem;
60
65template <class TypeTag>
67{
71
75 static void registerParameters()
76 {}
77};
78
84template <class TypeTag>
86{};
87
92template <class TypeTag>
94{
96
97protected:
98 void update_(const ElementContext&,
99 unsigned,
100 unsigned)
101 {}
102};
103
120template <class TypeTag>
122{
130
131 enum { dimWorld = GridView::dimensionworld };
132 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
133
134 using Toolbox = MathToolbox<Evaluation>;
135 using EvalDimVector = Dune::FieldVector<Evaluation, dimWorld>;
136 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
137 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
138
139public:
144 const DimMatrix& intrinsicPermability() const
145 { return K_; }
146
153 const EvalDimVector& potentialGrad(unsigned phaseIdx) const
154 { return potentialGrad_[phaseIdx]; }
155
162 const EvalDimVector& filterVelocity(unsigned phaseIdx) const
163 { return filterVelocity_[phaseIdx]; }
164
174 const Evaluation& volumeFlux(unsigned phaseIdx) const
175 { return volumeFlux_[phaseIdx]; }
176
177protected:
178 short upstreamIndex_(unsigned phaseIdx) const
179 { return upstreamDofIdx_[phaseIdx]; }
180
181 short downstreamIndex_(unsigned phaseIdx) const
182 { return downstreamDofIdx_[phaseIdx]; }
183
189 void calculateGradients_(const ElementContext& elemCtx,
190 unsigned faceIdx,
191 unsigned timeIdx)
192 {
193 const auto& gradCalc = elemCtx.gradientCalculator();
194 PressureCallback<TypeTag> pressureCallback(elemCtx);
195
196 const auto& scvf = elemCtx.stencil(timeIdx).interiorFace(faceIdx);
197 const auto& faceNormal = scvf.normal();
198
199 const unsigned i = scvf.interiorIndex();
200 const unsigned j = scvf.exteriorIndex();
201 interiorDofIdx_ = static_cast<short>(i);
202 exteriorDofIdx_ = static_cast<short>(j);
203 const unsigned focusDofIdx = elemCtx.focusDofIndex();
204
205 // calculate the "raw" pressure gradient
206 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
207 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
208 Valgrind::SetUndefined(potentialGrad_[phaseIdx]);
209 continue;
210 }
211
212 pressureCallback.setPhaseIndex(phaseIdx);
213 gradCalc.calculateGradient(potentialGrad_[phaseIdx],
214 elemCtx,
215 faceIdx,
216 pressureCallback);
217 Valgrind::CheckDefined(potentialGrad_[phaseIdx]);
218 }
219
220 // correct the pressure gradients by the gravitational acceleration
221 if (Parameters::Get<Parameters::EnableGravity>()) {
222 // estimate the gravitational acceleration at a given SCV face
223 // using the arithmetic mean
224 const auto& gIn = elemCtx.problem().gravity(elemCtx, i, timeIdx);
225 const auto& gEx = elemCtx.problem().gravity(elemCtx, j, timeIdx);
226
227 const auto& intQuantsIn = elemCtx.intensiveQuantities(i, timeIdx);
228 const auto& intQuantsEx = elemCtx.intensiveQuantities(j, timeIdx);
229
230 const auto& posIn = elemCtx.pos(i, timeIdx);
231 const auto& posEx = elemCtx.pos(j, timeIdx);
232 const auto& posFace = scvf.integrationPos();
233
234 // the distance between the centers of the control volumes
235 DimVector distVecIn(posIn);
236 DimVector distVecEx(posEx);
237 DimVector distVecTotal(posEx);
238
239 distVecIn -= posFace;
240 distVecEx -= posFace;
241 distVecTotal -= posIn;
242 const Scalar absDistTotalSquared = distVecTotal.two_norm2();
243 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
244 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
245 continue;
246 }
247
248 // calculate the hydrostatic pressure at the integration point of the face
249 Evaluation pStatIn;
250
251 if (std::is_same<Scalar, Evaluation>::value ||
252 interiorDofIdx_ == static_cast<int>(focusDofIdx))
253 {
254 const Evaluation& rhoIn = intQuantsIn.fluidState().density(phaseIdx);
255 pStatIn = -rhoIn * (gIn * distVecIn);
256 }
257 else {
258 const Scalar rhoIn = Toolbox::value(intQuantsIn.fluidState().density(phaseIdx));
259 pStatIn = -rhoIn * (gIn * distVecIn);
260 }
261
262 // the quantities on the exterior side of the face do not influence the
263 // result for the TPFA scheme, so they can be treated as scalar values.
264 Evaluation pStatEx;
265
266 if (std::is_same<Scalar, Evaluation>::value ||
267 exteriorDofIdx_ == static_cast<int>(focusDofIdx))
268 {
269 const Evaluation& rhoEx = intQuantsEx.fluidState().density(phaseIdx);
270 pStatEx = -rhoEx * (gEx * distVecEx);
271 }
272 else {
273 const Scalar rhoEx = Toolbox::value(intQuantsEx.fluidState().density(phaseIdx));
274 pStatEx = -rhoEx * (gEx * distVecEx);
275 }
276
277 // compute the hydrostatic gradient between the two control volumes (this
278 // gradient exhibitis the same direction as the vector between the two
279 // control volume centers and the length (pStaticExterior -
280 // pStaticInterior)/distanceInteriorToExterior
281 Dune::FieldVector<Evaluation, dimWorld> f(distVecTotal);
282 f *= (pStatEx - pStatIn) / absDistTotalSquared;
283
284 // calculate the final potential gradient
285 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
286 potentialGrad_[phaseIdx][dimIdx] += f[dimIdx];
287 }
288
289 for (unsigned dimIdx = 0; dimIdx < potentialGrad_[phaseIdx].size(); ++dimIdx) {
290 if (!isfinite(potentialGrad_[phaseIdx][dimIdx])) {
291 throw NumericalProblem("Non-finite potential gradient for phase '"
292 + std::string(FluidSystem::phaseName(phaseIdx)) + "'");
293 }
294 }
295 }
296 }
297
298 Valgrind::SetUndefined(K_);
299 elemCtx.problem().intersectionIntrinsicPermeability(K_, elemCtx, faceIdx, timeIdx);
300 Valgrind::CheckDefined(K_);
301
302 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
303 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
304 Valgrind::SetUndefined(potentialGrad_[phaseIdx]);
305 continue;
306 }
307
308 // determine the upstream and downstream DOFs
309 Evaluation tmp = 0.0;
310 for (unsigned dimIdx = 0; dimIdx < faceNormal.size(); ++dimIdx) {
311 tmp += potentialGrad_[phaseIdx][dimIdx]*faceNormal[dimIdx];
312 }
313
314 if (tmp > 0) {
317 }
318 else {
321 }
322
323 // we only carry the derivatives along if the upstream DOF is the one which
324 // we currently focus on
325 const auto& up = elemCtx.intensiveQuantities(upstreamDofIdx_[phaseIdx], timeIdx);
326 if (upstreamDofIdx_[phaseIdx] == static_cast<int>(focusDofIdx)) {
327 mobility_[phaseIdx] = up.mobility(phaseIdx);
328 }
329 else {
330 mobility_[phaseIdx] = Toolbox::value(up.mobility(phaseIdx));
331 }
332 }
333 }
334
341 template <class FluidState>
342 void calculateBoundaryGradients_(const ElementContext& elemCtx,
343 unsigned boundaryFaceIdx,
344 unsigned timeIdx,
345 const FluidState& fluidState)
346 {
347 const auto& gradCalc = elemCtx.gradientCalculator();
348 BoundaryPressureCallback<TypeTag, FluidState> pressureCallback(elemCtx, fluidState);
349
350 // calculate the pressure gradient
351 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
352 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
353 Valgrind::SetUndefined(potentialGrad_[phaseIdx]);
354 continue;
355 }
356
357 pressureCallback.setPhaseIndex(phaseIdx);
358 gradCalc.calculateBoundaryGradient(potentialGrad_[phaseIdx],
359 elemCtx,
360 boundaryFaceIdx,
361 pressureCallback);
362 Valgrind::CheckDefined(potentialGrad_[phaseIdx]);
363 }
364
365 const auto& scvf = elemCtx.stencil(timeIdx).boundaryFace(boundaryFaceIdx);
366 const auto i = scvf.interiorIndex();
367 interiorDofIdx_ = static_cast<short>(i);
368 exteriorDofIdx_ = -1;
369 int focusDofIdx = elemCtx.focusDofIndex();
370
371 // calculate the intrinsic permeability
372 const auto& intQuantsIn = elemCtx.intensiveQuantities(i, timeIdx);
373 K_ = intQuantsIn.intrinsicPermeability();
374
375 // correct the pressure gradients by the gravitational acceleration
376 if (Parameters::Get<Parameters::EnableGravity>()) {
377 // estimate the gravitational acceleration at a given SCV face
378 // using the arithmetic mean
379 const auto& gIn = elemCtx.problem().gravity(elemCtx, i, timeIdx);
380 const auto& posIn = elemCtx.pos(i, timeIdx);
381 const auto& posFace = scvf.integrationPos();
382
383 // the distance between the face center and the center of the control volume
384 DimVector distVecIn(posIn);
385 distVecIn -= posFace;
386 const Scalar absDistSquared = distVecIn.two_norm2();
387 const Scalar gTimesDist = gIn * distVecIn;
388
389 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
390 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
391 continue;
392 }
393
394 // calculate the hydrostatic pressure at the integration point of the face
395 const Evaluation rhoIn = intQuantsIn.fluidState().density(phaseIdx);
396 const Evaluation pStatIn = -gTimesDist * rhoIn;
397
398 Valgrind::CheckDefined(pStatIn);
399 // compute the hydrostatic gradient between the control volume and face integration
400 // point. This gradient exhibitis the same direction as the vector between the
401 // control volume center and face integration point (-distVecIn / absDist) and the
402 // length of the vector is -pStaticIn / absDist. Note that the two negatives become
403 // + and the 1 / (absDist * absDist) -> absDistSquared.
404 EvalDimVector f(distVecIn);
405 f *= pStatIn / absDistSquared;
406
407 // calculate the final potential gradient
408 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
409 potentialGrad_[phaseIdx][dimIdx] += f[dimIdx];
410 }
411
412 Valgrind::CheckDefined(potentialGrad_[phaseIdx]);
413 for (unsigned dimIdx = 0; dimIdx < potentialGrad_[phaseIdx].size(); ++dimIdx) {
414 if (!isfinite(potentialGrad_[phaseIdx][dimIdx])) {
415 throw NumericalProblem("Non-finite potential gradient for phase '"
416 + std::string(FluidSystem::phaseName(phaseIdx)) + "'");
417 }
418 }
419 }
420 }
421
422 // determine the upstream and downstream DOFs
423 const auto& faceNormal = scvf.normal();
424
425 const auto& matParams = elemCtx.problem().materialLawParams(elemCtx, i, timeIdx);
426
427 std::array<Scalar, numPhases> kr;
428 MaterialLaw::relativePermeabilities(kr, matParams, fluidState);
429 Valgrind::CheckDefined(kr);
430
431 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
432 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
433 continue;
434 }
435
436 Evaluation tmp = 0.0;
437 for (unsigned dimIdx = 0; dimIdx < faceNormal.size(); ++dimIdx) {
438 tmp += potentialGrad_[phaseIdx][dimIdx] * faceNormal[dimIdx];
439 }
440
441 if (tmp > 0) {
444 }
445 else {
448 }
449
450 // take the phase mobility from the DOF in upstream direction
451 if (upstreamDofIdx_[phaseIdx] < 0) {
452 if (interiorDofIdx_ == focusDofIdx) {
453 mobility_[phaseIdx] = kr[phaseIdx] / fluidState.viscosity(phaseIdx);
454 }
455 else {
456 mobility_[phaseIdx] = Toolbox::value(kr[phaseIdx]) /
457 Toolbox::value(fluidState.viscosity(phaseIdx));
458 }
459 }
460 else if (upstreamDofIdx_[phaseIdx] != focusDofIdx) {
461 mobility_[phaseIdx] = Toolbox::value(intQuantsIn.mobility(phaseIdx));
462 }
463 else {
464 mobility_[phaseIdx] = intQuantsIn.mobility(phaseIdx);
465 }
466 Valgrind::CheckDefined(mobility_[phaseIdx]);
467 }
468 }
469
476 void calculateFluxes_(const ElementContext& elemCtx, unsigned scvfIdx, unsigned timeIdx)
477 {
478 const auto& scvf = elemCtx.stencil(timeIdx).interiorFace(scvfIdx);
479 const DimVector& normal = scvf.normal();
480 Valgrind::CheckDefined(normal);
481
482 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
483 filterVelocity_[phaseIdx] = 0.0;
484 volumeFlux_[phaseIdx] = 0.0;
485 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
486 continue;
487 }
488
489 asImp_().calculateFilterVelocity_(phaseIdx);
490 Valgrind::CheckDefined(filterVelocity_[phaseIdx]);
491
492 volumeFlux_[phaseIdx] = 0.0;
493 for (unsigned i = 0; i < normal.size(); ++i) {
494 volumeFlux_[phaseIdx] += filterVelocity_[phaseIdx][i] * normal[i];
495 }
496 }
497 }
498
505 void calculateBoundaryFluxes_(const ElementContext& elemCtx,
506 unsigned boundaryFaceIdx,
507 unsigned timeIdx)
508 {
509 const auto& scvf = elemCtx.stencil(timeIdx).boundaryFace(boundaryFaceIdx);
510 const DimVector& normal = scvf.normal();
511 Valgrind::CheckDefined(normal);
512
513 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
514 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
515 filterVelocity_[phaseIdx] = 0.0;
516 volumeFlux_[phaseIdx] = 0.0;
517 continue;
518 }
519
520 asImp_().calculateFilterVelocity_(phaseIdx);
521 Valgrind::CheckDefined(filterVelocity_[phaseIdx]);
522 volumeFlux_[phaseIdx] = 0.0;
523 for (unsigned i = 0; i < normal.size(); ++i) {
524 volumeFlux_[phaseIdx] += filterVelocity_[phaseIdx][i] * normal[i];
525 }
526 }
527 }
528
529 void calculateFilterVelocity_(unsigned phaseIdx)
530 {
531#ifndef NDEBUG
532 assert(isfinite(mobility_[phaseIdx]));
533 for (unsigned i = 0; i < K_.M(); ++i) {
534 for (unsigned j = 0; j < K_.N(); ++j) {
535 assert(std::isfinite(K_[i][j]));
536 }
537 }
538#endif
539
540 K_.mv(potentialGrad_[phaseIdx], filterVelocity_[phaseIdx]);
541 filterVelocity_[phaseIdx] *= - mobility_[phaseIdx];
542
543#ifndef NDEBUG
544 for (unsigned i = 0; i < filterVelocity_[phaseIdx].size(); ++i) {
545 assert(isfinite(filterVelocity_[phaseIdx][i]));
546 }
547#endif
548 }
549
550private:
551 Implementation& asImp_()
552 { return *static_cast<Implementation*>(this); }
553
554 const Implementation& asImp_() const
555 { return *static_cast<const Implementation*>(this); }
556
557protected:
558 // intrinsic permeability tensor and its square root
559 DimMatrix K_;
560
561 // mobilities of all fluid phases [1 / (Pa s)]
562 std::array<Evaluation, numPhases> mobility_;
563
564 // filter velocities of all phases [m/s]
565 std::array<EvalDimVector, numPhases> filterVelocity_;
566
567 // the volumetric flux of all fluid phases over the control
568 // volume's face [m^3/s / m^2]
569 std::array<Evaluation, numPhases> volumeFlux_;
570
571 // pressure potential gradients of all phases [Pa / m]
572 std::array<EvalDimVector, numPhases> potentialGrad_;
573
574 // upstream, downstream, interior and exterior DOFs
575 std::array<short, numPhases> upstreamDofIdx_;
576 std::array<short, numPhases> downstreamDofIdx_;
579};
580
581} // namespace Opm
582
583#endif
Callback class for a phase pressure.
Definition: quantitycallbacks.hh:134
void setPhaseIndex(unsigned phaseIdx)
Set the index of the fluid phase for which the pressure should be returned.
Definition: quantitycallbacks.hh:163
Provides the defaults for the parameters required by the Darcy velocity approach.
Definition: darcyfluxmodule.hh:86
Provides the Darcy flux module.
Definition: darcyfluxmodule.hh:122
std::array< EvalDimVector, numPhases > filterVelocity_
Definition: darcyfluxmodule.hh:565
std::array< EvalDimVector, numPhases > potentialGrad_
Definition: darcyfluxmodule.hh:572
void calculateBoundaryGradients_(const ElementContext &elemCtx, unsigned boundaryFaceIdx, unsigned timeIdx, const FluidState &fluidState)
Calculate the gradients at the grid boundary which are required to determine the volumetric fluxes.
Definition: darcyfluxmodule.hh:342
short interiorDofIdx_
Definition: darcyfluxmodule.hh:577
const EvalDimVector & filterVelocity(unsigned phaseIdx) const
Return the filter velocity of a fluid phase at the face's integration point [m/s].
Definition: darcyfluxmodule.hh:162
std::array< Evaluation, numPhases > volumeFlux_
Definition: darcyfluxmodule.hh:569
void calculateBoundaryFluxes_(const ElementContext &elemCtx, unsigned boundaryFaceIdx, unsigned timeIdx)
Calculate the volumetric fluxes at a boundary face of all fluid phases.
Definition: darcyfluxmodule.hh:505
std::array< Evaluation, numPhases > mobility_
Definition: darcyfluxmodule.hh:562
short upstreamIndex_(unsigned phaseIdx) const
Definition: darcyfluxmodule.hh:178
std::array< short, numPhases > upstreamDofIdx_
Definition: darcyfluxmodule.hh:575
const EvalDimVector & potentialGrad(unsigned phaseIdx) const
Return the pressure potential gradient of a fluid phase at the face's integration point [Pa/m].
Definition: darcyfluxmodule.hh:153
void calculateFluxes_(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Calculate the volumetric fluxes of all phases.
Definition: darcyfluxmodule.hh:476
const Evaluation & volumeFlux(unsigned phaseIdx) const
Return the volume flux of a fluid phase at the face's integration point .
Definition: darcyfluxmodule.hh:174
short downstreamIndex_(unsigned phaseIdx) const
Definition: darcyfluxmodule.hh:181
std::array< short, numPhases > downstreamDofIdx_
Definition: darcyfluxmodule.hh:576
void calculateFilterVelocity_(unsigned phaseIdx)
Definition: darcyfluxmodule.hh:529
void calculateGradients_(const ElementContext &elemCtx, unsigned faceIdx, unsigned timeIdx)
Calculate the gradients which are required to determine the volumetric fluxes.
Definition: darcyfluxmodule.hh:189
short exteriorDofIdx_
Definition: darcyfluxmodule.hh:578
const DimMatrix & intrinsicPermability() const
Returns the intrinsic permeability tensor for a given sub-control volume face.
Definition: darcyfluxmodule.hh:144
DimMatrix K_
Definition: darcyfluxmodule.hh:559
Provides the intensive quantities for the Darcy flux module.
Definition: darcyfluxmodule.hh:94
void update_(const ElementContext &, unsigned, unsigned)
Definition: darcyfluxmodule.hh:98
Callback class for a phase pressure.
Definition: quantitycallbacks.hh:85
void setPhaseIndex(unsigned phaseIdx)
Set the index of the fluid phase for which the pressure should be returned.
Definition: quantitycallbacks.hh:109
Defines the common parameters for the porous medium multi-phase models.
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...
Specifies a flux module which uses the Darcy relation.
Definition: darcyfluxmodule.hh:67
static void registerParameters()
Register all run-time parameters for the flux module.
Definition: darcyfluxmodule.hh:75