discretefractureintensivequantities.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_INTENSIVE_QUANTITIES_HH
29#define EWOMS_DISCRETE_FRACTURE_INTENSIVE_QUANTITIES_HH
30
31#include <opm/material/common/Valgrind.hpp>
32
34
35#include <algorithm>
36#include <array>
37
38namespace Opm {
39
48template <class TypeTag>
50{
57
58 enum { numPhases = FluidSystem::numPhases };
59 enum { dimWorld = GridView::dimensionworld };
60
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!");
65
66 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
67 enum { wettingPhaseIdx = MaterialLaw::wettingPhaseIdx };
68 enum { nonWettingPhaseIdx = MaterialLaw::nonWettingPhaseIdx };
69 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
70 using FluidState = ImmiscibleFluidState<Scalar, FluidSystem,
71 /*storeEnthalpy=*/enableEnergy>;
72
73public:
77
81 void update(const ElementContext& elemCtx, unsigned vertexIdx, unsigned timeIdx)
82 {
83 ParentType::update(elemCtx, vertexIdx, timeIdx);
84
85 const auto& problem = elemCtx.problem();
86 const auto& fractureMapper = problem.fractureMapper();
87 const unsigned globalVertexIdx = elemCtx.globalSpaceIndex(vertexIdx, timeIdx);
88
89 Valgrind::SetUndefined(fractureFluidState_);
90 Valgrind::SetUndefined(fractureVolume_);
91 Valgrind::SetUndefined(fracturePorosity_);
92 Valgrind::SetUndefined(fractureIntrinsicPermeability_);
93 Valgrind::SetUndefined(fractureRelativePermeabilities_);
94
95 // do nothing if there is no fracture within the current degree of freedom
96 if (!fractureMapper.isFractureVertex(globalVertexIdx)) {
98 return;
99 }
100
101 // Make sure that the wetting saturation in the matrix fluid
102 // state does not get larger than 1
103 const Scalar SwMatrix =
104 std::min<Scalar>(1.0, this->fluidState_.saturation(wettingPhaseIdx));
105 this->fluidState_.setSaturation(wettingPhaseIdx, SwMatrix);
106 this->fluidState_.setSaturation(nonWettingPhaseIdx, 1 - SwMatrix);
107
108 // retrieve the facture width and intrinsic permeability from
109 // the problem
111 problem.fracturePorosity(elemCtx, vertexIdx, timeIdx);
113 problem.fractureIntrinsicPermeability(elemCtx, vertexIdx, timeIdx);
114
115 // compute the fracture volume for the current sub-control
116 // volume. note, that we don't take overlaps of fractures into
117 // account for this.
118 fractureVolume_ = 0;
119 const auto& vertexPos = elemCtx.pos(vertexIdx, timeIdx);
120 for (unsigned vertex2Idx = 0; vertex2Idx < elemCtx.numDof(/*timeIdx=*/0); ++vertex2Idx) {
121 const unsigned globalVertex2Idx = elemCtx.globalSpaceIndex(vertex2Idx, timeIdx);
122
123 if (vertexIdx == vertex2Idx ||
124 !fractureMapper.isFractureEdge(globalVertexIdx, globalVertex2Idx))
125 {
126 continue;
127 }
128
129 const Scalar fractureWidth =
130 problem.fractureWidth(elemCtx, vertexIdx, vertex2Idx, timeIdx);
131
132 auto distVec = elemCtx.pos(vertex2Idx, timeIdx);
133 distVec -= vertexPos;
134
135 const Scalar edgeLength = distVec.two_norm();
136
137 // the fracture is always adjacent to two sub-control
138 // volumes of the control volume, so when calculating the
139 // volume of the fracture which gets attributed to one
140 // SCV, the fracture width needs to divided by 2. Also,
141 // only half of the edge is located in the current control
142 // volume, so its length also needs to divided by 2.
143 fractureVolume_ += (fractureWidth / 2) * (edgeLength / 2);
144 }
145
147 // set the fluid state for the fracture.
149
150 // start with the same fluid state as in the matrix. This
151 // implies equal saturations, pressures, temperatures,
152 // enthalpies, etc.
153 fractureFluidState_.assign(this->fluidState_);
154
155 // ask the problem for the material law parameters of the
156 // fracture.
157 const auto& fractureMatParams =
158 problem.fractureMaterialLawParams(elemCtx, vertexIdx, timeIdx);
159
160 // calculate the fracture saturations which would be required
161 // to be consistent with the pressures
162 std::array<Scalar, numPhases> saturations;
163 MaterialLaw::saturations(saturations, fractureMatParams,
165 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
166 fractureFluidState_.setSaturation(phaseIdx, saturations[phaseIdx]);
167 }
168
169 // Make sure that the wetting saturation in the fracture does
170 // not get negative
171 const Scalar SwFracture =
172 std::max<Scalar>(0.0, fractureFluidState_.saturation(wettingPhaseIdx));
173 fractureFluidState_.setSaturation(wettingPhaseIdx, SwFracture);
174 fractureFluidState_.setSaturation(nonWettingPhaseIdx, 1 - SwFracture);
175
176 // calculate the relative permeabilities of the fracture
177 MaterialLaw::relativePermeabilities(fractureRelativePermeabilities_,
178 fractureMatParams,
180
181 // make sure that valgrind complains if the fluid state is not
182 // fully defined.
183 fractureFluidState_.checkDefined();
184 }
185
186public:
193 Scalar fractureRelativePermeability(unsigned phaseIdx) const
194 { return fractureRelativePermeabilities_[phaseIdx]; }
195
202 Scalar fractureMobility(unsigned phaseIdx) const
203 {
204 return fractureRelativePermeabilities_[phaseIdx] /
205 fractureFluidState_.viscosity(phaseIdx);
206 }
207
211 Scalar fracturePorosity() const
212 { return fracturePorosity_; }
213
218 const DimMatrix& fractureIntrinsicPermeability() const
220
225 Scalar fractureVolume() const
226 { return fractureVolume_; }
227
232 const FluidState& fractureFluidState() const
233 { return fractureFluidState_; }
234
235protected:
240 std::array<Scalar, numPhases> fractureRelativePermeabilities_;
241};
242
243} // namespace Opm
244
245#endif
Contains the quantities which are are constant within a finite volume in the discret fracture immisci...
Definition: discretefractureintensivequantities.hh:50
std::array< Scalar, numPhases > fractureRelativePermeabilities_
Definition: discretefractureintensivequantities.hh:240
Scalar fracturePorosity() const
Returns the average porosity within the fracture.
Definition: discretefractureintensivequantities.hh:211
Scalar fractureVolume_
Definition: discretefractureintensivequantities.hh:237
DimMatrix fractureIntrinsicPermeability_
Definition: discretefractureintensivequantities.hh:239
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
DiscreteFractureIntensiveQuantities(const DiscreteFractureIntensiveQuantities &)=default
DiscreteFractureIntensiveQuantities & operator=(const DiscreteFractureIntensiveQuantities &)=default
FluidState fractureFluidState_
Definition: discretefractureintensivequantities.hh:236
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
Scalar fracturePorosity_
Definition: discretefractureintensivequantities.hh:238
Contains the quantities which are are constant within a finite volume for the immiscible multi-phase ...
Definition: immiscibleintensivequantities.hh:54
FluidState fluidState_
Definition: immiscibleintensivequantities.hh:190
void update(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: immiscibleintensivequantities.hh:93
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