28#ifndef EWOMS_DISCRETE_FRACTURE_LOCAL_RESIDUAL_BASE_HH
29#define EWOMS_DISCRETE_FRACTURE_LOCAL_RESIDUAL_BASE_HH
41template <
class TypeTag>
52 enum { conti0EqIdx = Indices::conti0EqIdx };
53 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
54 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
68 const ElementContext& elemCtx,
71 unsigned phaseIdx)
const
73 EqVector phaseStorage(0.0);
76 const auto& problem = elemCtx.problem();
77 const auto& fractureMapper = problem.fractureMapper();
78 unsigned globalIdx = elemCtx.globalSpaceIndex(dofIdx, timeIdx);
80 if (!fractureMapper.isFractureVertex(globalIdx)) {
83 storage += phaseStorage;
87 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
88 const auto& scv = elemCtx.stencil(timeIdx).subControlVolume(dofIdx);
91 phaseStorage *= 1 - intQuants.fractureVolume()/scv.volume();
94 const auto& fsFracture = intQuants.fractureFluidState();
96 phaseStorage[conti0EqIdx + phaseIdx] +=
97 intQuants.fracturePorosity()*
98 fsFracture.saturation(phaseIdx) *
99 fsFracture.density(phaseIdx) *
100 intQuants.fractureVolume()/scv.volume();
102 EnergyModule::addFracturePhaseStorage(phaseStorage, intQuants, scv,
106 storage += phaseStorage;
113 const ElementContext& elemCtx,
115 unsigned timeIdx)
const
119 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
121 unsigned i = extQuants.interiorIndex();
122 unsigned j = extQuants.exteriorIndex();
123 unsigned I = elemCtx.globalSpaceIndex(i, timeIdx);
124 unsigned J = elemCtx.globalSpaceIndex(j, timeIdx);
125 const auto& fractureMapper = elemCtx.problem().fractureMapper();
126 if (!fractureMapper.isFractureEdge(I, J))
131 const auto& scvf = elemCtx.stencil(timeIdx).interiorFace(scvfIdx);
132 Scalar scvfArea = scvf.area();
135 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
136 if (!elemCtx.model().phaseIsConsidered(phaseIdx))
143 flux[conti0EqIdx + phaseIdx] *=
144 1 - extQuants.fractureWidth() / (2 * scvfArea);
147 unsigned upIdx =
static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
148 const auto& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
149 flux[conti0EqIdx + phaseIdx] +=
150 extQuants.fractureVolumeFlux(phaseIdx) * up.fractureFluidState().density(phaseIdx);
153 EnergyModule::handleFractureFlux(flux, elemCtx, scvfIdx, timeIdx);
Calculates the local residual of the discrete fracture immiscible multi-phase model.
Definition: discretefracturelocalresidual.hh:43
void addPhaseStorage(EqVector &storage, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx, unsigned phaseIdx) const
Adds the amount all conservation quantities (e.g. phase mass) within a single fluid phase.
Definition: discretefracturelocalresidual.hh:67
void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Evaluates the total mass flux of all conservation quantities over a face of a sub-control volume.
Definition: discretefracturelocalresidual.hh:112
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:50
Calculates the local residual of the immiscible multi-phase model.
Definition: immisciblelocalresidual.hh:46
void addPhaseStorage(Dune::FieldVector< LhsEval, numEq > &storage, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx, unsigned phaseIdx) const
Adds the amount all conservation quantities (e.g. phase mass) within a single fluid phase.
Definition: immisciblelocalresidual.hh:76
void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Evaluates the total mass flux of all conservation quantities over a face of a sub-control volume.
Definition: immisciblelocalresidual.hh:114
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