forchheimerfluxmodule.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_FORCHHEIMER_FLUX_MODULE_HH
31#define EWOMS_FORCHHEIMER_FLUX_MODULE_HH
32
33#include <opm/common/Exceptions.hpp>
34
35#include <opm/material/common/Valgrind.hpp>
36
39
40#include <dune/common/fvector.hh>
41#include <dune/common/fmatrix.hh>
42
43#include <algorithm>
44#include <array>
45#include <cassert>
46#include <cmath>
47#include <stdexcept>
48#include <string>
49
50namespace Opm {
51template <class TypeTag>
52class ForchheimerIntensiveQuantities;
53
54template <class TypeTag>
55class ForchheimerExtensiveQuantities;
56
57template <class TypeTag>
58class ForchheimerBaseProblem;
59
64template <class TypeTag>
66{
70
74 static void registerParameters()
75 {}
76};
77
83template <class TypeTag>
85{
88
89public:
99 template <class Context>
100 Scalar ergunCoefficient(const Context&,
101 unsigned,
102 unsigned) const
103 {
104 throw std::logic_error("Not implemented: Problem::ergunCoefficient()");
105 }
106
118 template <class Context>
119 Evaluation mobilityPassabilityRatio(Context& context,
120 unsigned spaceIdx,
121 unsigned timeIdx,
122 unsigned phaseIdx) const
123 {
124 return 1.0 / context.intensiveQuantities(spaceIdx, timeIdx).fluidState().viscosity(phaseIdx);
125 }
126};
127
132template <class TypeTag>
134{
137
138 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
139
140public:
148 const Evaluation& ergunCoefficient() const
149 { return ergunCoefficient_; }
150
154 const Evaluation& mobilityPassabilityRatio(unsigned phaseIdx) const
155 { return mobilityPassabilityRatio_[phaseIdx]; }
156
157protected:
158 void update_(const ElementContext& elemCtx, unsigned dofIdx, unsigned timeIdx)
159 {
160 const auto& problem = elemCtx.problem();
161 ergunCoefficient_ = problem.ergunCoefficient(elemCtx, dofIdx, timeIdx);
162
163 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
164 mobilityPassabilityRatio_[phaseIdx] =
165 problem.mobilityPassabilityRatio(elemCtx,
166 dofIdx,
167 timeIdx,
168 phaseIdx);
169 }
170 }
171
172private:
173 Evaluation ergunCoefficient_;
174 std::array<Evaluation, numPhases> mobilityPassabilityRatio_;
175};
176
216template <class TypeTag>
218 : public DarcyExtensiveQuantities<TypeTag>
219{
226
227 enum { dimWorld = GridView::dimensionworld };
228 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
229
230 using Toolbox = MathToolbox<Evaluation>;
231
232 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
233 using DimEvalVector = Dune::FieldVector<Evaluation, dimWorld>;
234 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
235 using DimEvalMatrix = Dune::FieldMatrix<Evaluation, dimWorld, dimWorld>;
236
237public:
241 const Evaluation& ergunCoefficient() const
242 { return ergunCoefficient_; }
243
250 Evaluation& mobilityPassabilityRatio(unsigned phaseIdx) const
251 { return mobilityPassabilityRatio_[phaseIdx]; }
252
253protected:
254 void calculateGradients_(const ElementContext& elemCtx,
255 unsigned faceIdx,
256 unsigned timeIdx)
257 {
258 DarcyExtQuants::calculateGradients_(elemCtx, faceIdx, timeIdx);
259
260 const auto focusDofIdx = elemCtx.focusDofIndex();
261 const unsigned i = static_cast<unsigned>(this->interiorDofIdx_);
262 const unsigned j = static_cast<unsigned>(this->exteriorDofIdx_);
263 const auto& intQuantsIn = elemCtx.intensiveQuantities(i, timeIdx);
264 const auto& intQuantsEx = elemCtx.intensiveQuantities(j, timeIdx);
265
266 // calculate the square root of the intrinsic permeability
267 assert(isDiagonal_(this->K_));
268 sqrtK_ = 0.0;
269 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
270 sqrtK_[dimIdx] = std::sqrt(this->K_[dimIdx][dimIdx]);
271 }
272
273 // obtain the Ergun coefficient. Lacking better ideas, we use its the arithmetic mean.
274 if (focusDofIdx == i) {
276 (intQuantsIn.ergunCoefficient() +
277 getValue(intQuantsEx.ergunCoefficient())) / 2;
278 }
279 else if (focusDofIdx == j) {
281 (getValue(intQuantsIn.ergunCoefficient()) +
282 intQuantsEx.ergunCoefficient()) / 2;
283 }
284 else {
286 (getValue(intQuantsIn.ergunCoefficient()) +
287 getValue(intQuantsEx.ergunCoefficient())) / 2;
288 }
289
290 // obtain the mobility to passability ratio for each phase.
291 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
292 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
293 continue;
294 }
295
296 const unsigned upIdx = static_cast<unsigned>(this->upstreamIndex_(phaseIdx));
297 const auto& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
298
299 if (focusDofIdx == upIdx) {
300 density_[phaseIdx] = up.fluidState().density(phaseIdx);
301 mobilityPassabilityRatio_[phaseIdx] = up.mobilityPassabilityRatio(phaseIdx);
302 }
303 else {
304 density_[phaseIdx] = getValue(up.fluidState().density(phaseIdx));
305 mobilityPassabilityRatio_[phaseIdx] = getValue(up.mobilityPassabilityRatio(phaseIdx));
306 }
307 }
308 }
309
310 template <class FluidState>
311 void calculateBoundaryGradients_(const ElementContext& elemCtx,
312 unsigned boundaryFaceIdx,
313 unsigned timeIdx,
314 const FluidState& fluidState)
315 {
317 boundaryFaceIdx,
318 timeIdx,
319 fluidState);
320
321 const auto focusDofIdx = elemCtx.focusDofIndex();
322 const unsigned i = static_cast<unsigned>(this->interiorDofIdx_);
323 const auto& intQuantsIn = elemCtx.intensiveQuantities(i, timeIdx);
324
325 // obtain the Ergun coefficient. Because we are on the boundary here, we will
326 // take the Ergun coefficient of the interior
327 if (focusDofIdx == i) {
328 ergunCoefficient_ = intQuantsIn.ergunCoefficient();
329 }
330 else {
331 ergunCoefficient_ = getValue(intQuantsIn.ergunCoefficient());
332 }
333
334 // calculate the square root of the intrinsic permeability
335 assert(isDiagonal_(this->K_));
336 sqrtK_ = 0.0;
337 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
338 sqrtK_[dimIdx] = std::sqrt(this->K_[dimIdx][dimIdx]);
339 }
340
341 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
342 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
343 continue;
344 }
345
346 if (focusDofIdx == i) {
347 density_[phaseIdx] = intQuantsIn.fluidState().density(phaseIdx);
348 mobilityPassabilityRatio_[phaseIdx] = intQuantsIn.mobilityPassabilityRatio(phaseIdx);
349 }
350 else {
351 density_[phaseIdx] = getValue(intQuantsIn.fluidState().density(phaseIdx));
352 mobilityPassabilityRatio_[phaseIdx] =
353 getValue(intQuantsIn.mobilityPassabilityRatio(phaseIdx));
354 }
355 }
356 }
357
364 void calculateFluxes_(const ElementContext& elemCtx, unsigned scvfIdx, unsigned timeIdx)
365 {
366 const auto focusDofIdx = elemCtx.focusDofIndex();
367 const auto i = asImp_().interiorIndex();
368 const auto j = asImp_().exteriorIndex();
369 const auto& intQuantsI = elemCtx.intensiveQuantities(i, timeIdx);
370 const auto& intQuantsJ = elemCtx.intensiveQuantities(j, timeIdx);
371
372 const auto& scvf = elemCtx.stencil(timeIdx).interiorFace(scvfIdx);
373 const auto& normal = scvf.normal();
374 Valgrind::CheckDefined(normal);
375
376 // obtain the Ergun coefficient from the intensive quantity object. Until a
377 // better method comes along, we use arithmetic averaging.
378 if (focusDofIdx == i) {
379 ergunCoefficient_ = (intQuantsI.ergunCoefficient() +
380 getValue(intQuantsJ.ergunCoefficient())) / 2;
381 }
382 else if (focusDofIdx == j) {
383 ergunCoefficient_ = (getValue(intQuantsI.ergunCoefficient()) +
384 intQuantsJ.ergunCoefficient()) / 2;
385 }
386 else {
387 ergunCoefficient_ = (getValue(intQuantsI.ergunCoefficient()) +
388 getValue(intQuantsJ.ergunCoefficient())) / 2;
389 }
390
392 // calculate the weights of the upstream and the downstream control volumes
394 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
395 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
396 this->filterVelocity_[phaseIdx] = 0.0;
397 this->volumeFlux_[phaseIdx] = 0.0;
398 continue;
399 }
400
402
403 this->volumeFlux_[phaseIdx] = 0.0;
404 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
405 this->volumeFlux_[phaseIdx] +=
406 this->filterVelocity_[phaseIdx][dimIdx]*normal[dimIdx];
407 }
408 }
409 }
410
414 void calculateBoundaryFluxes_(const ElementContext& elemCtx,
415 unsigned bfIdx,
416 unsigned timeIdx)
417 {
418 const auto& boundaryFace = elemCtx.stencil(timeIdx).boundaryFace(bfIdx);
419 const auto& normal = boundaryFace.normal();
420
422 // calculate the weights of the upstream and the downstream degrees of freedom
424 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
425 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
426 this->filterVelocity_[phaseIdx] = 0.0;
427 this->volumeFlux_[phaseIdx] = 0.0;
428 continue;
429 }
430
432
433 this->volumeFlux_[phaseIdx] = 0.0;
434 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
435 this->volumeFlux_[phaseIdx] +=
436 this->filterVelocity_[phaseIdx][dimIdx]*normal[dimIdx];
437 }
438 }
439 }
440
441 void calculateForchheimerFlux_(unsigned phaseIdx)
442 {
443 // initial guess: filter velocity is zero
444 DimEvalVector& velocity = this->filterVelocity_[phaseIdx];
445 velocity = 0.0;
446
447 // the change of velocity between two consecutive Newton iterations
448 DimEvalVector deltaV(1e5);
449 // the function value that is to be minimized of the equation that is to be
450 // fulfilled
451 DimEvalVector residual;
452 // derivative of equation that is to be solved
453 DimEvalMatrix gradResid;
454
455 // search by means of the Newton method for a root of Forchheimer equation
456 unsigned newtonIter = 0;
457 while (deltaV.one_norm() > 1e-11) {
458 if (newtonIter >= 50) {
459 throw NumericalProblem("Could not determine Forchheimer velocity within "
460 + std::to_string(newtonIter) + " iterations");
461 }
462 ++newtonIter;
463
464 // calculate the residual and its Jacobian matrix
465 gradForchheimerResid_(residual, gradResid, phaseIdx);
466
467 // newton method
468 gradResid.solve(deltaV, residual);
469 velocity -= deltaV;
470 }
471 }
472
473 void forchheimerResid_(DimEvalVector& residual, unsigned phaseIdx) const
474 {
475 const DimEvalVector& velocity = this->filterVelocity_[phaseIdx];
476
477 // Obtaining the upstreamed quantities
478 const auto& mobility = this->mobility_[phaseIdx];
479 const auto& density = density_[phaseIdx];
480 const auto& mobPassabilityRatio = mobilityPassabilityRatio_[phaseIdx];
481
482 // optain the quantites for the integration point
483 const auto& pGrad = this->potentialGrad_[phaseIdx];
484
485 // residual = v_\alpha
486 residual = velocity;
487
488 // residual += mobility_\alpha K(\grad p_\alpha - \rho_\alpha g)
489 // -> this->K_.usmv(mobility, pGrad, residual);
490 assert(isDiagonal_(this->K_));
491 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++ dimIdx) {
492 residual[dimIdx] += mobility*pGrad[dimIdx]*this->K_[dimIdx][dimIdx];
493 }
494
495 // Forchheimer turbulence correction:
496 //
497 // residual +=
498 // \rho_\alpha
499 // * mobility_\alpha
500 // * C_E / \eta_{r,\alpha}
501 // * abs(v_\alpha) * sqrt(K)*v_\alpha
502 //
503 // -> sqrtK_.usmv(density*mobilityPassabilityRatio*ergunCoefficient_*velocity.two_norm(),
504 // velocity,
505 // residual);
506 Evaluation absVel = 0.0;
507 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
508 absVel += velocity[dimIdx]*velocity[dimIdx];
509 }
510 // the derivatives of the square root of 0 are undefined, so we must guard
511 // against this case
512 if (absVel <= 0.0) {
513 absVel = 0.0;
514 }
515 else {
516 absVel = Toolbox::sqrt(absVel);
517 }
518 const auto& alpha = density * mobPassabilityRatio * ergunCoefficient_ * absVel;
519 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
520 residual[dimIdx] += sqrtK_[dimIdx] * alpha * velocity[dimIdx];
521 }
522 Valgrind::CheckDefined(residual);
523 }
524
525 void gradForchheimerResid_(DimEvalVector& residual,
526 DimEvalMatrix& gradResid,
527 unsigned phaseIdx)
528 {
529 // TODO (?) use AD for this.
530 DimEvalVector& velocity = this->filterVelocity_[phaseIdx];
531 forchheimerResid_(residual, phaseIdx);
532
533 constexpr Scalar eps = 1e-11;
534 DimEvalVector tmp;
535 for (unsigned i = 0; i < dimWorld; ++i) {
536 const Scalar coordEps = std::max(eps, Toolbox::scalarValue(velocity[i]) * (1 + eps));
537 velocity[i] += coordEps;
538 forchheimerResid_(tmp, phaseIdx);
539 tmp -= residual;
540 tmp /= coordEps;
541 gradResid[i] = tmp;
542 velocity[i] -= coordEps;
543 }
544 }
545
553 bool isDiagonal_(const DimMatrix& K) const
554 {
555 for (unsigned i = 0; i < dimWorld; i++) {
556 for (unsigned j = 0; j < dimWorld; j++) {
557 if (i == j) {
558 continue;
559 }
560
561 if (std::abs(K[i][j]) > 1e-25) {
562 return false;
563 }
564 }
565 }
566 return true;
567 }
568
569private:
570 Implementation& asImp_()
571 { return *static_cast<Implementation *>(this); }
572
573 const Implementation& asImp_() const
574 { return *static_cast<const Implementation *>(this); }
575
576protected:
577 // intrinsic permeability tensor and its square root
578 DimVector sqrtK_;
579
580 // Ergun coefficient of all phases at the integration point
582
583 // Passability of all phases at the integration point
584 std::array<Evaluation, numPhases> mobilityPassabilityRatio_;
585
586 // Density of all phases at the integration point
587 std::array<Evaluation, numPhases> density_;
588};
589
590} // namespace Opm
591
592#endif
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
std::array< Evaluation, numPhases > volumeFlux_
Definition: darcyfluxmodule.hh:569
std::array< Evaluation, numPhases > mobility_
Definition: darcyfluxmodule.hh:562
short upstreamIndex_(unsigned phaseIdx) const
Definition: darcyfluxmodule.hh:178
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
DimMatrix K_
Definition: darcyfluxmodule.hh:559
Provides the defaults for the parameters required by the Forchheimer velocity approach.
Definition: forchheimerfluxmodule.hh:85
Evaluation mobilityPassabilityRatio(Context &context, unsigned spaceIdx, unsigned timeIdx, unsigned phaseIdx) const
Returns the ratio between the phase mobility and its passability for a given fluid phase .
Definition: forchheimerfluxmodule.hh:119
Scalar ergunCoefficient(const Context &, unsigned, unsigned) const
Returns the Ergun coefficient.
Definition: forchheimerfluxmodule.hh:100
Provides the Forchheimer flux module.
Definition: forchheimerfluxmodule.hh:219
void calculateGradients_(const ElementContext &elemCtx, unsigned faceIdx, unsigned timeIdx)
Definition: forchheimerfluxmodule.hh:254
void calculateForchheimerFlux_(unsigned phaseIdx)
Definition: forchheimerfluxmodule.hh:441
const Evaluation & ergunCoefficient() const
Return the Ergun coefficent at the face's integration point.
Definition: forchheimerfluxmodule.hh:241
void calculateBoundaryGradients_(const ElementContext &elemCtx, unsigned boundaryFaceIdx, unsigned timeIdx, const FluidState &fluidState)
Definition: forchheimerfluxmodule.hh:311
void forchheimerResid_(DimEvalVector &residual, unsigned phaseIdx) const
Definition: forchheimerfluxmodule.hh:473
std::array< Evaluation, numPhases > mobilityPassabilityRatio_
Definition: forchheimerfluxmodule.hh:584
DimVector sqrtK_
Definition: forchheimerfluxmodule.hh:578
void gradForchheimerResid_(DimEvalVector &residual, DimEvalMatrix &gradResid, unsigned phaseIdx)
Definition: forchheimerfluxmodule.hh:525
void calculateBoundaryFluxes_(const ElementContext &elemCtx, unsigned bfIdx, unsigned timeIdx)
Calculate the volumetric flux rates of all phases at the domain boundary.
Definition: forchheimerfluxmodule.hh:414
Evaluation ergunCoefficient_
Definition: forchheimerfluxmodule.hh:581
void calculateFluxes_(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Calculate the volumetric fluxes of all phases.
Definition: forchheimerfluxmodule.hh:364
std::array< Evaluation, numPhases > density_
Definition: forchheimerfluxmodule.hh:587
bool isDiagonal_(const DimMatrix &K) const
Check whether all off-diagonal entries of a tensor are zero.
Definition: forchheimerfluxmodule.hh:553
Evaluation & mobilityPassabilityRatio(unsigned phaseIdx) const
Return the ratio of the mobility divided by the passability at the face's integration point for a giv...
Definition: forchheimerfluxmodule.hh:250
Provides the intensive quantities for the Forchheimer module.
Definition: forchheimerfluxmodule.hh:134
const Evaluation & ergunCoefficient() const
Returns the Ergun coefficient.
Definition: forchheimerfluxmodule.hh:148
void update_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: forchheimerfluxmodule.hh:158
const Evaluation & mobilityPassabilityRatio(unsigned phaseIdx) const
Returns the passability of a phase.
Definition: forchheimerfluxmodule.hh:154
This file contains the necessary classes to calculate the volumetric fluxes out of a pressure potenti...
Declare the properties used by the infrastructure code of the finite volume discretizations.
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
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Specifies a flux module which uses the Forchheimer relation.
Definition: forchheimerfluxmodule.hh:66
static void registerParameters()
Register all run-time parameters for the flux module.
Definition: forchheimerfluxmodule.hh:74