28#ifndef EWOMS_DISCRETE_FRACTURE_EXTENSIVE_QUANTITIES_HH
29#define EWOMS_DISCRETE_FRACTURE_EXTENSIVE_QUANTITIES_HH
33#include <dune/common/fvector.hh>
34#include <dune/common/fmatrix.hh>
44template <
class TypeTag>
54 enum { dimWorld = GridView::dimensionworld };
55 enum { numPhases = FluidSystem::numPhases };
57 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
58 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
64 void update(
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
68 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
69 const auto& stencil = elemCtx.stencil(timeIdx);
70 const auto& scvf = stencil.interiorFace(scvfIdx);
71 unsigned insideScvIdx = scvf.interiorIndex();
72 unsigned outsideScvIdx = scvf.exteriorIndex();
74 unsigned globalI = elemCtx.globalSpaceIndex(insideScvIdx, timeIdx);
75 unsigned globalJ = elemCtx.globalSpaceIndex(outsideScvIdx, timeIdx);
76 const auto& fractureMapper = elemCtx.problem().fractureMapper();
77 if (!fractureMapper.isFractureEdge(globalI, globalJ))
82 elemCtx.problem().fractureFaceIntrinsicPermeability(fractureIntrinsicPermeability_,
83 elemCtx, scvfIdx, timeIdx);
85 auto distDirection = elemCtx.pos(outsideScvIdx, timeIdx);
86 distDirection -= elemCtx.pos(insideScvIdx, timeIdx);
87 distDirection /= distDirection.two_norm();
89 const auto& problem = elemCtx.problem();
90 fractureWidth_ = problem.fractureWidth(elemCtx, insideScvIdx,
91 outsideScvIdx, timeIdx);
92 assert(fractureWidth_ < scvf.area());
94 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
95 const auto& pGrad = extQuants.potentialGrad(phaseIdx);
97 unsigned upstreamIdx =
static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
98 const auto& up = elemCtx.intensiveQuantities(upstreamIdx, timeIdx);
101 fractureIntrinsicPermeability_.mv(pGrad,
102 fractureFilterVelocity_[phaseIdx]);
103 fractureFilterVelocity_[phaseIdx] *= -up.fractureMobility(phaseIdx);
108 fractureVolumeFlux_[phaseIdx] = 0;
109 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx)
110 fractureVolumeFlux_[phaseIdx] +=
111 (fractureFilterVelocity_[phaseIdx][dimIdx] * distDirection[dimIdx])
112 * (fractureWidth_ / 2.0) / scvf.area();
118 {
return fractureIntrinsicPermeability_; }
121 {
return fractureVolumeFlux_[phaseIdx]; }
124 {
return fractureWidth_; }
127 {
return fractureFilterVelocity_[phaseIdx]; }
130 DimMatrix fractureIntrinsicPermeability_;
131 DimVector fractureFilterVelocity_[numPhases];
132 Scalar fractureVolumeFlux_[numPhases];
133 Scalar fractureWidth_;
This class expresses all intensive quantities of the discrete fracture model.
Definition: discretefractureextensivequantities.hh:46
Scalar fractureVolumeFlux(unsigned phaseIdx) const
Definition: discretefractureextensivequantities.hh:120
Scalar fractureWidth() const
Definition: discretefractureextensivequantities.hh:123
const DimVector & fractureFilterVelocity(unsigned phaseIdx) const
Definition: discretefractureextensivequantities.hh:126
const DimMatrix & fractureIntrinsicPermeability() const
Definition: discretefractureextensivequantities.hh:117
void update(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Update the extensive quantities for a given sub-control-volume-face.
Definition: discretefractureextensivequantities.hh:64
This class provides the data all quantities that are required to calculate the fluxes of the fluid ph...
Definition: immiscibleextensivequantities.hh:52
void update(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Update the extensive quantities for a given sub-control-volume-face.
Definition: immiscibleextensivequantities.hh:75
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