opm-common
Loading...
Searching...
No Matches
BlackOilFunctions.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_MATERIAL_FLUIDSYSTEMS_BLACKOILFUNCTIONS_HEADER_INCLUDED
28#define OPM_MATERIAL_FLUIDSYSTEMS_BLACKOILFUNCTIONS_HEADER_INCLUDED
29
30#include <opm/common/TimingMacros.hpp>
31
34
39
40#include <opm/common/utility/gpuDecorators.hpp>
41
42#include <array>
43#include <cstddef>
44#include <memory>
45#include <stdexcept>
46#include <string>
47#include <string_view>
48#include <vector>
49
50namespace Opm::BlackOil
51{
52OPM_GENERATE_HAS_MEMBER(Rs, ) // Creates 'HasMember_Rs<T>'.
53OPM_GENERATE_HAS_MEMBER(Rv, ) // Creates 'HasMember_Rv<T>'.
54OPM_GENERATE_HAS_MEMBER(Rvw, ) // Creates 'HasMember_Rvw<T>'.
55OPM_GENERATE_HAS_MEMBER(Rsw, ) // Creates 'HasMember_Rsw<T>'.
56OPM_GENERATE_HAS_MEMBER(saltConcentration, )
57OPM_GENERATE_HAS_MEMBER(saltSaturation, )
58OPM_GENERATE_HAS_MEMBER(solventSaturation, )
59OPM_GENERATE_HAS_MEMBER(solventDensity, )
60OPM_GENERATE_HAS_MEMBER(solventInvB, )
62
63template <class FluidSystem, class FluidState, class LhsEval>
64OPM_HOST_DEVICE LhsEval
65getRs_(typename std::enable_if<!HasMember_Rs<FluidState>::value, const FluidState&>::type fluidState,
66 unsigned regionIdx)
67{
68 const auto& XoG = decay<LhsEval>(fluidState.massFraction(FluidSystem::oilPhaseIdx, FluidSystem::gasCompIdx));
69 return FluidSystem::convertXoGToRs(XoG, regionIdx);
70}
71
72template <class FluidSystem, class FluidState, class LhsEval>
73OPM_HOST_DEVICE auto
74getRs_(typename std::enable_if<HasMember_Rs<FluidState>::value, const FluidState&>::type fluidState, unsigned)
75 -> decltype(decay<LhsEval>(fluidState.Rs()))
76{
77 return decay<LhsEval>(fluidState.Rs());
78}
79
80template <class FluidSystem, class FluidState, class LhsEval>
81OPM_HOST_DEVICE LhsEval
82getRv_(typename std::enable_if<!HasMember_Rv<FluidState>::value, const FluidState&>::type fluidState,
83 unsigned regionIdx)
84{
85 const auto& XgO = decay<LhsEval>(fluidState.massFraction(FluidSystem::gasPhaseIdx, FluidSystem::oilCompIdx));
86 return FluidSystem::convertXgOToRv(XgO, regionIdx);
87}
88
89template <class FluidSystem, class FluidState, class LhsEval>
90OPM_HOST_DEVICE auto
91getRv_(typename std::enable_if<HasMember_Rv<FluidState>::value, const FluidState&>::type fluidState, unsigned)
92 -> decltype(decay<LhsEval>(fluidState.Rv()))
93{
94 return decay<LhsEval>(fluidState.Rv());
95}
96
97template <class FluidSystem, class FluidState, class LhsEval>
98OPM_HOST_DEVICE LhsEval
99getRvw_(typename std::enable_if<!HasMember_Rvw<FluidState>::value, const FluidState&>::type fluidState,
100 unsigned regionIdx)
101{
102 const auto& XgW = decay<LhsEval>(fluidState.massFraction(FluidSystem::gasPhaseIdx, FluidSystem::waterCompIdx));
103 return FluidSystem::convertXgWToRvw(XgW, regionIdx);
104}
105
106template <class FluidSystem, class FluidState, class LhsEval>
107OPM_HOST_DEVICE auto
108getRvw_(typename std::enable_if<HasMember_Rvw<FluidState>::value, const FluidState&>::type fluidState, unsigned)
109 -> decltype(decay<LhsEval>(fluidState.Rvw()))
110{
111 return decay<LhsEval>(fluidState.Rvw());
112}
113
114template <class FluidSystem, class FluidState, class LhsEval>
115OPM_HOST_DEVICE LhsEval
116getRsw_(typename std::enable_if<!HasMember_Rsw<FluidState>::value, const FluidState&>::type fluidState,
117 unsigned regionIdx)
118{
119 const auto& XwG = decay<LhsEval>(fluidState.massFraction(FluidSystem::waterPhaseIdx, FluidSystem::gasCompIdx));
120 return FluidSystem::convertXwGToRsw(XwG, regionIdx);
121}
122
123template <class FluidSystem, class FluidState, class LhsEval>
124OPM_HOST_DEVICE auto
125getRsw_(typename std::enable_if<HasMember_Rsw<FluidState>::value, const FluidState&>::type fluidState, unsigned)
126 -> decltype(decay<LhsEval>(fluidState.Rsw()))
127{
128 return decay<LhsEval>(fluidState.Rsw());
129}
130
131template <class FluidState, class LhsEval>
132OPM_HOST_DEVICE LhsEval
133getSaltConcentration_(typename std::enable_if<!HasMember_saltConcentration<FluidState>::value, const FluidState&>::type,
134 unsigned)
135{
136 return 0.0;
137}
138
139template <class FluidState, class LhsEval>
140OPM_HOST_DEVICE auto
141getSaltConcentration_(
142 typename std::enable_if<HasMember_saltConcentration<FluidState>::value, const FluidState&>::type fluidState,
143 unsigned) -> decltype(decay<LhsEval>(fluidState.saltConcentration()))
144{
145 return decay<LhsEval>(fluidState.saltConcentration());
146}
147
148template <class FluidSystem, class FluidState, class LhsEval>
149OPM_HOST_DEVICE LhsEval
150getSaltSaturation_(typename std::enable_if<!HasMember_saltSaturation<FluidState>::value, const FluidState&>::type,
151 unsigned)
152{
153 return 0.0;
154}
155
156template <class FluidSystem, class FluidState, class LhsEval>
157OPM_HOST_DEVICE auto
158getSaltSaturation_(
159 typename std::enable_if<HasMember_saltSaturation<FluidState>::value, const FluidState&>::type fluidState, unsigned)
160 -> decltype(decay<LhsEval>(fluidState.saltSaturation()))
161{
162 return decay<LhsEval>(fluidState.saltSaturation());
163}
164
165template <class FluidState, class LhsEval>
166OPM_HOST_DEVICE LhsEval
167getSolventSaturation_(typename std::enable_if<!HasMember_solventSaturation<FluidState>::value, const FluidState&>::type,
168 unsigned)
169{
170 return 0.0;
171}
172
173template <class FluidState, class LhsEval>
174OPM_HOST_DEVICE auto
175getSolventSaturation_(
176 typename std::enable_if<HasMember_solventSaturation<FluidState>::value, const FluidState&>::type fluidState,
177 unsigned) -> decltype(decay<LhsEval>(fluidState.solventSaturation()))
178{
179 return decay<LhsEval>(fluidState.solventSaturation());
180}
181
182template <class FluidState, class LhsEval>
183OPM_HOST_DEVICE LhsEval
184getSolventDensity_(typename std::enable_if<!HasMember_solventDensity<FluidState>::value, const FluidState&>::type,
185 unsigned)
186{
187 return 0.0;
188}
189
190template <class FluidState, class LhsEval>
191OPM_HOST_DEVICE auto
192getSolventDensity_(
193 typename std::enable_if<HasMember_solventDensity<FluidState>::value, const FluidState&>::type fluidState,
194 unsigned) -> decltype(decay<LhsEval>(fluidState.solventDensity()))
195{
196 return decay<LhsEval>(fluidState.solventDensity());
197}
198
199template <class FluidState, class LhsEval>
200OPM_HOST_DEVICE LhsEval
201getSolventInvB_(typename std::enable_if<!HasMember_solventInvB<FluidState>::value, const FluidState&>::type,
202 unsigned)
203{
204 return 0.0;
205}
206
207template <class FluidState, class LhsEval>
208OPM_HOST_DEVICE auto
209getSolventInvB_(
210 typename std::enable_if<HasMember_solventInvB<FluidState>::value, const FluidState&>::type fluidState,
211 unsigned) -> decltype(decay<LhsEval>(fluidState.solventInvB()))
212{
213 return decay<LhsEval>(fluidState.solventInvB());
214}
215
216template <class FluidState, class LhsEval>
217OPM_HOST_DEVICE LhsEval
218getRsSolw_(typename std::enable_if<!HasMember_rsSolw<FluidState>::value, const FluidState&>::type,
219 unsigned)
220{
221 return 0.0;
222}
223
224template <class FluidState, class LhsEval>
225OPM_HOST_DEVICE auto
226getRsSolw_(
227 typename std::enable_if<HasMember_rsSolw<FluidState>::value, const FluidState&>::type fluidState,
228 unsigned) -> decltype(decay<LhsEval>(fluidState.rsSolw()))
229{
230 return decay<LhsEval>(fluidState.rsSolw());
231}
232
233} // namespace Opm::BlackOil
234#endif
The base class for all fluid systems.
This macro generates a class HasMember_${MEMBER_NAME} which can be used for template specialization.
#define OPM_GENERATE_HAS_MEMBER(MEMBER_NAME,...)
This macro generates a class HasMember_${MEMBER_NAME} which can be used for template specialization.
Definition HasMemberGeneratorMacros.hpp:49
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...
A parameter cache which does nothing.
Some templates to wrap the valgrind client request macros.