30#ifndef EWOMS_FORCHHEIMER_FLUX_MODULE_HH
31#define EWOMS_FORCHHEIMER_FLUX_MODULE_HH
33#include <opm/common/Exceptions.hpp>
35#include <opm/material/common/Valgrind.hpp>
40#include <dune/common/fvector.hh>
41#include <dune/common/fmatrix.hh>
51template <
class TypeTag>
52class ForchheimerIntensiveQuantities;
54template <
class TypeTag>
55class ForchheimerExtensiveQuantities;
57template <
class TypeTag>
58class ForchheimerBaseProblem;
64template <
class TypeTag>
83template <
class TypeTag>
99 template <
class Context>
104 throw std::logic_error(
"Not implemented: Problem::ergunCoefficient()");
118 template <
class Context>
122 unsigned phaseIdx)
const
124 return 1.0 / context.intensiveQuantities(spaceIdx, timeIdx).fluidState().viscosity(phaseIdx);
132template <
class TypeTag>
138 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
149 {
return ergunCoefficient_; }
155 {
return mobilityPassabilityRatio_[phaseIdx]; }
158 void update_(
const ElementContext& elemCtx,
unsigned dofIdx,
unsigned timeIdx)
160 const auto& problem = elemCtx.problem();
161 ergunCoefficient_ = problem.ergunCoefficient(elemCtx, dofIdx, timeIdx);
163 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
164 mobilityPassabilityRatio_[phaseIdx] =
165 problem.mobilityPassabilityRatio(elemCtx,
173 Evaluation ergunCoefficient_;
174 std::array<Evaluation, numPhases> mobilityPassabilityRatio_;
216template <
class TypeTag>
227 enum { dimWorld = GridView::dimensionworld };
228 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
230 using Toolbox = MathToolbox<Evaluation>;
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>;
260 const auto focusDofIdx = elemCtx.focusDofIndex();
263 const auto& intQuantsIn = elemCtx.intensiveQuantities(i, timeIdx);
264 const auto& intQuantsEx = elemCtx.intensiveQuantities(j, timeIdx);
269 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
270 sqrtK_[dimIdx] = std::sqrt(this->
K_[dimIdx][dimIdx]);
274 if (focusDofIdx == i) {
276 (intQuantsIn.ergunCoefficient() +
277 getValue(intQuantsEx.ergunCoefficient())) / 2;
279 else if (focusDofIdx == j) {
281 (getValue(intQuantsIn.ergunCoefficient()) +
282 intQuantsEx.ergunCoefficient()) / 2;
286 (getValue(intQuantsIn.ergunCoefficient()) +
287 getValue(intQuantsEx.ergunCoefficient())) / 2;
291 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
292 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
296 const unsigned upIdx =
static_cast<unsigned>(this->
upstreamIndex_(phaseIdx));
297 const auto& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
299 if (focusDofIdx == upIdx) {
300 density_[phaseIdx] = up.fluidState().density(phaseIdx);
304 density_[phaseIdx] = getValue(up.fluidState().density(phaseIdx));
310 template <
class Flu
idState>
312 unsigned boundaryFaceIdx,
314 const FluidState& fluidState)
321 const auto focusDofIdx = elemCtx.focusDofIndex();
323 const auto& intQuantsIn = elemCtx.intensiveQuantities(i, timeIdx);
327 if (focusDofIdx == i) {
337 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
338 sqrtK_[dimIdx] = std::sqrt(this->
K_[dimIdx][dimIdx]);
341 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
342 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
346 if (focusDofIdx == i) {
347 density_[phaseIdx] = intQuantsIn.fluidState().density(phaseIdx);
351 density_[phaseIdx] = getValue(intQuantsIn.fluidState().density(phaseIdx));
353 getValue(intQuantsIn.mobilityPassabilityRatio(phaseIdx));
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);
372 const auto& scvf = elemCtx.stencil(timeIdx).interiorFace(scvfIdx);
373 const auto& normal = scvf.normal();
374 Valgrind::CheckDefined(normal);
378 if (focusDofIdx == i) {
380 getValue(intQuantsJ.ergunCoefficient())) / 2;
382 else if (focusDofIdx == j) {
384 intQuantsJ.ergunCoefficient()) / 2;
388 getValue(intQuantsJ.ergunCoefficient())) / 2;
394 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
395 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
404 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
418 const auto& boundaryFace = elemCtx.stencil(timeIdx).boundaryFace(bfIdx);
419 const auto& normal = boundaryFace.normal();
424 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
425 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
434 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
448 DimEvalVector deltaV(1e5);
451 DimEvalVector residual;
453 DimEvalMatrix gradResid;
456 unsigned newtonIter = 0;
457 while (deltaV.one_norm() > 1e-11) {
458 if (newtonIter >= 50) {
459 throw NumericalProblem(
"Could not determine Forchheimer velocity within "
468 gradResid.solve(deltaV, residual);
478 const auto& mobility = this->
mobility_[phaseIdx];
479 const auto& density =
density_[phaseIdx];
491 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++ dimIdx) {
492 residual[dimIdx] += mobility*pGrad[dimIdx]*this->
K_[dimIdx][dimIdx];
506 Evaluation absVel = 0.0;
507 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
508 absVel += velocity[dimIdx]*velocity[dimIdx];
516 absVel = Toolbox::sqrt(absVel);
519 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
520 residual[dimIdx] +=
sqrtK_[dimIdx] * alpha * velocity[dimIdx];
522 Valgrind::CheckDefined(residual);
526 DimEvalMatrix& gradResid,
533 constexpr Scalar eps = 1e-11;
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;
542 velocity[i] -= coordEps;
555 for (
unsigned i = 0; i < dimWorld; i++) {
556 for (
unsigned j = 0; j < dimWorld; j++) {
561 if (std::abs(K[i][j]) > 1e-25) {
570 Implementation& asImp_()
571 {
return *
static_cast<Implementation *
>(
this); }
573 const Implementation& asImp_()
const
574 {
return *
static_cast<const Implementation *
>(
this); }
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