30#ifndef EWOMS_DARCY_FLUX_MODULE_HH
31#define EWOMS_DARCY_FLUX_MODULE_HH
33#include <dune/common/fvector.hh>
34#include <dune/common/fmatrix.hh>
36#include <opm/common/Exceptions.hpp>
38#include <opm/material/common/Valgrind.hpp>
46template <
class TypeTag>
47class DarcyIntensiveQuantities;
49template <
class TypeTag>
50class DarcyExtensiveQuantities;
52template <
class TypeTag>
53class DarcyBaseProblem;
59template <
class TypeTag>
78template <
class TypeTag>
86template <
class TypeTag>
113template <
class TypeTag>
124 enum { dimWorld = GridView::dimensionworld };
125 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
127 using Toolbox = MathToolbox<Evaluation>;
128 using ParameterCache =
typename FluidSystem::template ParameterCache<Evaluation>;
129 using EvalDimVector = Dune::FieldVector<Evaluation, dimWorld>;
130 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
131 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
187 const auto& gradCalc = elemCtx.gradientCalculator();
190 const auto& scvf = elemCtx.stencil(timeIdx).interiorFace(faceIdx);
191 const auto& faceNormal = scvf.normal();
193 unsigned i = scvf.interiorIndex();
194 unsigned j = scvf.exteriorIndex();
197 unsigned focusDofIdx = elemCtx.focusDofIndex();
200 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
201 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
215 if (Parameters::Get<Parameters::EnableGravity>()) {
218 const auto& gIn = elemCtx.problem().gravity(elemCtx, i, timeIdx);
219 const auto& gEx = elemCtx.problem().gravity(elemCtx, j, timeIdx);
221 const auto& intQuantsIn = elemCtx.intensiveQuantities(i, timeIdx);
222 const auto& intQuantsEx = elemCtx.intensiveQuantities(j, timeIdx);
224 const auto& posIn = elemCtx.pos(i, timeIdx);
225 const auto& posEx = elemCtx.pos(j, timeIdx);
226 const auto& posFace = scvf.integrationPos();
229 DimVector distVecIn(posIn);
230 DimVector distVecEx(posEx);
231 DimVector distVecTotal(posEx);
233 distVecIn -= posFace;
234 distVecEx -= posFace;
235 distVecTotal -= posIn;
236 Scalar absDistTotalSquared = distVecTotal.two_norm2();
237 for (
unsigned phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
238 if (!elemCtx.model().phaseIsConsidered(phaseIdx))
244 if (std::is_same<Scalar, Evaluation>::value ||
247 const Evaluation& rhoIn = intQuantsIn.fluidState().density(phaseIdx);
248 pStatIn = - rhoIn*(gIn*distVecIn);
251 Scalar rhoIn = Toolbox::value(intQuantsIn.fluidState().density(phaseIdx));
252 pStatIn = - rhoIn*(gIn*distVecIn);
259 if (std::is_same<Scalar, Evaluation>::value ||
262 const Evaluation& rhoEx = intQuantsEx.fluidState().density(phaseIdx);
263 pStatEx = - rhoEx*(gEx*distVecEx);
266 Scalar rhoEx = Toolbox::value(intQuantsEx.fluidState().density(phaseIdx));
267 pStatEx = - rhoEx*(gEx*distVecEx);
274 Dune::FieldVector<Evaluation, dimWorld> f(distVecTotal);
275 f *= (pStatEx - pStatIn)/absDistTotalSquared;
278 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx)
281 for (
unsigned dimIdx = 0; dimIdx <
potentialGrad_[phaseIdx].size(); ++dimIdx) {
283 throw NumericalProblem(
"Non-finite potential gradient for phase '"
284 + std::string(FluidSystem::phaseName(phaseIdx))+
"'");
290 Valgrind::SetUndefined(
K_);
291 elemCtx.problem().intersectionIntrinsicPermeability(
K_, elemCtx, faceIdx, timeIdx);
292 Valgrind::CheckDefined(
K_);
294 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
295 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
301 Evaluation tmp = 0.0;
302 for (
unsigned dimIdx = 0; dimIdx < faceNormal.size(); ++dimIdx)
316 const auto& up = elemCtx.intensiveQuantities(
upstreamDofIdx_[phaseIdx], timeIdx);
318 mobility_[phaseIdx] = up.mobility(phaseIdx);
320 mobility_[phaseIdx] = Toolbox::value(up.mobility(phaseIdx));
330 template <
class Flu
idState>
332 unsigned boundaryFaceIdx,
334 const FluidState& fluidState)
336 const auto& gradCalc = elemCtx.gradientCalculator();
340 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
341 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
354 const auto& scvf = elemCtx.stencil(timeIdx).boundaryFace(boundaryFaceIdx);
355 auto i = scvf.interiorIndex();
358 int focusDofIdx = elemCtx.focusDofIndex();
361 const auto& intQuantsIn = elemCtx.intensiveQuantities(i, timeIdx);
362 K_ = intQuantsIn.intrinsicPermeability();
365 if (Parameters::Get<Parameters::EnableGravity>()) {
368 const auto& gIn = elemCtx.problem().gravity(elemCtx, i, timeIdx);
369 const auto& posIn = elemCtx.pos(i, timeIdx);
370 const auto& posFace = scvf.integrationPos();
373 DimVector distVecIn(posIn);
374 distVecIn -= posFace;
375 Scalar absDistSquared = distVecIn.two_norm2();
376 Scalar gTimesDist = gIn*distVecIn;
378 for (
unsigned phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
379 if (!elemCtx.model().phaseIsConsidered(phaseIdx))
383 Evaluation rhoIn = intQuantsIn.fluidState().density(phaseIdx);
384 Evaluation pStatIn = - gTimesDist*rhoIn;
386 Valgrind::CheckDefined(pStatIn);
392 EvalDimVector f(distVecIn);
393 f *= pStatIn / absDistSquared;
396 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx)
400 for (
unsigned dimIdx = 0; dimIdx <
potentialGrad_[phaseIdx].size(); ++dimIdx) {
402 throw NumericalProblem(
"Non-finite potential gradient for phase '"
403 + std::string(FluidSystem::phaseName(phaseIdx))+
"'");
410 const auto& faceNormal = scvf.normal();
412 const auto& matParams = elemCtx.problem().materialLawParams(elemCtx, i, timeIdx);
414 Scalar kr[numPhases];
415 MaterialLaw::relativePermeabilities(kr, matParams, fluidState);
416 Valgrind::CheckDefined(kr);
418 for (
unsigned phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
419 if (!elemCtx.model().phaseIsConsidered(phaseIdx))
422 Evaluation tmp = 0.0;
423 for (
unsigned dimIdx = 0; dimIdx < faceNormal.size(); ++dimIdx)
439 kr[phaseIdx] / fluidState.viscosity(phaseIdx);
442 Toolbox::value(kr[phaseIdx])
443 / Toolbox::value(fluidState.viscosity(phaseIdx));
446 mobility_[phaseIdx] = Toolbox::value(intQuantsIn.mobility(phaseIdx));
448 mobility_[phaseIdx] = intQuantsIn.mobility(phaseIdx);
449 Valgrind::CheckDefined(
mobility_[phaseIdx]);
461 const auto& scvf = elemCtx.stencil(timeIdx).interiorFace(scvfIdx);
462 const DimVector& normal = scvf.normal();
463 Valgrind::CheckDefined(normal);
465 for (
unsigned phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
468 if (!elemCtx.model().phaseIsConsidered(phaseIdx))
471 asImp_().calculateFilterVelocity_(phaseIdx);
475 for (
unsigned i = 0; i < normal.size(); ++i)
487 unsigned boundaryFaceIdx,
490 const auto& scvf = elemCtx.stencil(timeIdx).boundaryFace(boundaryFaceIdx);
491 const DimVector& normal = scvf.normal();
492 Valgrind::CheckDefined(normal);
494 for (
unsigned phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
495 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
501 asImp_().calculateFilterVelocity_(phaseIdx);
504 for (
unsigned i = 0; i < normal.size(); ++i)
513 for (
unsigned i = 0; i <
K_.M(); ++ i)
514 for (
unsigned j = 0; j <
K_.N(); ++ j)
515 assert(std::isfinite(
K_[i][j]));
528 Implementation& asImp_()
529 {
return *
static_cast<Implementation*
>(
this); }
531 const Implementation& asImp_()
const
532 {
return *
static_cast<const Implementation*
>(
this); }
Callback class for a phase pressure.
Definition: quantitycallbacks.hh:133
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:80
Provides the Darcy flux module.
Definition: darcyfluxmodule.hh:115
Evaluation mobility_[numPhases]
Definition: darcyfluxmodule.hh:539
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:331
short interiorDofIdx_
Definition: darcyfluxmodule.hh:554
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:156
void calculateBoundaryFluxes_(const ElementContext &elemCtx, unsigned boundaryFaceIdx, unsigned timeIdx)
Calculate the volumetric fluxes at a boundary face of all fluid phases.
Definition: darcyfluxmodule.hh:486
short upstreamIndex_(unsigned phaseIdx) const
Definition: darcyfluxmodule.hh:172
EvalDimVector potentialGrad_[numPhases]
Definition: darcyfluxmodule.hh:549
EvalDimVector filterVelocity_[numPhases]
Definition: darcyfluxmodule.hh:542
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:147
short upstreamDofIdx_[numPhases]
Definition: darcyfluxmodule.hh:552
void calculateFluxes_(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Calculate the volumetric fluxes of all phases.
Definition: darcyfluxmodule.hh:459
short downstreamDofIdx_[numPhases]
Definition: darcyfluxmodule.hh:553
const Evaluation & volumeFlux(unsigned phaseIdx) const
Return the volume flux of a fluid phase at the face's integration point .
Definition: darcyfluxmodule.hh:168
short downstreamIndex_(unsigned phaseIdx) const
Definition: darcyfluxmodule.hh:175
Evaluation volumeFlux_[numPhases]
Definition: darcyfluxmodule.hh:546
void calculateFilterVelocity_(unsigned phaseIdx)
Definition: darcyfluxmodule.hh:509
void calculateGradients_(const ElementContext &elemCtx, unsigned faceIdx, unsigned timeIdx)
Calculate the gradients which are required to determine the volumetric fluxes.
Definition: darcyfluxmodule.hh:183
short exteriorDofIdx_
Definition: darcyfluxmodule.hh:555
const DimMatrix & intrinsicPermability() const
Returns the intrinsic permeability tensor for a given sub-control volume face.
Definition: darcyfluxmodule.hh:138
DimMatrix K_
Definition: darcyfluxmodule.hh:536
Provides the intensive quantities for the Darcy flux module.
Definition: darcyfluxmodule.hh:88
void update_(const ElementContext &, unsigned, unsigned)
Definition: darcyfluxmodule.hh:91
Callback class for a phase pressure.
Definition: quantitycallbacks.hh:84
void setPhaseIndex(unsigned phaseIdx)
Set the index of the fluid phase for which the pressure should be returned.
Definition: quantitycallbacks.hh:108
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:37
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:235
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:61
static void registerParameters()
Register all run-time parameters for the flux module.
Definition: darcyfluxmodule.hh:69