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>
120template <
class TypeTag>
131 enum { dimWorld = GridView::dimensionworld };
132 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
134 using Toolbox = MathToolbox<Evaluation>;
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>;
194 const auto& gradCalc = elemCtx.gradientCalculator();
197 const auto& scvf = elemCtx.stencil(timeIdx).interiorFace(faceIdx);
198 const auto& faceNormal = scvf.normal();
200 const unsigned i = scvf.interiorIndex();
201 const unsigned j = scvf.exteriorIndex();
204 const unsigned focusDofIdx = elemCtx.focusDofIndex();
207 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
208 if (!elemCtx.model().phaseIsConsidered(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);
228 const auto& intQuantsIn = elemCtx.intensiveQuantities(i, timeIdx);
229 const auto& intQuantsEx = elemCtx.intensiveQuantities(j, timeIdx);
231 const auto& posIn = elemCtx.pos(i, timeIdx);
232 const auto& posEx = elemCtx.pos(j, timeIdx);
233 const auto& posFace = scvf.integrationPos();
236 DimVector distVecIn(posIn);
237 DimVector distVecEx(posEx);
238 DimVector distVecTotal(posEx);
240 distVecIn -= posFace;
241 distVecEx -= posFace;
242 distVecTotal -= posIn;
243 const Scalar absDistTotalSquared = distVecTotal.two_norm2();
244 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
245 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
252 if (std::is_same<Scalar, Evaluation>::value ||
255 const Evaluation& rhoIn = intQuantsIn.fluidState().density(phaseIdx);
256 pStatIn = -rhoIn * (gIn * distVecIn);
259 const Scalar rhoIn = Toolbox::value(intQuantsIn.fluidState().density(phaseIdx));
260 pStatIn = -rhoIn * (gIn * distVecIn);
267 if (std::is_same<Scalar, Evaluation>::value ||
270 const Evaluation& rhoEx = intQuantsEx.fluidState().density(phaseIdx);
271 pStatEx = -rhoEx * (gEx * distVecEx);
274 const Scalar rhoEx = Toolbox::value(intQuantsEx.fluidState().density(phaseIdx));
275 pStatEx = -rhoEx * (gEx * distVecEx);
282 Dune::FieldVector<Evaluation, dimWorld> f(distVecTotal);
283 f *= (pStatEx - pStatIn) / absDistTotalSquared;
286 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
290 for (
unsigned dimIdx = 0; dimIdx <
potentialGrad_[phaseIdx].size(); ++dimIdx) {
292 throw NumericalProblem(
"Non-finite potential gradient for phase '"
293 + std::string(FluidSystem::phaseName(phaseIdx)) +
"'");
299 Valgrind::SetUndefined(
K_);
300 elemCtx.problem().intersectionIntrinsicPermeability(
K_, elemCtx, faceIdx, timeIdx);
301 Valgrind::CheckDefined(
K_);
303 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
304 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
310 Evaluation tmp = 0.0;
311 for (
unsigned dimIdx = 0; dimIdx < faceNormal.size(); ++dimIdx) {
326 const auto& up = elemCtx.intensiveQuantities(
upstreamDofIdx_[phaseIdx], timeIdx);
328 mobility_[phaseIdx] = up.mobility(phaseIdx);
331 mobility_[phaseIdx] = Toolbox::value(up.mobility(phaseIdx));
342 template <
class Flu
idState>
344 unsigned boundaryFaceIdx,
346 const FluidState& fluidState)
348 const auto& gradCalc = elemCtx.gradientCalculator();
352 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
353 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
366 const auto& scvf = elemCtx.stencil(timeIdx).boundaryFace(boundaryFaceIdx);
367 const auto i = scvf.interiorIndex();
370 int focusDofIdx = elemCtx.focusDofIndex();
373 const auto& intQuantsIn = elemCtx.intensiveQuantities(i, timeIdx);
374 K_ = intQuantsIn.intrinsicPermeability();
377 if (Parameters::Get<Parameters::EnableGravity>()) {
380 const auto& gIn = elemCtx.problem().gravity(elemCtx, i, timeIdx);
381 const auto& posIn = elemCtx.pos(i, timeIdx);
382 const auto& posFace = scvf.integrationPos();
385 DimVector distVecIn(posIn);
386 distVecIn -= posFace;
387 const Scalar absDistSquared = distVecIn.two_norm2();
388 const Scalar gTimesDist = gIn * distVecIn;
390 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
391 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
396 const Evaluation rhoIn = intQuantsIn.fluidState().density(phaseIdx);
397 const Evaluation pStatIn = -gTimesDist * rhoIn;
399 Valgrind::CheckDefined(pStatIn);
405 EvalDimVector f(distVecIn);
406 f *= pStatIn / absDistSquared;
409 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
414 for (
unsigned dimIdx = 0; dimIdx <
potentialGrad_[phaseIdx].size(); ++dimIdx) {
416 throw NumericalProblem(
"Non-finite potential gradient for phase '"
417 + std::string(FluidSystem::phaseName(phaseIdx)) +
"'");
424 const auto& faceNormal = scvf.normal();
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);
432 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
433 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
437 Evaluation tmp = 0.0;
438 for (
unsigned dimIdx = 0; dimIdx < faceNormal.size(); ++dimIdx) {
454 mobility_[phaseIdx] = kr[phaseIdx] / fluidState.viscosity(phaseIdx);
457 mobility_[phaseIdx] = Toolbox::value(kr[phaseIdx]) /
458 Toolbox::value(fluidState.viscosity(phaseIdx));
462 mobility_[phaseIdx] = Toolbox::value(intQuantsIn.mobility(phaseIdx));
465 mobility_[phaseIdx] = intQuantsIn.mobility(phaseIdx);
467 Valgrind::CheckDefined(
mobility_[phaseIdx]);
479 const auto& scvf = elemCtx.stencil(timeIdx).interiorFace(scvfIdx);
480 const DimVector& normal = scvf.normal();
481 Valgrind::CheckDefined(normal);
483 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
486 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
490 asImp_().calculateFilterVelocity_(phaseIdx);
494 for (
unsigned i = 0; i < normal.size(); ++i) {
507 unsigned boundaryFaceIdx,
510 const auto& scvf = elemCtx.stencil(timeIdx).boundaryFace(boundaryFaceIdx);
511 const DimVector& normal = scvf.normal();
512 Valgrind::CheckDefined(normal);
514 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
515 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
521 asImp_().calculateFilterVelocity_(phaseIdx);
524 for (
unsigned i = 0; i < normal.size(); ++i) {
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]));
552 Implementation& asImp_()
553 {
return *
static_cast<Implementation*
>(
this); }
555 const Implementation& asImp_()
const
556 {
return *
static_cast<const Implementation*
>(
this); }
Callback class for a phase pressure.
Definition: quantitycallbacks.hh:134
void setPhaseIndex(unsigned phaseIdx)
Set the index of the fluid phase for which the pressure should be returned.
Definition: quantitycallbacks.hh:164
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
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
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
std::array< Evaluation, numPhases > volumeFlux_
Definition: darcyfluxmodule.hh:570
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
std::array< Evaluation, numPhases > mobility_
Definition: darcyfluxmodule.hh:563
short upstreamIndex_(unsigned phaseIdx) const
Definition: darcyfluxmodule.hh:179
std::array< short, numPhases > upstreamDofIdx_
Definition: darcyfluxmodule.hh:576
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
short downstreamIndex_(unsigned phaseIdx) const
Definition: darcyfluxmodule.hh:182
std::array< short, numPhases > downstreamDofIdx_
Definition: darcyfluxmodule.hh:577
void calculateFilterVelocity_(unsigned phaseIdx)
Definition: darcyfluxmodule.hh:530
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
const DimMatrix & intrinsicPermability() const
Returns the intrinsic permeability tensor for a given sub-control volume face.
Definition: darcyfluxmodule.hh:145
DimMatrix K_
Definition: darcyfluxmodule.hh:560
Provides the intensive quantities for the Darcy flux module.
Definition: darcyfluxmodule.hh:94
void update_(const ElementContext &, unsigned, unsigned)
Definition: darcyfluxmodule.hh:98
Callback class for a phase pressure.
Definition: quantitycallbacks.hh:85
void setPhaseIndex(unsigned phaseIdx)
Set the index of the fluid phase for which the pressure should be returned.
Definition: quantitycallbacks.hh:109
Defines the common parameters for the porous medium multi-phase models.
Defines the common properties required by the porous medium multi-phase models.
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
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