28#ifndef EWOMS_DISCRETE_FRACTURE_EXTENSIVE_QUANTITIES_HH
29#define EWOMS_DISCRETE_FRACTURE_EXTENSIVE_QUANTITIES_HH
31#include <dune/common/fmatrix.hh>
32#include <dune/common/fvector.hh>
47template <
class TypeTag>
57 enum { dimWorld = GridView::dimensionworld };
58 enum { numPhases = FluidSystem::numPhases };
60 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
61 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
67 void update(
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
71 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
72 const auto& stencil = elemCtx.stencil(timeIdx);
73 const auto& scvf = stencil.interiorFace(scvfIdx);
74 const unsigned insideScvIdx = scvf.interiorIndex();
75 const unsigned outsideScvIdx = scvf.exteriorIndex();
77 const unsigned globalI = elemCtx.globalSpaceIndex(insideScvIdx, timeIdx);
78 const unsigned globalJ = elemCtx.globalSpaceIndex(outsideScvIdx, timeIdx);
79 const auto& fractureMapper = elemCtx.problem().fractureMapper();
80 if (!fractureMapper.isFractureEdge(globalI, globalJ)) {
86 elemCtx.problem().fractureFaceIntrinsicPermeability(fractureIntrinsicPermeability_,
87 elemCtx, scvfIdx, timeIdx);
89 auto distDirection = elemCtx.pos(outsideScvIdx, timeIdx);
90 distDirection -= elemCtx.pos(insideScvIdx, timeIdx);
91 distDirection /= distDirection.two_norm();
93 const auto& problem = elemCtx.problem();
94 fractureWidth_ = problem.fractureWidth(elemCtx, insideScvIdx,
95 outsideScvIdx, timeIdx);
96 assert(fractureWidth_ < scvf.area());
98 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
99 const auto& pGrad = extQuants.potentialGrad(phaseIdx);
101 const unsigned upstreamIdx =
static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
102 const auto& up = elemCtx.intensiveQuantities(upstreamIdx, timeIdx);
105 fractureIntrinsicPermeability_.mv(pGrad,
106 fractureFilterVelocity_[phaseIdx]);
107 fractureFilterVelocity_[phaseIdx] *= -up.fractureMobility(phaseIdx);
112 fractureVolumeFlux_[phaseIdx] = 0;
113 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
114 fractureVolumeFlux_[phaseIdx] +=
115 (fractureFilterVelocity_[phaseIdx][dimIdx] * distDirection[dimIdx])
116 * (fractureWidth_ / 2.0) / scvf.area();
123 {
return fractureIntrinsicPermeability_; }
126 {
return fractureVolumeFlux_[phaseIdx]; }
129 {
return fractureWidth_; }
132 {
return fractureFilterVelocity_[phaseIdx]; }
135 DimMatrix fractureIntrinsicPermeability_;
136 std::array<DimVector, numPhases> fractureFilterVelocity_;
137 std::array<Scalar, numPhases> fractureVolumeFlux_;
138 Scalar fractureWidth_;
This class expresses all intensive quantities of the discrete fracture model.
Definition: discretefractureextensivequantities.hh:49
Scalar fractureVolumeFlux(unsigned phaseIdx) const
Definition: discretefractureextensivequantities.hh:125
Scalar fractureWidth() const
Definition: discretefractureextensivequantities.hh:128
const DimVector & fractureFilterVelocity(unsigned phaseIdx) const
Definition: discretefractureextensivequantities.hh:131
const DimMatrix & fractureIntrinsicPermeability() const
Definition: discretefractureextensivequantities.hh:122
void update(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Update the extensive quantities for a given sub-control-volume-face.
Definition: discretefractureextensivequantities.hh:67
This class provides the data all quantities that are required to calculate the fluxes of the fluid ph...
Definition: immiscibleextensivequantities.hh:53
void update(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Update the extensive quantities for a given sub-control-volume-face.
Definition: immiscibleextensivequantities.hh:76
Definition: blackoilboundaryratevector.hh:39
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