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>
139 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
150 {
return ergunCoefficient_; }
156 {
return mobilityPassabilityRatio_[phaseIdx]; }
159 void update_(
const ElementContext& elemCtx,
unsigned dofIdx,
unsigned timeIdx)
161 const auto& problem = elemCtx.problem();
162 ergunCoefficient_ = problem.ergunCoefficient(elemCtx, dofIdx, timeIdx);
164 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
165 mobilityPassabilityRatio_[phaseIdx] =
166 problem.mobilityPassabilityRatio(elemCtx,
174 Evaluation ergunCoefficient_;
175 std::array<Evaluation, numPhases> mobilityPassabilityRatio_;
217template <
class TypeTag>
230 enum { dimWorld = GridView::dimensionworld };
231 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
233 using Toolbox = MathToolbox<Evaluation>;
235 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
236 using DimEvalVector = Dune::FieldVector<Evaluation, dimWorld>;
237 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
238 using DimEvalMatrix = Dune::FieldMatrix<Evaluation, dimWorld, dimWorld>;
263 const auto focusDofIdx = elemCtx.focusDofIndex();
266 const auto& intQuantsIn = elemCtx.intensiveQuantities(i, timeIdx);
267 const auto& intQuantsEx = elemCtx.intensiveQuantities(j, timeIdx);
272 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
273 sqrtK_[dimIdx] = std::sqrt(this->
K_[dimIdx][dimIdx]);
277 if (focusDofIdx == i) {
279 (intQuantsIn.ergunCoefficient() +
280 getValue(intQuantsEx.ergunCoefficient())) / 2;
282 else if (focusDofIdx == j) {
284 (getValue(intQuantsIn.ergunCoefficient()) +
285 intQuantsEx.ergunCoefficient()) / 2;
289 (getValue(intQuantsIn.ergunCoefficient()) +
290 getValue(intQuantsEx.ergunCoefficient())) / 2;
294 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
295 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
299 const unsigned upIdx =
static_cast<unsigned>(this->
upstreamIndex_(phaseIdx));
300 const auto& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
302 if (focusDofIdx == upIdx) {
303 density_[phaseIdx] = up.fluidState().density(phaseIdx);
307 density_[phaseIdx] = getValue(up.fluidState().density(phaseIdx));
313 template <
class Flu
idState>
315 unsigned boundaryFaceIdx,
317 const FluidState& fluidState)
324 const auto focusDofIdx = elemCtx.focusDofIndex();
326 const auto& intQuantsIn = elemCtx.intensiveQuantities(i, timeIdx);
330 if (focusDofIdx == i) {
340 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
341 sqrtK_[dimIdx] = std::sqrt(this->
K_[dimIdx][dimIdx]);
344 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
345 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
349 if (focusDofIdx == i) {
350 density_[phaseIdx] = intQuantsIn.fluidState().density(phaseIdx);
354 density_[phaseIdx] = getValue(intQuantsIn.fluidState().density(phaseIdx));
356 getValue(intQuantsIn.mobilityPassabilityRatio(phaseIdx));
369 const auto focusDofIdx = elemCtx.focusDofIndex();
370 const auto i = asImp_().interiorIndex();
371 const auto j = asImp_().exteriorIndex();
372 const auto& intQuantsI = elemCtx.intensiveQuantities(i, timeIdx);
373 const auto& intQuantsJ = elemCtx.intensiveQuantities(j, timeIdx);
375 const auto& scvf = elemCtx.stencil(timeIdx).interiorFace(scvfIdx);
376 const auto& normal = scvf.normal();
377 Valgrind::CheckDefined(normal);
381 if (focusDofIdx == i) {
383 getValue(intQuantsJ.ergunCoefficient())) / 2;
385 else if (focusDofIdx == j) {
387 intQuantsJ.ergunCoefficient()) / 2;
391 getValue(intQuantsJ.ergunCoefficient())) / 2;
397 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
398 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
407 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
421 const auto& boundaryFace = elemCtx.stencil(timeIdx).boundaryFace(bfIdx);
422 const auto& normal = boundaryFace.normal();
427 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
428 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
437 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
451 DimEvalVector deltaV(1e5);
454 DimEvalVector residual;
456 DimEvalMatrix gradResid;
459 unsigned newtonIter = 0;
460 while (deltaV.one_norm() > 1e-11) {
461 if (newtonIter >= 50) {
462 throw NumericalProblem(
"Could not determine Forchheimer velocity within "
471 gradResid.solve(deltaV, residual);
481 const auto& mobility = this->
mobility_[phaseIdx];
482 const auto& density =
density_[phaseIdx];
494 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++ dimIdx) {
495 residual[dimIdx] += mobility*pGrad[dimIdx]*this->
K_[dimIdx][dimIdx];
509 Evaluation absVel = 0.0;
510 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
511 absVel += velocity[dimIdx]*velocity[dimIdx];
519 absVel = Toolbox::sqrt(absVel);
522 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
523 residual[dimIdx] +=
sqrtK_[dimIdx] * alpha * velocity[dimIdx];
525 Valgrind::CheckDefined(residual);
529 DimEvalMatrix& gradResid,
536 constexpr Scalar eps = 1e-11;
538 for (
unsigned i = 0; i < dimWorld; ++i) {
539 const Scalar coordEps = std::max(eps, Toolbox::scalarValue(velocity[i]) * (1 + eps));
540 velocity[i] += coordEps;
545 velocity[i] -= coordEps;
558 for (
unsigned i = 0; i < dimWorld; i++) {
559 for (
unsigned j = 0; j < dimWorld; j++) {
564 if (std::abs(K[i][j]) > 1e-25) {
573 Implementation& asImp_()
574 {
return *
static_cast<Implementation *
>(
this); }
576 const Implementation& asImp_()
const
577 {
return *
static_cast<const Implementation *
>(
this); }
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
std::array< Evaluation, numPhases > volumeFlux_
Definition: darcyfluxmodule.hh:570
std::array< Evaluation, numPhases > mobility_
Definition: darcyfluxmodule.hh:563
short upstreamIndex_(unsigned phaseIdx) const
Definition: darcyfluxmodule.hh:179
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
DimMatrix K_
Definition: darcyfluxmodule.hh:560
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:220
void calculateGradients_(const ElementContext &elemCtx, unsigned faceIdx, unsigned timeIdx)
Definition: forchheimerfluxmodule.hh:257
void calculateForchheimerFlux_(unsigned phaseIdx)
Definition: forchheimerfluxmodule.hh:444
const Evaluation & ergunCoefficient() const
Return the Ergun coefficent at the face's integration point.
Definition: forchheimerfluxmodule.hh:244
void calculateBoundaryGradients_(const ElementContext &elemCtx, unsigned boundaryFaceIdx, unsigned timeIdx, const FluidState &fluidState)
Definition: forchheimerfluxmodule.hh:314
void forchheimerResid_(DimEvalVector &residual, unsigned phaseIdx) const
Definition: forchheimerfluxmodule.hh:476
std::array< Evaluation, numPhases > mobilityPassabilityRatio_
Definition: forchheimerfluxmodule.hh:587
DimVector sqrtK_
Definition: forchheimerfluxmodule.hh:581
void gradForchheimerResid_(DimEvalVector &residual, DimEvalMatrix &gradResid, unsigned phaseIdx)
Definition: forchheimerfluxmodule.hh:528
void calculateBoundaryFluxes_(const ElementContext &elemCtx, unsigned bfIdx, unsigned timeIdx)
Calculate the volumetric flux rates of all phases at the domain boundary.
Definition: forchheimerfluxmodule.hh:417
Evaluation ergunCoefficient_
Definition: forchheimerfluxmodule.hh:584
void calculateFluxes_(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Calculate the volumetric fluxes of all phases.
Definition: forchheimerfluxmodule.hh:367
std::array< Evaluation, numPhases > density_
Definition: forchheimerfluxmodule.hh:590
bool isDiagonal_(const DimMatrix &K) const
Check whether all off-diagonal entries of a tensor are zero.
Definition: forchheimerfluxmodule.hh:556
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:253
Provides the intensive quantities for the Forchheimer module.
Definition: forchheimerfluxmodule.hh:134
const Evaluation & ergunCoefficient() const
Returns the Ergun coefficient.
Definition: forchheimerfluxmodule.hh:149
void update_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: forchheimerfluxmodule.hh:159
const Evaluation & mobilityPassabilityRatio(unsigned phaseIdx) const
Returns the passability of a phase.
Definition: forchheimerfluxmodule.hh:155
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: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
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