My Project
Loading...
Searching...
No Matches
LinearMaterial.hpp
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*/
27#ifndef OPM_LINEAR_MATERIAL_HPP
28#define OPM_LINEAR_MATERIAL_HPP
29
31
34
35#include <type_traits>
36
37namespace Opm {
38
49template <class TraitsT, class ParamsT = LinearMaterialParams<TraitsT> >
50class LinearMaterial : public TraitsT
51{
52public:
53 typedef TraitsT Traits;
54 typedef ParamsT Params;
55 typedef typename Traits::Scalar Scalar;
56
58 using Traits::numPhases;
59
62 static const bool implementsTwoPhaseApi = (numPhases == 2);
63
66 static const bool implementsTwoPhaseSatApi = (numPhases == 2);
67
70 static const bool isSaturationDependent = true;
71
74 static const bool isPressureDependent = false;
75
78 static const bool isTemperatureDependent = false;
79
82 static const bool isCompositionDependent = false;
83
96 template <class ContainerT, class FluidState>
97 static void capillaryPressures(ContainerT& values,
98 const Params& params,
99 const FluidState& state)
100 {
101 typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
102
103 for (unsigned phaseIdx = 0; phaseIdx < Traits::numPhases; ++phaseIdx) {
104 const Evaluation& S =
105 decay<Evaluation>(state.saturation(phaseIdx));
106 Valgrind::CheckDefined(S);
107
108 values[phaseIdx] =
109 S*params.pcMaxSat(phaseIdx) +
110 (1.0 - S)*params.pcMinSat(phaseIdx);
111 }
112 }
113
117 template <class ContainerT, class FluidState>
118 static void saturations(ContainerT& /*values*/,
119 const Params& /*params*/,
120 const FluidState& /*state*/)
121 {
122 throw std::runtime_error("Not implemented: LinearMaterial::saturations()");
123 }
124
128 template <class ContainerT, class FluidState>
129 static void relativePermeabilities(ContainerT& values,
130 const Params& /*params*/,
131 const FluidState& state)
132 {
133 typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
134
135 for (unsigned phaseIdx = 0; phaseIdx < Traits::numPhases; ++phaseIdx) {
136 const Evaluation& S =
137 decay<Evaluation>(state.saturation(phaseIdx));
138 Valgrind::CheckDefined(S);
139
140 values[phaseIdx] = max(min(S,1.0),0.0);
141 }
142 }
143
147 template <class FluidState, class Evaluation = typename FluidState::Scalar>
148 static Evaluation pcnw(const Params& params, const FluidState& fs)
149 {
150 const Evaluation& Sw =
151 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
152 Valgrind::CheckDefined(Sw);
153
154 const Evaluation& wPhasePressure =
155 Sw*params.pcMaxSat(Traits::wettingPhaseIdx) +
156 (1.0 - Sw)*params.pcMinSat(Traits::wettingPhaseIdx);
157
158 const Evaluation& Sn =
159 decay<Evaluation>(fs.saturation(Traits::nonWettingPhaseIdx));
160 Valgrind::CheckDefined(Sn);
161
162 const Evaluation& nPhasePressure =
163 Sn*params.pcMaxSat(Traits::nonWettingPhaseIdx) +
164 (1.0 - Sn)*params.pcMinSat(Traits::nonWettingPhaseIdx);
165
166 return nPhasePressure - wPhasePressure;
167 }
168
169
170 template <class Evaluation = Scalar>
171 static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
172 twoPhaseSatPcnw(const Params& params, const Evaluation& Sw)
173 {
174 const Evaluation& wPhasePressure =
175 Sw*params.pcMaxSat(Traits::wettingPhaseIdx) +
176 (1.0 - Sw)*params.pcMinSat(Traits::wettingPhaseIdx);
177
178 const Evaluation& nPhasePressure =
179 (1.0 - Sw)*params.pcMaxSat(Traits::nonWettingPhaseIdx) +
180 Sw*params.pcMinSat(Traits::nonWettingPhaseIdx);
181
182 return nPhasePressure - wPhasePressure;
183 }
184
189 template <class FluidState, class Evaluation = typename FluidState::Scalar>
190 static Evaluation Sw(const Params& /*params*/, const FluidState& /*fs*/)
191 { throw std::runtime_error("Not implemented: Sw()"); }
192
193 template <class Evaluation = Scalar>
194 static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
195 twoPhaseSatSw(const Params& /*params*/, const Evaluation& /*Sw*/)
196 { throw std::runtime_error("Not implemented: twoPhaseSatSw()"); }
197
202 template <class FluidState, class Evaluation = typename FluidState::Scalar>
203 static Evaluation Sn(const Params& /*params*/, const FluidState& /*fs*/)
204 { throw std::runtime_error("Not implemented: Sn()"); }
205
206 template <class Evaluation = Scalar>
207 static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
208 twoPhaseSatSn(const Params& /*params*/, const Evaluation& /*Sw*/)
209 { throw std::runtime_error("Not implemented: twoPhaseSatSn()"); }
210
217 template <class FluidState, class Evaluation = typename FluidState::Scalar>
218 static typename std::enable_if<Traits::numPhases == 3, Evaluation>::type
219 Sg(const Params& /*params*/, const FluidState& /*fs*/)
220 { throw std::runtime_error("Not implemented: Sg()"); }
221
225 template <class FluidState, class Evaluation = typename FluidState::Scalar>
226 static Evaluation krw(const Params& /*params*/, const FluidState& fs)
227 {
228 const Evaluation& sw =
229 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
230 return max(0.0, min(1.0, sw));
231 }
232
233 template <class Evaluation = Scalar>
234 static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
235 twoPhaseSatKrw(const Params& /*params*/, const Evaluation& Sw)
236 { return max(0.0, min(1.0, Sw)); }
237
241 template <class FluidState, class Evaluation = typename FluidState::Scalar>
242 static Evaluation krn(const Params& /*params*/, const FluidState& fs)
243 {
244 const Evaluation& sn =
245 decay<Evaluation>(fs.saturation(Traits::nonWettingPhaseIdx));
246 return max(0.0, min(1.0, sn));
247 }
248
249 template <class Evaluation = Scalar>
250 static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
251 twoPhaseSatKrn(const Params& /*params*/, const Evaluation& Sw)
252 {
253 return max(0.0, min(1.0, Sw));
254 }
255
261 template <class FluidState, class Evaluation = typename FluidState::Scalar>
262 static typename std::enable_if<Traits::numPhases == 3, Evaluation>::type
263 krg(const Params& /*params*/, const FluidState& fs)
264 {
265 const Evaluation& sg =
266 decay<Evaluation>(fs.saturation(Traits::gasPhaseIdx));
267 return max(0.0, min(1.0, sg));
268 }
269
275 template <class FluidState, class Evaluation = Scalar>
276 static typename std::enable_if<Traits::numPhases == 3, Evaluation>::type
277 pcgn(const Params& params, const FluidState& fs)
278 {
279 const Evaluation& sn =
280 decay<Evaluation>(fs.saturation(Traits::nonWettingPhaseIdx));
281 Valgrind::CheckDefined(sn);
282
283 const Evaluation& nPhasePressure =
284 Sn*params.pcMaxSat(Traits::nonWettingPhaseIdx) +
285 (1.0 - Sn)*params.pcMinSat(Traits::nonWettingPhaseIdx);
286
287 const Evaluation& sg =
288 decay<Evaluation>(fs.saturation(Traits::gasPhaseIdx));
289 Valgrind::CheckDefined(sg);
290
291 const Evaluation& gPhasePressure =
292 Sg*params.pcMaxSat(Traits::gasPhaseIdx) +
293 (1.0 - Sg)*params.pcMinSat(Traits::gasPhaseIdx);
294
295 return gPhasePressure - nPhasePressure;
296 }
297};
298} // namespace Opm
299
300#endif
Reference implementation of params for the linear M-phase material material.
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...
Some templates to wrap the valgrind client request macros.
Implements a linear saturation-capillary pressure relation.
Definition LinearMaterial.hpp:51
static Evaluation pcnw(const Params &params, const FluidState &fs)
The difference between the pressures of the non-wetting and wetting phase.
Definition LinearMaterial.hpp:148
static Evaluation krn(const Params &, const FluidState &fs)
The relative permability of the liquid non-wetting phase.
Definition LinearMaterial.hpp:242
static void relativePermeabilities(ContainerT &values, const Params &, const FluidState &state)
The relative permeability of all phases.
Definition LinearMaterial.hpp:129
static const bool implementsTwoPhaseSatApi
Specify whether this material law implements the two-phase convenience API which only depends on the ...
Definition LinearMaterial.hpp:66
static const bool isPressureDependent
Specify whether the quantities defined by this material law are dependent on the absolute pressure.
Definition LinearMaterial.hpp:74
static std::enable_if< Traits::numPhases==3, Evaluation >::type pcgn(const Params &params, const FluidState &fs)
The difference between the pressures of the gas and the non-wetting phase.
Definition LinearMaterial.hpp:277
static const bool isSaturationDependent
Specify whether the quantities defined by this material law are saturation dependent.
Definition LinearMaterial.hpp:70
static std::enable_if< Traits::numPhases==3, Evaluation >::type Sg(const Params &, const FluidState &)
Calculate gas phase saturation given that the rest of the fluid state has been initialized.
Definition LinearMaterial.hpp:219
static Evaluation Sw(const Params &, const FluidState &)
Calculate wetting phase saturation given that the rest of the fluid state has been initialized.
Definition LinearMaterial.hpp:190
static Evaluation Sn(const Params &, const FluidState &)
Calculate non-wetting liquid phase saturation given that the rest of the fluid state has been initial...
Definition LinearMaterial.hpp:203
static const bool isTemperatureDependent
Specify whether the quantities defined by this material law are temperature dependent.
Definition LinearMaterial.hpp:78
static void saturations(ContainerT &, const Params &, const FluidState &)
The inverse of the capillary pressure.
Definition LinearMaterial.hpp:118
static const bool isCompositionDependent
Specify whether the quantities defined by this material law are dependent on the phase composition.
Definition LinearMaterial.hpp:82
static const bool implementsTwoPhaseApi
Specify whether this material law implements the two-phase convenience API.
Definition LinearMaterial.hpp:62
static void capillaryPressures(ContainerT &values, const Params &params, const FluidState &state)
The linear capillary pressure-saturation curve.
Definition LinearMaterial.hpp:97
static std::enable_if< Traits::numPhases==3, Evaluation >::type krg(const Params &, const FluidState &fs)
The relative permability of the gas phase.
Definition LinearMaterial.hpp:263
static Evaluation krw(const Params &, const FluidState &fs)
The relative permability of the wetting phase.
Definition LinearMaterial.hpp:226
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30