discretefractureextensivequantities.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef EWOMS_DISCRETE_FRACTURE_EXTENSIVE_QUANTITIES_HH
29#define EWOMS_DISCRETE_FRACTURE_EXTENSIVE_QUANTITIES_HH
30
31#include <dune/common/fmatrix.hh>
32#include <dune/common/fvector.hh>
33
35
36#include <array>
37#include <cassert>
38
39namespace Opm {
40
47template <class TypeTag>
49{
51
56
57 enum { dimWorld = GridView::dimensionworld };
58 enum { numPhases = FluidSystem::numPhases };
59
60 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
61 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
62
63public:
67 void update(const ElementContext& elemCtx, unsigned scvfIdx, unsigned timeIdx)
68 {
69 ParentType::update(elemCtx, scvfIdx, timeIdx);
70
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();
76
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)) {
81 // do nothing if no fracture goes though the current edge
82 return;
83 }
84
85 // average the intrinsic permeability of the fracture
86 elemCtx.problem().fractureFaceIntrinsicPermeability(fractureIntrinsicPermeability_,
87 elemCtx, scvfIdx, timeIdx);
88
89 auto distDirection = elemCtx.pos(outsideScvIdx, timeIdx);
90 distDirection -= elemCtx.pos(insideScvIdx, timeIdx);
91 distDirection /= distDirection.two_norm();
92
93 const auto& problem = elemCtx.problem();
94 fractureWidth_ = problem.fractureWidth(elemCtx, insideScvIdx,
95 outsideScvIdx, timeIdx);
96 assert(fractureWidth_ < scvf.area());
97
98 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
99 const auto& pGrad = extQuants.potentialGrad(phaseIdx);
100
101 const unsigned upstreamIdx = static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
102 const auto& up = elemCtx.intensiveQuantities(upstreamIdx, timeIdx);
103
104 // multiply with the fracture mobility of the upstream vertex
105 fractureIntrinsicPermeability_.mv(pGrad,
106 fractureFilterVelocity_[phaseIdx]);
107 fractureFilterVelocity_[phaseIdx] *= -up.fractureMobility(phaseIdx);
108
109 // divide the volume flux by two. This is required because
110 // a fracture is always shared by two sub-control-volume
111 // faces.
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();
117 }
118 }
119 }
120
121public:
122 const DimMatrix& fractureIntrinsicPermeability() const
123 { return fractureIntrinsicPermeability_; }
124
125 Scalar fractureVolumeFlux(unsigned phaseIdx) const
126 { return fractureVolumeFlux_[phaseIdx]; }
127
128 Scalar fractureWidth() const
129 { return fractureWidth_; }
130
131 const DimVector& fractureFilterVelocity(unsigned phaseIdx) const
132 { return fractureFilterVelocity_[phaseIdx]; }
133
134private:
135 DimMatrix fractureIntrinsicPermeability_;
136 std::array<DimVector, numPhases> fractureFilterVelocity_;
137 std::array<Scalar, numPhases> fractureVolumeFlux_;
138 Scalar fractureWidth_;
139};
140
141} // namespace Opm
142
143#endif
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