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