opm-common
Loading...
Searching...
No Matches
TwoPhaseLETCurvesParams.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 3 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_TWO_PHASE_LET_CURVES_PARAMS_HPP
28
#define OPM_TWO_PHASE_LET_CURVES_PARAMS_HPP
29
30
#include <
opm/material/common/Valgrind.hpp
>
31
#include <
opm/material/common/EnsureFinalized.hpp
>
32
33
namespace
Opm
{
34
43
template
<
class
TraitsT>
44
class
TwoPhaseLETCurvesParams :
public
EnsureFinalized
45
{
46
using
Scalar =
typename
TraitsT::Scalar;
47
48
public
:
49
using
Traits = TraitsT;
50
51
static
constexpr
int
wIdx = 0;
//wetting phase index for two phase let
52
static
constexpr
int
nwIdx = 1;
//non-wetting phase index for two phase let
53
54
TwoPhaseLETCurvesParams()
55
{
56
Valgrind::SetUndefined
(*
this
);
57
}
58
59
virtual
~TwoPhaseLETCurvesParams() {}
60
65
void
finalize
()
66
{
67
EnsureFinalized :: finalize ();
68
69
// printLETCoeffs();
70
}
71
75
Scalar
Smin
(
const
unsigned
phaseIdx)
const
76
{ EnsureFinalized::check();
if
(phaseIdx<Traits::numPhases)
return
Smin_[phaseIdx];
return
0.0;}
77
81
Scalar
dS
(
const
unsigned
phaseIdx)
const
82
{ EnsureFinalized::check();
if
(phaseIdx<Traits::numPhases)
return
dS_[phaseIdx];
return
0.0;}
83
87
Scalar
Sminpc
()
const
88
{ EnsureFinalized::check();
return
Sminpc_; }
89
93
Scalar
dSpc
()
const
94
{ EnsureFinalized::check();
return
dSpc_; }
95
99
Scalar
L
(
const
unsigned
phaseIdx)
const
100
{ EnsureFinalized::check();
if
(phaseIdx<Traits::numPhases)
return
L_[phaseIdx];
return
0.0;}
101
105
Scalar
E
(
const
unsigned
phaseIdx)
const
106
{ EnsureFinalized::check();
if
(phaseIdx<Traits::numPhases)
return
E_[phaseIdx];
return
0.0;}
107
111
Scalar
T
(
const
unsigned
phaseIdx)
const
112
{ EnsureFinalized::check();
if
(phaseIdx<Traits::numPhases)
return
T_[phaseIdx];
return
0.0;}
113
117
Scalar
Krt
(
const
unsigned
phaseIdx)
const
118
{ EnsureFinalized::check();
if
(phaseIdx<Traits::numPhases)
return
Krt_[phaseIdx];
return
0.0;}
119
123
Scalar
Lpc
()
const
124
{ EnsureFinalized::check();
return
Lpc_; }
125
129
Scalar
Epc
()
const
130
{ EnsureFinalized::check();
return
Epc_; }
131
135
Scalar
Tpc
()
const
136
{ EnsureFinalized::check();
return
Tpc_; }
137
141
Scalar
Pcir
()
const
142
{ EnsureFinalized::check();
return
Pcir_; }
143
147
Scalar
Pct
()
const
148
{ EnsureFinalized::check();
return
Pct_; }
149
156
template
<
class
Container>
157
void
setKrwSamples
(
const
Container& letProp,
const
Container& )
158
{
159
setLETCoeffs(wIdx, letProp[2], letProp[3], letProp[4], letProp[5]);
160
//Smin_[wIdx] = letProp[1];
161
//dS_[nwIdx] = 1.0 - letProp[0];
162
Smin_[wIdx] = letProp[0];
163
dS_[wIdx] = letProp[1] - letProp[0];
164
}
165
172
template
<
class
Container>
173
void
setKrnSamples
(
const
Container& letProp,
const
Container& )
174
{
175
setLETCoeffs(nwIdx, letProp[2], letProp[3], letProp[4], letProp[5]);
176
//Smin_[nwIdx] = letProp[1];
177
//dS_[wIdx] = 1.0 - letProp[0];
178
Smin_[nwIdx] = letProp[0];
179
dS_[nwIdx] = letProp[1] - letProp[0];
180
}
181
182
189
template
<
class
Container>
190
void
setPcnwSamples
(
const
Container& letProp,
const
Container& )
191
{
192
setLETPcCoeffs(letProp[2], letProp[3], letProp[4], letProp[5], letProp[6]);