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>
132template <
class TypeTag>
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);
164 mobilityPassabilityRatio_[
phaseIdx] =
165 problem.mobilityPassabilityRatio(elemCtx,
173 Evaluation ergunCoefficient_;
174 std::array<Evaluation, numPhases> mobilityPassabilityRatio_;
216template <
class TypeTag>
227 enum { dimWorld = GridView::dimensionworld };
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>;
242 {
return ergunCoefficient_; }
251 {
return mobilityPassabilityRatio_[
phaseIdx]; }
254 void calculateGradients_(
const ElementContext& elemCtx,
261 const unsigned i =
static_cast<unsigned>(this->interiorDofIdx_);
262 const unsigned j =
static_cast<unsigned>(this->exteriorDofIdx_);
292 if (!elemCtx.model().phaseIsConsidered(
phaseIdx)) {
296 const unsigned upIdx =
static_cast<unsigned>(this->upstreamIndex_(
phaseIdx));
310 template <
class Flu
idState>
311 void calculateBoundaryGradients_(
const ElementContext& elemCtx,
314 const FluidState& fluidState)
322 const unsigned i =
static_cast<unsigned>(this->interiorDofIdx_);
328 ergunCoefficient_ =
intQuantsIn.ergunCoefficient();
342 if (!elemCtx.model().phaseIsConsidered(
phaseIdx)) {
352 mobilityPassabilityRatio_[
phaseIdx] =
367 const auto i = asImp_().interiorIndex();
368 const auto j = asImp_().exteriorIndex();
373 const auto& normal =
scvf.normal();
374 Valgrind::CheckDefined(normal);
379 ergunCoefficient_ = (
intQuantsI.ergunCoefficient() +
395 if (!elemCtx.model().phaseIsConsidered(
phaseIdx)) {
396 this->filterVelocity_[
phaseIdx] = 0.0;
401 calculateForchheimerFlux_(
phaseIdx);
418 const auto& boundaryFace = elemCtx.stencil(
timeIdx).boundaryFace(
bfIdx);
419 const auto& normal = boundaryFace.normal();
425 if (!elemCtx.model().phaseIsConsidered(
phaseIdx)) {
426 this->filterVelocity_[
phaseIdx] = 0.0;
431 calculateForchheimerFlux_(
phaseIdx);
441 void calculateForchheimerFlux_(
unsigned phaseIdx)
444 DimEvalVector& velocity = this->filterVelocity_[
phaseIdx];
451 DimEvalVector residual;
457 while (
deltaV.one_norm() > 1
e-11) {
460 + std::to_string(
newtonIter) +
" iterations");
473 void forchheimerResid_(DimEvalVector& residual,
unsigned phaseIdx)
const
475 const DimEvalVector& velocity = this->filterVelocity_[
phaseIdx];
478 const auto& mobility = this->mobility_[
phaseIdx];
479 const auto& density = density_[
phaseIdx];
522 Valgrind::CheckDefined(residual);
525 void gradForchheimerResid_(DimEvalVector& residual,
530 DimEvalVector& velocity = this->filterVelocity_[
phaseIdx];
531 forchheimerResid_(residual,
phaseIdx);
533 constexpr Scalar eps = 1
e-11;
535 for (
unsigned i = 0; i < dimWorld; ++i) {
536 const Scalar
coordEps = std::max(eps, Toolbox::scalarValue(velocity[i]) * (1 + eps));
555 for (
unsigned i = 0; i < dimWorld; i++) {
556 for (
unsigned j = 0; j < dimWorld; j++) {
561 if (std::abs(K[i][j]) > 1
e-25) {
570 Implementation& asImp_()
571 {
return *
static_cast<Implementation *
>(
this); }
573 const Implementation& asImp_()
const
574 {
return *
static_cast<const Implementation *
>(
this); }
581 Evaluation ergunCoefficient_;
584 std::array<Evaluation, numPhases> mobilityPassabilityRatio_;
587 std::array<Evaluation, numPhases> density_;
Provides the Darcy flux module.
Definition darcyfluxmodule.hh:122
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
void calculateGradients_(const ElementContext &elemCtx, unsigned faceIdx, unsigned timeIdx)
Calculate the gradients which are required to determine the volumetric fluxes.
Definition darcyfluxmodule.hh:189
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
const Evaluation & ergunCoefficient() const
Return the Ergun coefficent at the face's integration point.
Definition forchheimerfluxmodule.hh:241
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
void calculateFluxes_(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Calculate the volumetric fluxes of all phases.
Definition forchheimerfluxmodule.hh:364
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
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.
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
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
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