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 EvalDimVector = Dune::FieldVector<Evaluation, dimWorld>;
136 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
137 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
154 {
return potentialGrad_[
phaseIdx]; }
163 {
return filterVelocity_[
phaseIdx]; }
178 short upstreamIndex_(
unsigned phaseIdx)
const
179 {
return upstreamDofIdx_[
phaseIdx]; }
181 short downstreamIndex_(
unsigned phaseIdx)
const
182 {
return downstreamDofIdx_[
phaseIdx]; }
193 const auto&
gradCalc = elemCtx.gradientCalculator();
196 const auto&
scvf = elemCtx.stencil(
timeIdx).interiorFace(faceIdx);
199 const unsigned i =
scvf.interiorIndex();
200 const unsigned j =
scvf.exteriorIndex();
201 interiorDofIdx_ =
static_cast<short>(i);
202 exteriorDofIdx_ =
static_cast<short>(j);
203 const unsigned focusDofIdx = elemCtx.focusDofIndex();
207 if (!elemCtx.model().phaseIsConsidered(
phaseIdx)) {
208 Valgrind::SetUndefined(potentialGrad_[
phaseIdx]);
217 Valgrind::CheckDefined(potentialGrad_[
phaseIdx]);
221 if (Parameters::Get<Parameters::EnableGravity>()) {
224 const auto&
gIn = elemCtx.problem().gravity(elemCtx, i,
timeIdx);
225 const auto&
gEx = elemCtx.problem().gravity(elemCtx, j,
timeIdx);
244 if (!elemCtx.model().phaseIsConsidered(
phaseIdx)) {
251 if (std::is_same<Scalar, Evaluation>::value ||
266 if (std::is_same<Scalar, Evaluation>::value ||
292 + std::string(FluidSystem::phaseName(
phaseIdx)) +
"'");
298 Valgrind::SetUndefined(K_);
299 elemCtx.problem().intersectionIntrinsicPermeability(K_, elemCtx, faceIdx,
timeIdx);
300 Valgrind::CheckDefined(K_);
303 if (!elemCtx.model().phaseIsConsidered(
phaseIdx)) {
304 Valgrind::SetUndefined(potentialGrad_[
phaseIdx]);
309 Evaluation tmp = 0.0;
315 upstreamDofIdx_[
phaseIdx] = exteriorDofIdx_;
316 downstreamDofIdx_[
phaseIdx] = interiorDofIdx_;
319 upstreamDofIdx_[
phaseIdx] = interiorDofIdx_;
320 downstreamDofIdx_[
phaseIdx] = exteriorDofIdx_;
325 const auto&
up = elemCtx.intensiveQuantities(upstreamDofIdx_[
phaseIdx],
timeIdx);
341 template <
class Flu
idState>
345 const FluidState& fluidState)
347 const auto&
gradCalc = elemCtx.gradientCalculator();
352 if (!elemCtx.model().phaseIsConsidered(
phaseIdx)) {
353 Valgrind::SetUndefined(potentialGrad_[
phaseIdx]);
362 Valgrind::CheckDefined(potentialGrad_[
phaseIdx]);
366 const auto i =
scvf.interiorIndex();
367 interiorDofIdx_ =
static_cast<short>(i);
368 exteriorDofIdx_ = -1;
376 if (Parameters::Get<Parameters::EnableGravity>()) {
379 const auto&
gIn = elemCtx.problem().gravity(elemCtx, i,
timeIdx);
390 if (!elemCtx.model().phaseIsConsidered(
phaseIdx)) {
398 Valgrind::CheckDefined(
pStatIn);
412 Valgrind::CheckDefined(potentialGrad_[
phaseIdx]);
416 + std::string(FluidSystem::phaseName(
phaseIdx)) +
"'");
425 const auto&
matParams = elemCtx.problem().materialLawParams(elemCtx, i,
timeIdx);
427 std::array<Scalar, numPhases>
kr;
428 MaterialLaw::relativePermeabilities(
kr,
matParams, fluidState);
429 Valgrind::CheckDefined(
kr);
432 if (!elemCtx.model().phaseIsConsidered(
phaseIdx)) {
436 Evaluation tmp = 0.0;
442 upstreamDofIdx_[
phaseIdx] = exteriorDofIdx_;
443 downstreamDofIdx_[
phaseIdx] = interiorDofIdx_;
446 upstreamDofIdx_[
phaseIdx] = interiorDofIdx_;
447 downstreamDofIdx_[
phaseIdx] = exteriorDofIdx_;
451 if (upstreamDofIdx_[
phaseIdx] < 0) {
457 Toolbox::value(fluidState.viscosity(
phaseIdx));
466 Valgrind::CheckDefined(mobility_[
phaseIdx]);
479 const DimVector& normal =
scvf.normal();
480 Valgrind::CheckDefined(normal);
485 if (!elemCtx.model().phaseIsConsidered(
phaseIdx)) {
489 asImp_().calculateFilterVelocity_(
phaseIdx);
490 Valgrind::CheckDefined(filterVelocity_[
phaseIdx]);
493 for (
unsigned i = 0; i < normal.size(); ++i) {
510 const DimVector& normal =
scvf.normal();
511 Valgrind::CheckDefined(normal);
514 if (!elemCtx.model().phaseIsConsidered(
phaseIdx)) {
520 asImp_().calculateFilterVelocity_(
phaseIdx);
521 Valgrind::CheckDefined(filterVelocity_[
phaseIdx]);
523 for (
unsigned i = 0; i < normal.size(); ++i) {
529 void calculateFilterVelocity_(
unsigned phaseIdx)
533 for (
unsigned i = 0; i < K_.M(); ++i) {
534 for (
unsigned j = 0; j < K_.N(); ++j) {
535 assert(std::isfinite(K_[i][j]));
544 for (
unsigned i = 0; i < filterVelocity_[
phaseIdx].size(); ++i) {
551 Implementation& asImp_()
552 {
return *
static_cast<Implementation*
>(
this); }
554 const Implementation& asImp_()
const
555 {
return *
static_cast<const Implementation*
>(
this); }
562 std::array<Evaluation, numPhases> mobility_;
565 std::array<EvalDimVector, numPhases> filterVelocity_;
569 std::array<Evaluation, numPhases> volumeFlux_;
572 std::array<EvalDimVector, numPhases> potentialGrad_;
575 std::array<short, numPhases> upstreamDofIdx_;
576 std::array<short, numPhases> downstreamDofIdx_;
577 short interiorDofIdx_;
578 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:342
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:162
void calculateBoundaryFluxes_(const ElementContext &elemCtx, unsigned boundaryFaceIdx, unsigned timeIdx)
Calculate the volumetric fluxes at a boundary face of all fluid phases.
Definition darcyfluxmodule.hh:505
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:153
void calculateFluxes_(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Calculate the volumetric fluxes of all phases.
Definition darcyfluxmodule.hh:476
const Evaluation & volumeFlux(unsigned phaseIdx) const
Return the volume flux of a fluid phase at the face's integration point .
Definition darcyfluxmodule.hh:174
void calculateGradients_(const ElementContext &elemCtx, unsigned faceIdx, unsigned timeIdx)
Calculate the gradients which are required to determine the volumetric fluxes.
Definition darcyfluxmodule.hh:189
const DimMatrix & intrinsicPermability() const
Returns the intrinsic permeability tensor for a given sub-control volume face.
Definition darcyfluxmodule.hh:144
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