28#ifndef EWOMS_DISCRETE_FRACTURE_INTENSIVE_QUANTITIES_HH
29#define EWOMS_DISCRETE_FRACTURE_INTENSIVE_QUANTITIES_HH
31#include <opm/material/common/Valgrind.hpp>
48template <
class TypeTag>
58 enum { numPhases = FluidSystem::numPhases };
59 enum { dimWorld = GridView::dimensionworld };
61 static_assert(dimWorld == 2,
"The fracture module currently is only "
62 "implemented for the 2D case!");
63 static_assert(numPhases == 2,
"The fracture module currently is only "
64 "implemented for two fluid phases!");
67 enum { wettingPhaseIdx = MaterialLaw::wettingPhaseIdx };
68 enum { nonWettingPhaseIdx = MaterialLaw::nonWettingPhaseIdx };
69 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
85 const auto& problem = elemCtx.problem();
86 const auto& fractureMapper = problem.fractureMapper();
89 Valgrind::SetUndefined(fractureFluidState_);
90 Valgrind::SetUndefined(fractureVolume_);
91 Valgrind::SetUndefined(fracturePorosity_);
92 Valgrind::SetUndefined(fractureIntrinsicPermeability_);
93 Valgrind::SetUndefined(fractureRelativePermeabilities_);
104 std::min<Scalar>(1.0, this->fluidState_.saturation(wettingPhaseIdx));
105 this->fluidState_.setSaturation(wettingPhaseIdx,
SwMatrix);
106 this->fluidState_.setSaturation(nonWettingPhaseIdx, 1 -
SwMatrix);
112 fractureIntrinsicPermeability_ =
129 const Scalar fractureWidth =
143 fractureVolume_ += (fractureWidth / 2) * (
edgeLength / 2);
153 fractureFluidState_.assign(this->fluidState_);
164 fractureFluidState_);
172 std::max<Scalar>(0.0, fractureFluidState_.saturation(wettingPhaseIdx));
173 fractureFluidState_.setSaturation(wettingPhaseIdx,
SwFracture);
174 fractureFluidState_.setSaturation(nonWettingPhaseIdx, 1 -
SwFracture);
177 MaterialLaw::relativePermeabilities(fractureRelativePermeabilities_,
179 fractureFluidState_);
183 fractureFluidState_.checkDefined();
194 {
return fractureRelativePermeabilities_[
phaseIdx]; }
204 return fractureRelativePermeabilities_[
phaseIdx] /
205 fractureFluidState_.viscosity(
phaseIdx);
212 {
return fracturePorosity_; }
219 {
return fractureIntrinsicPermeability_; }
226 {
return fractureVolume_; }
233 {
return fractureFluidState_; }
236 FluidState fractureFluidState_;
237 Scalar fractureVolume_;
238 Scalar fracturePorosity_;
239 DimMatrix fractureIntrinsicPermeability_;
240 std::array<Scalar, numPhases> fractureRelativePermeabilities_;
Contains the quantities which are are constant within a finite volume in the discret fracture immisci...
Definition discretefractureintensivequantities.hh:50
Scalar fracturePorosity() const
Returns the average porosity within the fracture.
Definition discretefractureintensivequantities.hh:211
const FluidState & fractureFluidState() const
Returns a fluid state object which represents the thermodynamic state of the fluids within the fractu...
Definition discretefractureintensivequantities.hh:232
Scalar fractureMobility(unsigned phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition discretefractureintensivequantities.hh:202
Scalar fractureRelativePermeability(unsigned phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition discretefractureintensivequantities.hh:193
Scalar fractureVolume() const
Return the volume [m^2] occupied by fractures within the given sub-control volume.
Definition discretefractureintensivequantities.hh:225
const DimMatrix & fractureIntrinsicPermeability() const
Returns the average intrinsic permeability within the fracture.
Definition discretefractureintensivequantities.hh:218
void update(const ElementContext &elemCtx, unsigned vertexIdx, unsigned timeIdx)
Definition discretefractureintensivequantities.hh:81
Contains the quantities which are are constant within a finite volume for the immiscible multi-phase ...
Definition immiscibleintensivequantities.hh:54
void update(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition immiscibleintensivequantities.hh:93
Contains the quantities which are are constant within a finite volume for the immiscible multi-phase ...
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