SplineTwoPhaseMaterial.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_SPLINE_TWO_PHASE_MATERIAL_HPP
28#define OPM_SPLINE_TWO_PHASE_MATERIAL_HPP
29
31
32#include <algorithm>
33#include <cmath>
34#include <cassert>
35
36namespace Opm {
43template <class TraitsT, class ParamsT = SplineTwoPhaseMaterialParams<TraitsT> >
44class SplineTwoPhaseMaterial : public TraitsT
45{
46 typedef typename ParamsT::SamplePoints SamplePoints;
47
48public:
50 typedef TraitsT Traits;
51
53 typedef ParamsT Params;
54
56 typedef typename Traits::Scalar Scalar;
57
59 static const int numPhases = Traits::numPhases;
60 static_assert(numPhases == 2,
61 "The piecewise linear two-phase capillary pressure law only"
62 "applies to the case of two fluid phases");
63
66 static const bool implementsTwoPhaseApi = true;
67
70 static const bool implementsTwoPhaseSatApi = true;
71
74 static const bool isSaturationDependent = true;
75
78 static const bool isPressureDependent = false;
79
82 static const bool isTemperatureDependent = false;
83
86 static const bool isCompositionDependent = false;
87
91 template <class Container, class FluidState>
92 static void capillaryPressures(Container& values, const Params& params, const FluidState& fluidState)
93 {
94 typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
95
96 values[Traits::wettingPhaseIdx] = 0.0; // reference phase
97 values[Traits::nonWettingPhaseIdx] = pcnw<FluidState, Evaluation>(params, fluidState);
98 }
99
104 template <class Container, class FluidState>
105 static void saturations(Container& /*values*/, const Params& /*params*/, const FluidState& /*fluidState*/)
106 { throw std::logic_error("Not implemented: saturations()"); }
107
111 template <class Container, class FluidState>
112 static void relativePermeabilities(Container& values, const Params& params, const FluidState& fluidState)
113 {
114 typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
115
116 values[Traits::wettingPhaseIdx] = krw<FluidState, Evaluation>(params, fluidState);
117 values[Traits::nonWettingPhaseIdx] = krn<FluidState, Evaluation>(params, fluidState);
118 }
119
123 template <class FluidState, class Evaluation = typename FluidState::Scalar>
124 static Evaluation pcnw(const Params& params, const FluidState& fluidState)
125 {
126 const Evaluation& Sw =
127 decay<Evaluation>(fluidState.saturation(Traits::wettingPhaseIdx));
128
129 return twoPhaseSatPcnw(params, Sw);
130 }
131
135 template <class Evaluation>
136 static Evaluation twoPhaseSatPcnw(const Params& params, const Evaluation& Sw)
137 {
138 // this assumes that the capillary pressure is monotonically decreasing
139 const auto& pcnwSpline = params.pcnwSpline();
140 if (Sw <= pcnwSpline.xAt(0))
141 return Evaluation(pcnwSpline.valueAt(0));
142 if (Sw >= pcnwSpline.xAt(pcnwSpline.numSamples() - 1))
143 return Evaluation(pcnwSpline.valueAt(pcnwSpline.numSamples() - 1));
144
145 return pcnwSpline.eval(Sw);
146 }
147
148 template <class Evaluation>
149 static Evaluation twoPhaseSatPcnwInv(const Params& params, const Evaluation& pcnw)
150 {
151 static const Evaluation nil(0.0);
152
153 // this assumes that the capillary pressure is monotonically decreasing
154 const auto& pcnwSpline = params.pcnwSpline();
155 if (pcnw >= pcnwSpline.valueAt(0))
156 return Evaluation(pcnwSpline.xAt(0));
157 if (pcnw <= pcnwSpline.y(pcnwSpline.numSamples() - 1))
158 return Evaluation(pcnwSpline.xAt(pcnwSpline.numSamples() - 1));
159
160 // the intersect() method of splines is a bit slow, but this code path is not too
161 // time critical...
162 return pcnwSpline.intersect(/*a=*/nil, /*b=*/nil, /*c=*/nil, /*d=*/pcnw);
163 }
164
168 template <class FluidState, class Evaluation = typename FluidState::Scalar>
169 static Evaluation Sw(const Params& /*params*/, const FluidState& /*fluidState*/)
170 { throw std::logic_error("Not implemented: Sw()"); }
171
172 template <class Evaluation>
173 static Evaluation twoPhaseSatSw(const Params& /*params*/, const Evaluation& /*pC*/)
174 { throw std::logic_error("Not implemented: twoPhaseSatSw()"); }
175
180 template <class FluidState, class Evaluation = typename FluidState::Scalar>
181 static Evaluation Sn(const Params& params, const FluidState& fluidState)
182 { return 1 - Sw<FluidState, Evaluation>(params, fluidState); }
183
184 template <class Evaluation>
185 static Evaluation twoPhaseSatSn(const Params& params, const Evaluation& pC)
186 { return 1 - twoPhaseSatSw(params, pC); }
187
192 template <class FluidState, class Evaluation = typename FluidState::Scalar>
193 static Evaluation krw(const Params& params, const FluidState& fluidState)
194 {
195 const Evaluation& Sw =
196 decay<Evaluation>(fluidState.saturation(Traits::wettingPhaseIdx));
197
198 return twoPhaseSatKrw(params, Sw);
199 }
200
201 template <class Evaluation>
202 static Evaluation twoPhaseSatKrw(const Params& params, const Evaluation& Sw)
203 {
204 const auto& krwSpline = params.krwSpline();
205 if (Sw <= krwSpline.xAt(0))
206 return Evaluation(krwSpline.valueAt(0));
207 if (Sw >= krwSpline.xAt(krwSpline.numSamples() - 1))
208 return Evaluation(krwSpline.valueAt(krwSpline.numSamples() - 1));
209
210 return krwSpline.eval(Sw);
211 }
212
213 template <class Evaluation>
214 static Evaluation twoPhaseSatKrwInv(const Params& params, const Evaluation& krw)
215 {
216 static const Evaluation nil(0.0);
217
218 const auto& krwSpline = params.krwSpline();
219 if (krw <= krwSpline.valueAt(0))
220 return Evaluation(krwSpline.xAt(0));
221 if (krw >= krwSpline.valueAt(krwSpline.numSamples() - 1))
222 return Evaluation(krwSpline.xAt(krwSpline.numSamples() - 1));
223
224 return krwSpline.intersect(/*a=*/nil, /*b=*/nil, /*c=*/nil, /*d=*/krw);
225 }
226
231 template <class FluidState, class Evaluation = typename FluidState::Scalar>
232 static Evaluation krn(const Params& params, const FluidState& fluidState)
233 {
234 const Evaluation& Sn =
235 decay<Evaluation>(fluidState.saturation(Traits::nonWettingPhaseIdx));
236
237 return twoPhaseSatKrn(params, 1.0 - Sn);
238 }
239
240 template <class Evaluation>
241 static Evaluation twoPhaseSatKrn(const Params& params, const Evaluation& Sw)
242 {
243 const auto& krnSpline = params.krnSpline();
244 if (Sw <= krnSpline.xAt(0))
245 return Evaluation(krnSpline.valueAt(0));
246 if (Sw >= krnSpline.xAt(krnSpline.numSamples() - 1))
247 return Evaluation(krnSpline.valueAt(krnSpline.numSamples() - 1));
248
249 return krnSpline.eval(Sw);
250 }
251
252 template <class Evaluation>
253 static Evaluation twoPhaseSatKrnInv(const Params& params, const Evaluation& krn)
254 {
255 static const Evaluation nil(0.0);
256
257 const auto& krnSpline = params.krnSpline();
258 if (krn >= krnSpline.valueAt(0))
259 return Evaluation(krnSpline.xAt(0));
260 if (krn <= krnSpline.valueAt(krnSpline.numSamples() - 1))
261 return Evaluation(krnSpline.xAt(krnSpline.numSamples() - 1));
262
263 return krnSpline.intersect(/*a=*/nil, /*b=*/nil, /*c=*/nil, /*d=*/krn);
264 }
265};
266} // namespace Opm
267
268#endif
Implementation of a tabulated capillary pressure and relperm law which uses spline curves as interpol...
Definition: SplineTwoPhaseMaterial.hpp:45
static const bool implementsTwoPhaseApi
Definition: SplineTwoPhaseMaterial.hpp:66
static Evaluation krw(const Params &params, const FluidState &fluidState)
The relative permeability for the wetting phase of the porous medium.
Definition: SplineTwoPhaseMaterial.hpp:193
static const bool isTemperatureDependent
Definition: SplineTwoPhaseMaterial.hpp:82
Traits::Scalar Scalar
The type of the scalar values for this law.
Definition: SplineTwoPhaseMaterial.hpp:56
static void capillaryPressures(Container &values, const Params &params, const FluidState &fluidState)
The capillary pressure-saturation curve.
Definition: SplineTwoPhaseMaterial.hpp:92
static Evaluation twoPhaseSatKrw(const Params &params, const Evaluation &Sw)
Definition: SplineTwoPhaseMaterial.hpp:202
static Evaluation Sw(const Params &, const FluidState &)
The saturation-capillary pressure curve.
Definition: SplineTwoPhaseMaterial.hpp:169
static Evaluation twoPhaseSatKrnInv(const Params &params, const Evaluation &krn)
Definition: SplineTwoPhaseMaterial.hpp:253
static void relativePermeabilities(Container &values, const Params &params, const FluidState &fluidState)
The relative permeabilities.
Definition: SplineTwoPhaseMaterial.hpp:112
static void saturations(Container &, const Params &, const FluidState &)
The saturations of the fluid phases starting from their pressure differences.
Definition: SplineTwoPhaseMaterial.hpp:105
static Evaluation krn(const Params &params, const FluidState &fluidState)
The relative permeability for the non-wetting phase of the porous medium.
Definition: SplineTwoPhaseMaterial.hpp:232
static Evaluation twoPhaseSatPcnwInv(const Params &params, const Evaluation &pcnw)
Definition: SplineTwoPhaseMaterial.hpp:149
static Evaluation twoPhaseSatKrwInv(const Params &params, const Evaluation &krw)
Definition: SplineTwoPhaseMaterial.hpp:214
static const bool implementsTwoPhaseSatApi
Definition: SplineTwoPhaseMaterial.hpp:70
static const bool isPressureDependent
Definition: SplineTwoPhaseMaterial.hpp:78
static Evaluation twoPhaseSatKrn(const Params &params, const Evaluation &Sw)
Definition: SplineTwoPhaseMaterial.hpp:241
static Evaluation pcnw(const Params &params, const FluidState &fluidState)
The capillary pressure-saturation curve.
Definition: SplineTwoPhaseMaterial.hpp:124
ParamsT Params
The type of the parameter objects for this law.
Definition: SplineTwoPhaseMaterial.hpp:53
static Evaluation Sn(const Params &params, const FluidState &fluidState)
Calculate the non-wetting phase saturations depending on the phase pressures.
Definition: SplineTwoPhaseMaterial.hpp:181
static const bool isCompositionDependent
Definition: SplineTwoPhaseMaterial.hpp:86
static const bool isSaturationDependent
Definition: SplineTwoPhaseMaterial.hpp:74
TraitsT Traits
The traits class for this material law.
Definition: SplineTwoPhaseMaterial.hpp:50
static Evaluation twoPhaseSatSw(const Params &, const Evaluation &)
Definition: SplineTwoPhaseMaterial.hpp:173
static const int numPhases
The number of fluid phases.
Definition: SplineTwoPhaseMaterial.hpp:59
static Evaluation twoPhaseSatPcnw(const Params &params, const Evaluation &Sw)
The saturation-capillary pressure curve.
Definition: SplineTwoPhaseMaterial.hpp:136
static Evaluation twoPhaseSatSn(const Params &params, const Evaluation &pC)
Definition: SplineTwoPhaseMaterial.hpp:185
Definition: Air_Mesitylene.hpp:34