30#ifndef EWOMS_DARCY_FLUX_MODULE_HH
31#define EWOMS_DARCY_FLUX_MODULE_HH
33#include <dune/common/fmatrix.hh>
34#include <dune/common/fvector.hh>
36#include <opm/common/Exceptions.hpp>
38#include <opm/material/common/Valgrind.hpp>
52template <
class TypeTag>
53class DarcyIntensiveQuantities;
55template <
class TypeTag>
56class DarcyExtensiveQuantities;
58template <
class TypeTag>
59class DarcyBaseProblem;
65template <
class TypeTag>
84template <
class TypeTag>
92template <
class TypeTag>
98 void update_(
const ElementContext&,
120template <
class TypeTag>
131 enum { dimWorld = GridView::dimensionworld };
135 using ParameterCache =
typename FluidSystem::template ParameterCache<Evaluation>;
136 using EvalDimVector = Dune::FieldVector<Evaluation, dimWorld>;
137 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
138 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
155 {
return potentialGrad_[
phaseIdx]; }
164 {
return filterVelocity_[
phaseIdx]; }
179 short upstreamIndex_(
unsigned phaseIdx)
const
180 {
return upstreamDofIdx_[
phaseIdx]; }
182 short downstreamIndex_(
unsigned phaseIdx)
const
183 {
return downstreamDofIdx_[
phaseIdx]; }
194 const auto&
gradCalc = elemCtx.gradientCalculator();
197 const auto&
scvf = elemCtx.stencil(
timeIdx).interiorFace(faceIdx);
200 const unsigned i =
scvf.interiorIndex();
201 const unsigned j =
scvf.exteriorIndex();
202 interiorDofIdx_ =
static_cast<short>(i);
203 exteriorDofIdx_ =
static_cast<short>(j);
204 const unsigned focusDofIdx = elemCtx.focusDofIndex();
208 if (!elemCtx.model().phaseIsConsidered(
phaseIdx)) {
209 Valgrind::SetUndefined(potentialGrad_[
phaseIdx]);
218 Valgrind::CheckDefined(potentialGrad_[
phaseIdx]);
222 if (Parameters::Get<Parameters::EnableGravity>()) {
225 const auto&
gIn = elemCtx.problem().gravity(elemCtx, i,
timeIdx);
226 const auto&
gEx = elemCtx.problem().gravity(elemCtx, j,
timeIdx);
245 if (!elemCtx.model().phaseIsConsidered(
phaseIdx)) {
252 if (std::is_same<Scalar, Evaluation>::value ||
267 if (std::is_same<Scalar, Evaluation>::value ||
293 + std::string(FluidSystem::phaseName(
phaseIdx)) +
"'");
299 Valgrind::SetUndefined(K_);
300 elemCtx.problem().intersectionIntrinsicPermeability(K_, elemCtx, faceIdx,
timeIdx);
301 Valgrind::CheckDefined(K_);
304 if (!elemCtx.model().phaseIsConsidered(
phaseIdx)) {
305 Valgrind::SetUndefined(potentialGrad_[
phaseIdx]);
310 Evaluation tmp = 0.0;
316 upstreamDofIdx_[
phaseIdx] = exteriorDofIdx_;
317 downstreamDofIdx_[
phaseIdx] = interiorDofIdx_;
320 upstreamDofIdx_[
phaseIdx] = interiorDofIdx_;
321 downstreamDofIdx_[
phaseIdx] = exteriorDofIdx_;
326 const auto&
up = elemCtx.intensiveQuantities(upstreamDofIdx_[
phaseIdx],
timeIdx);
342 template <
class Flu
idState>
346 const FluidState& fluidState)
348 const auto&
gradCalc = elemCtx.gradientCalculator();
353 if (!elemCtx.model().phaseIsConsidered(
phaseIdx)) {
354 Valgrind::SetUndefined(potentialGrad_[
phaseIdx]);
363 Valgrind::CheckDefined(potentialGrad_[
phaseIdx]);
367 const auto i =
scvf.interiorIndex();
368 interiorDofIdx_ =
static_cast<short>(i);
369 exteriorDofIdx_ = -1;
377 if (Parameters::Get<Parameters::EnableGravity>()) {
380 const auto&
gIn = elemCtx.problem().gravity(elemCtx, i,
timeIdx);
391 if (!elemCtx.model().phaseIsConsidered(
phaseIdx)) {
399 Valgrind::CheckDefined(
pStatIn);
413 Valgrind::CheckDefined(potentialGrad_[
phaseIdx]);
417 + std::string(FluidSystem::phaseName(
phaseIdx)) +
"'");
426 const auto&
matParams = elemCtx.problem().materialLawParams(elemCtx, i,
timeIdx);
428 std::array<Scalar, numPhases>
kr;
429 MaterialLaw::relativePermeabilities(
kr,
matParams, fluidState);
430 Valgrind::CheckDefined(
kr);
433 if (!elemCtx.model().phaseIsConsidered(
phaseIdx)) {
437 Evaluation tmp = 0.0;
443 upstreamDofIdx_[
phaseIdx] = exteriorDofIdx_;
444 downstreamDofIdx_[
phaseIdx] = interiorDofIdx_;
447 upstreamDofIdx_[
phaseIdx] = interiorDofIdx_;
448 downstreamDofIdx_[
phaseIdx] = exteriorDofIdx_;
452 if (upstreamDofIdx_[
phaseIdx] < 0) {
458 Toolbox::value(fluidState.viscosity(
phaseIdx));
467 Valgrind::CheckDefined(mobility_[
phaseIdx]);
480 const DimVector& normal =
scvf.normal();
481 Valgrind::CheckDefined(normal);
486 if (!elemCtx.model().phaseIsConsidered(
phaseIdx)) {
490 asImp_().calculateFilterVelocity_(
phaseIdx);
491 Valgrind::CheckDefined(filterVelocity_[
phaseIdx]);
494 for (
unsigned i = 0; i < normal.size(); ++i) {
511 const DimVector& normal =
scvf.normal();
512 Valgrind::CheckDefined(normal);
515 if (!elemCtx.model().phaseIsConsidered(
phaseIdx)) {
521 asImp_().calculateFilterVelocity_(
phaseIdx);
522 Valgrind::CheckDefined(filterVelocity_[
phaseIdx]);
524 for (
unsigned i = 0; i < normal.size(); ++i) {
530 void calculateFilterVelocity_(
unsigned phaseIdx)
534 for (
unsigned i = 0; i < K_.M(); ++i) {
535 for (
unsigned j = 0; j < K_.N(); ++j) {
536 assert(std::isfinite(K_[i][j]));
545 for (
unsigned i = 0; i < filterVelocity_[
phaseIdx].size(); ++i) {
552 Implementation& asImp_()
553 {
return *
static_cast<Implementation*
>(
this); }
555 const Implementation& asImp_()
const
556 {
return *
static_cast<const Implementation*
>(
this); }
563 std::array<Evaluation, numPhases> mobility_;
566 std::array<EvalDimVector, numPhases> filterVelocity_;
570 std::array<Evaluation, numPhases> volumeFlux_;
573 std::array<EvalDimVector, numPhases> potentialGrad_;
576 std::array<short, numPhases> upstreamDofIdx_;
577 std::array<short, numPhases> downstreamDofIdx_;
578 short interiorDofIdx_;
579 short exteriorDofIdx_;
Callback class for a phase pressure.
Definition quantitycallbacks.hh:134
Provides the defaults for the parameters required by the Darcy velocity approach.
Definition darcyfluxmodule.hh:86
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:343
const EvalDimVector & filterVelocity(unsigned phaseIdx) const
Return the filter velocity of a fluid phase at the face's integration point [m/s].
Definition darcyfluxmodule.hh:163
void calculateBoundaryFluxes_(const ElementContext &elemCtx, unsigned boundaryFaceIdx, unsigned timeIdx)
Calculate the volumetric fluxes at a boundary face of all fluid phases.
Definition darcyfluxmodule.hh:506
const EvalDimVector & potentialGrad(unsigned phaseIdx) const
Return the pressure potential gradient of a fluid phase at the face's integration point [Pa/m].
Definition darcyfluxmodule.hh:154
void calculateFluxes_(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Calculate the volumetric fluxes of all phases.
Definition darcyfluxmodule.hh:477
const Evaluation & volumeFlux(unsigned phaseIdx) const
Return the volume flux of a fluid phase at the face's integration point .
Definition darcyfluxmodule.hh:175
void calculateGradients_(const ElementContext &elemCtx, unsigned faceIdx, unsigned timeIdx)
Calculate the gradients which are required to determine the volumetric fluxes.
Definition darcyfluxmodule.hh:190
const DimMatrix & intrinsicPermability() const
Returns the intrinsic permeability tensor for a given sub-control volume face.
Definition darcyfluxmodule.hh:145
Provides the intensive quantities for the Darcy flux module.
Definition darcyfluxmodule.hh:94
Callback class for a phase pressure.
Definition quantitycallbacks.hh:85
Defines the common parameters for the porous medium multi-phase models.
Defines the common properties required by the porous medium multi-phase models.
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
This method contains all callback classes for quantities that are required by some extensive quantiti...
Specifies a flux module which uses the Darcy relation.
Definition darcyfluxmodule.hh:67
static void registerParameters()
Register all run-time parameters for the flux module.
Definition darcyfluxmodule.hh:75