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 EvalDimVector = Dune::FieldVector<Evaluation, dimWorld>;
136 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
137 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
193 const auto& gradCalc = elemCtx.gradientCalculator();
196 const auto& scvf = elemCtx.stencil(timeIdx).interiorFace(faceIdx);
197 const auto& faceNormal = scvf.normal();
199 const unsigned i = scvf.interiorIndex();
200 const unsigned j = scvf.exteriorIndex();
203 const unsigned focusDofIdx = elemCtx.focusDofIndex();
206 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
207 if (!elemCtx.model().phaseIsConsidered(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);
227 const auto& intQuantsIn = elemCtx.intensiveQuantities(i, timeIdx);
228 const auto& intQuantsEx = elemCtx.intensiveQuantities(j, timeIdx);
230 const auto& posIn = elemCtx.pos(i, timeIdx);
231 const auto& posEx = elemCtx.pos(j, timeIdx);
232 const auto& posFace = scvf.integrationPos();
235 DimVector distVecIn(posIn);
236 DimVector distVecEx(posEx);
237 DimVector distVecTotal(posEx);
239 distVecIn -= posFace;
240 distVecEx -= posFace;
241 distVecTotal -= posIn;
242 const Scalar absDistTotalSquared = distVecTotal.two_norm2();
243 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
244 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
251 if (std::is_same<Scalar, Evaluation>::value ||
254 const Evaluation& rhoIn = intQuantsIn.fluidState().density(phaseIdx);
255 pStatIn = -rhoIn * (gIn * distVecIn);
258 const Scalar rhoIn = Toolbox::value(intQuantsIn.fluidState().density(phaseIdx));
259 pStatIn = -rhoIn * (gIn * distVecIn);
266 if (std::is_same<Scalar, Evaluation>::value ||
269 const Evaluation& rhoEx = intQuantsEx.fluidState().density(phaseIdx);
270 pStatEx = -rhoEx * (gEx * distVecEx);
273 const Scalar rhoEx = Toolbox::value(intQuantsEx.fluidState().density(phaseIdx));
274 pStatEx = -rhoEx * (gEx * distVecEx);
281 Dune::FieldVector<Evaluation, dimWorld> f(distVecTotal);
282 f *= (pStatEx - pStatIn) / absDistTotalSquared;
285 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
289 for (
unsigned dimIdx = 0; dimIdx <
potentialGrad_[phaseIdx].size(); ++dimIdx) {
291 throw NumericalProblem(
"Non-finite potential gradient for phase '"
292 + std::string(FluidSystem::phaseName(phaseIdx)) +
"'");
298 Valgrind::SetUndefined(
K_);
299 elemCtx.problem().intersectionIntrinsicPermeability(
K_, elemCtx, faceIdx, timeIdx);
300 Valgrind::CheckDefined(
K_);
302 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
303 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
309 Evaluation tmp = 0.0;
310 for (
unsigned dimIdx = 0; dimIdx < faceNormal.size(); ++dimIdx) {
325 const auto& up = elemCtx.intensiveQuantities(
upstreamDofIdx_[phaseIdx], timeIdx);
327 mobility_[phaseIdx] = up.mobility(phaseIdx);
330 mobility_[phaseIdx] = Toolbox::value(up.mobility(phaseIdx));
341 template <
class Flu
idState>
343 unsigned boundaryFaceIdx,
345 const FluidState& fluidState)
347 const auto& gradCalc = elemCtx.gradientCalculator();
351 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
352 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
365 const auto& scvf = elemCtx.stencil(timeIdx).boundaryFace(boundaryFaceIdx);
366 const auto i = scvf.interiorIndex();
369 int focusDofIdx = elemCtx.focusDofIndex();
372 const auto& intQuantsIn = elemCtx.intensiveQuantities(i, timeIdx);
373 K_ = intQuantsIn.intrinsicPermeability();
376 if (Parameters::Get<Parameters::EnableGravity>()) {
379 const auto& gIn = elemCtx.problem().gravity(elemCtx, i, timeIdx);
380 const auto& posIn = elemCtx.pos(i, timeIdx);
381 const auto& posFace = scvf.integrationPos();
384 DimVector distVecIn(posIn);
385 distVecIn -= posFace;
386 const Scalar absDistSquared = distVecIn.two_norm2();
387 const Scalar gTimesDist = gIn * distVecIn;
389 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
390 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
395 const Evaluation rhoIn = intQuantsIn.fluidState().density(phaseIdx);
396 const Evaluation pStatIn = -gTimesDist * rhoIn;
398 Valgrind::CheckDefined(pStatIn);
404 EvalDimVector f(distVecIn);
405 f *= pStatIn / absDistSquared;
408 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
413 for (
unsigned dimIdx = 0; dimIdx <
potentialGrad_[phaseIdx].size(); ++dimIdx) {
415 throw NumericalProblem(
"Non-finite potential gradient for phase '"
416 + std::string(FluidSystem::phaseName(phaseIdx)) +
"'");
423 const auto& faceNormal = scvf.normal();
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);
431 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
432 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
436 Evaluation tmp = 0.0;
437 for (
unsigned dimIdx = 0; dimIdx < faceNormal.size(); ++dimIdx) {
453 mobility_[phaseIdx] = kr[phaseIdx] / fluidState.viscosity(phaseIdx);
456 mobility_[phaseIdx] = Toolbox::value(kr[phaseIdx]) /
457 Toolbox::value(fluidState.viscosity(phaseIdx));
461 mobility_[phaseIdx] = Toolbox::value(intQuantsIn.mobility(phaseIdx));
464 mobility_[phaseIdx] = intQuantsIn.mobility(phaseIdx);
466 Valgrind::CheckDefined(
mobility_[phaseIdx]);
478 const auto& scvf = elemCtx.stencil(timeIdx).interiorFace(scvfIdx);
479 const DimVector& normal = scvf.normal();
480 Valgrind::CheckDefined(normal);
482 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
485 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
489 asImp_().calculateFilterVelocity_(phaseIdx);
493 for (
unsigned i = 0; i < normal.size(); ++i) {
506 unsigned boundaryFaceIdx,
509 const auto& scvf = elemCtx.stencil(timeIdx).boundaryFace(boundaryFaceIdx);
510 const DimVector& normal = scvf.normal();
511 Valgrind::CheckDefined(normal);
513 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
514 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
520 asImp_().calculateFilterVelocity_(phaseIdx);
523 for (
unsigned i = 0; i < normal.size(); ++i) {
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]));
551 Implementation& asImp_()
552 {
return *
static_cast<Implementation*
>(
this); }
554 const Implementation& asImp_()
const
555 {
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:163
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:565
std::array< EvalDimVector, numPhases > potentialGrad_
Definition: darcyfluxmodule.hh:572
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
short interiorDofIdx_
Definition: darcyfluxmodule.hh:577
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
std::array< Evaluation, numPhases > volumeFlux_
Definition: darcyfluxmodule.hh:569
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
std::array< Evaluation, numPhases > mobility_
Definition: darcyfluxmodule.hh:562
short upstreamIndex_(unsigned phaseIdx) const
Definition: darcyfluxmodule.hh:178
std::array< short, numPhases > upstreamDofIdx_
Definition: darcyfluxmodule.hh:575
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
short downstreamIndex_(unsigned phaseIdx) const
Definition: darcyfluxmodule.hh:181
std::array< short, numPhases > downstreamDofIdx_
Definition: darcyfluxmodule.hh:576
void calculateFilterVelocity_(unsigned phaseIdx)
Definition: darcyfluxmodule.hh:529
void calculateGradients_(const ElementContext &elemCtx, unsigned faceIdx, unsigned timeIdx)
Calculate the gradients which are required to determine the volumetric fluxes.
Definition: darcyfluxmodule.hh:189
short exteriorDofIdx_
Definition: darcyfluxmodule.hh:578
const DimMatrix & intrinsicPermability() const
Returns the intrinsic permeability tensor for a given sub-control volume face.
Definition: darcyfluxmodule.hh:144
DimMatrix K_
Definition: darcyfluxmodule.hh:559
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: blackoilbioeffectsmodules.hh:43
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