30#ifndef EWOMS_FORCHHEIMER_FLUX_MODULE_HH
31#define EWOMS_FORCHHEIMER_FLUX_MODULE_HH
35#include <opm/common/Exceptions.hpp>
39#include <opm/material/common/Valgrind.hpp>
41#include <dune/common/fvector.hh>
42#include <dune/common/fmatrix.hh>
47template <
class TypeTag>
48class ForchheimerIntensiveQuantities;
50template <
class TypeTag>
51class ForchheimerExtensiveQuantities;
53template <
class TypeTag>
54class ForchheimerBaseProblem;
60template <
class TypeTag>
79template <
class TypeTag>
95 template <
class Context>
100 throw std::logic_error(
"Not implemented: Problem::ergunCoefficient()");
114 template <
class Context>
118 unsigned phaseIdx)
const
120 return 1.0 / context.intensiveQuantities(spaceIdx, timeIdx).fluidState().viscosity(phaseIdx);
128template <
class TypeTag>
135 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
146 {
return ergunCoefficient_; }
152 {
return mobilityPassabilityRatio_[phaseIdx]; }
155 void update_(
const ElementContext& elemCtx,
unsigned dofIdx,
unsigned timeIdx)
157 const auto& problem = elemCtx.problem();
158 ergunCoefficient_ = problem.ergunCoefficient(elemCtx, dofIdx, timeIdx);
160 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
161 mobilityPassabilityRatio_[phaseIdx] =
162 problem.mobilityPassabilityRatio(elemCtx,
169 Evaluation ergunCoefficient_;
170 Evaluation mobilityPassabilityRatio_[numPhases];
212template <
class TypeTag>
225 enum { dimWorld = GridView::dimensionworld };
226 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
228 using Toolbox = MathToolbox<Evaluation>;
230 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
231 using DimEvalVector = Dune::FieldVector<Evaluation, dimWorld>;
232 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
233 using DimEvalMatrix = Dune::FieldMatrix<Evaluation, dimWorld, dimWorld>;
258 auto focusDofIdx = elemCtx.focusDofIndex();
261 const auto& intQuantsIn = elemCtx.intensiveQuantities(i, timeIdx);
262 const auto& intQuantsEx = elemCtx.intensiveQuantities(j, timeIdx);
267 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx)
268 sqrtK_[dimIdx] = std::sqrt(this->
K_[dimIdx][dimIdx]);
271 if (focusDofIdx == i) {
273 (intQuantsIn.ergunCoefficient() +
274 getValue(intQuantsEx.ergunCoefficient()))/2;
276 else if (focusDofIdx == j)
278 (getValue(intQuantsIn.ergunCoefficient()) +
279 intQuantsEx.ergunCoefficient())/2;
282 (getValue(intQuantsIn.ergunCoefficient()) +
283 getValue(intQuantsEx.ergunCoefficient()))/2;
286 for (
unsigned phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
287 if (!elemCtx.model().phaseIsConsidered(phaseIdx))
290 unsigned upIdx =
static_cast<unsigned>(this->
upstreamIndex_(phaseIdx));
291 const auto& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
293 if (focusDofIdx == upIdx) {
295 up.fluidState().density(phaseIdx);
297 up.mobilityPassabilityRatio(phaseIdx);
301 getValue(up.fluidState().density(phaseIdx));
303 getValue(up.mobilityPassabilityRatio(phaseIdx));
308 template <
class Flu
idState>
310 unsigned boundaryFaceIdx,
312 const FluidState& fluidState)
319 auto focusDofIdx = elemCtx.focusDofIndex();
321 const auto& intQuantsIn = elemCtx.intensiveQuantities(i, timeIdx);
325 if (focusDofIdx == i)
333 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx)
334 sqrtK_[dimIdx] = std::sqrt(this->
K_[dimIdx][dimIdx]);
336 for (
unsigned phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
337 if (!elemCtx.model().phaseIsConsidered(phaseIdx))
340 if (focusDofIdx == i) {
341 density_[phaseIdx] = intQuantsIn.fluidState().density(phaseIdx);
346 getValue(intQuantsIn.fluidState().density(phaseIdx));
348 getValue(intQuantsIn.mobilityPassabilityRatio(phaseIdx));
361 auto focusDofIdx = elemCtx.focusDofIndex();
362 auto i = asImp_().interiorIndex();
363 auto j = asImp_().exteriorIndex();
364 const auto& intQuantsI = elemCtx.intensiveQuantities(i, timeIdx);
365 const auto& intQuantsJ = elemCtx.intensiveQuantities(j, timeIdx);
367 const auto& scvf = elemCtx.stencil(timeIdx).interiorFace(scvfIdx);
368 const auto& normal = scvf.normal();
369 Valgrind::CheckDefined(normal);
373 if (focusDofIdx == i)
375 (intQuantsI.ergunCoefficient() +
376 getValue(intQuantsJ.ergunCoefficient())) / 2;
377 else if (focusDofIdx == j)
379 (getValue(intQuantsI.ergunCoefficient()) +
380 intQuantsJ.ergunCoefficient()) / 2;
383 (getValue(intQuantsI.ergunCoefficient()) +
384 getValue(intQuantsJ.ergunCoefficient())) / 2;
389 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; phaseIdx++) {
390 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
399 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++ dimIdx)
412 const auto& boundaryFace = elemCtx.stencil(timeIdx).boundaryFace(bfIdx);
413 const auto& normal = boundaryFace.normal();
418 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; phaseIdx++) {
419 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
428 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx)
441 DimEvalVector deltaV(1e5);
444 DimEvalVector residual;
446 DimEvalMatrix gradResid;
449 unsigned newtonIter = 0;
450 while (deltaV.one_norm() > 1e-11) {
451 if (newtonIter >= 50)
452 throw NumericalProblem(
"Could not determine Forchheimer velocity within "
453 + std::to_string(newtonIter)+
" iterations");
460 gradResid.solve(deltaV, residual);
470 const auto& mobility = this->
mobility_[phaseIdx];
471 const auto& density =
density_[phaseIdx];
483 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++ dimIdx)
484 residual[dimIdx] += mobility*pGrad[dimIdx]*this->
K_[dimIdx][dimIdx];
497 Evaluation absVel = 0.0;
498 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx)
499 absVel += velocity[dimIdx]*velocity[dimIdx];
505 absVel = Toolbox::sqrt(absVel);
507 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx)
508 residual[dimIdx] +=
sqrtK_[dimIdx]*alpha*velocity[dimIdx];
509 Valgrind::CheckDefined(residual);
513 DimEvalMatrix& gradResid,
522 for (
unsigned i = 0; i < dimWorld; ++i) {
523 Scalar coordEps = std::max(eps, Toolbox::scalarValue(velocity[i]) * (1 + eps));
524 velocity[i] += coordEps;
529 velocity[i] -= coordEps;
542 for (
unsigned i = 0; i < dimWorld; i++) {
543 for (
unsigned j = 0; j < dimWorld; j++) {
547 if (std::abs(K[i][j]) > 1e-25)
555 Implementation& asImp_()
556 {
return *
static_cast<Implementation *
>(
this); }
558 const Implementation& asImp_()
const
559 {
return *
static_cast<const Implementation *
>(
this); }
Provides the Darcy flux module.
Definition: darcyfluxmodule.hh:115
Evaluation mobility_[numPhases]
Definition: darcyfluxmodule.hh:539
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:331
short interiorDofIdx_
Definition: darcyfluxmodule.hh:554
short upstreamIndex_(unsigned phaseIdx) const
Definition: darcyfluxmodule.hh:172
EvalDimVector potentialGrad_[numPhases]
Definition: darcyfluxmodule.hh:549
EvalDimVector filterVelocity_[numPhases]
Definition: darcyfluxmodule.hh:542
Evaluation volumeFlux_[numPhases]
Definition: darcyfluxmodule.hh:546
void calculateGradients_(const ElementContext &elemCtx, unsigned faceIdx, unsigned timeIdx)
Calculate the gradients which are required to determine the volumetric fluxes.
Definition: darcyfluxmodule.hh:183
short exteriorDofIdx_
Definition: darcyfluxmodule.hh:555
DimMatrix K_
Definition: darcyfluxmodule.hh:536
Provides the defaults for the parameters required by the Forchheimer velocity approach.
Definition: forchheimerfluxmodule.hh:81
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:115
Scalar ergunCoefficient(const Context &, unsigned, unsigned) const
Returns the Ergun coefficient.
Definition: forchheimerfluxmodule.hh:96
Provides the Forchheimer flux module.
Definition: forchheimerfluxmodule.hh:215
void calculateGradients_(const ElementContext &elemCtx, unsigned faceIdx, unsigned timeIdx)
Definition: forchheimerfluxmodule.hh:252
void calculateForchheimerFlux_(unsigned phaseIdx)
Definition: forchheimerfluxmodule.hh:434
const Evaluation & ergunCoefficient() const
Return the Ergun coefficent at the face's integration point.
Definition: forchheimerfluxmodule.hh:239
void calculateBoundaryGradients_(const ElementContext &elemCtx, unsigned boundaryFaceIdx, unsigned timeIdx, const FluidState &fluidState)
Definition: forchheimerfluxmodule.hh:309
void forchheimerResid_(DimEvalVector &residual, unsigned phaseIdx) const
Definition: forchheimerfluxmodule.hh:465
Evaluation density_[numPhases]
Definition: forchheimerfluxmodule.hh:572
Evaluation mobilityPassabilityRatio_[numPhases]
Definition: forchheimerfluxmodule.hh:569
DimVector sqrtK_
Definition: forchheimerfluxmodule.hh:563
void gradForchheimerResid_(DimEvalVector &residual, DimEvalMatrix &gradResid, unsigned phaseIdx)
Definition: forchheimerfluxmodule.hh:512
void calculateBoundaryFluxes_(const ElementContext &elemCtx, unsigned bfIdx, unsigned timeIdx)
Calculate the volumetric flux rates of all phases at the domain boundary.
Definition: forchheimerfluxmodule.hh:408
Evaluation ergunCoefficient_
Definition: forchheimerfluxmodule.hh:566
void calculateFluxes_(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Calculate the volumetric fluxes of all phases.
Definition: forchheimerfluxmodule.hh:359
bool isDiagonal_(const DimMatrix &K) const
Check whether all off-diagonal entries of a tensor are zero.
Definition: forchheimerfluxmodule.hh:540
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:248
Provides the intensive quantities for the Forchheimer module.
Definition: forchheimerfluxmodule.hh:130
const Evaluation & ergunCoefficient() const
Returns the Ergun coefficient.
Definition: forchheimerfluxmodule.hh:145
void update_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: forchheimerfluxmodule.hh:155
const Evaluation & mobilityPassabilityRatio(unsigned phaseIdx) const
Returns the passability of a phase.
Definition: forchheimerfluxmodule.hh:151
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: blackoilboundaryratevector.hh:37
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:235
Specifies a flux module which uses the Forchheimer relation.
Definition: forchheimerfluxmodule.hh:62
static void registerParameters()
Register all run-time parameters for the flux module.
Definition: forchheimerfluxmodule.hh:70