opm-common
Loading...
Searching...
No Matches
SatCurveMultiplexerParams.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_SAT_CURVE_MULTIPLEXER_PARAMS_HPP
28#define OPM_SAT_CURVE_MULTIPLEXER_PARAMS_HPP
29
30
31#include "TwoPhaseLETCurves.hpp"
35
37
38#include <cassert>
39#include <memory>
40#include <stdexcept>
41#include <type_traits>
42
43namespace Opm {
44
45enum class SatCurveMultiplexerApproach {
46 PiecewiseLinear,
47 LET
48};
49
50template <SatCurveMultiplexerApproach ApproachArg>
52{
53 static constexpr auto approach = ApproachArg;
54};
55
56// Helper struct to tell if a parameter pack starts with SatCurveMultiplexerDispatch.
57template <typename ...Args>
58struct FrontIsSatCurveMultiplexerDispatch : public std::false_type {};
59template <SatCurveMultiplexerApproach Value, typename ...Args>
60struct FrontIsSatCurveMultiplexerDispatch<SatCurveMultiplexerDispatch<Value>, Args...> : public std::true_type {};
61template <typename ...Args>
62constexpr bool FrontIsSatCurveMultiplexerDispatchV = FrontIsSatCurveMultiplexerDispatch<Args...>::value;
63
64
73template <class TraitsT>
75{
76public:
77 using Traits = TraitsT;
78 using Scalar = typename TraitsT::Scalar;
79 enum { numPhases = 2 };
80
81private:
82 using LETTwoPhaseLaw = TwoPhaseLETCurves<Traits>;
83 using PLTwoPhaseLaw = PiecewiseLinearTwoPhaseMaterial<Traits>;
84
85 using LETParams = typename LETTwoPhaseLaw::Params;
86 using PLParams = typename PLTwoPhaseLaw::Params;
87
88 template <class ParamT>
89 struct Deleter
90 {
91 inline void operator () ( void* ptr )
92 {
93 delete static_cast< ParamT* > (ptr);
94 }
95 };
96
97 using ParamPointerType = std::shared_ptr<void>;
98
99public:
100
105 {
106 }
107
109 : realParams_()
110 {
111 setApproach( other.approach() );
112 }
113
114 SatCurveMultiplexerParams& operator= ( const SatCurveMultiplexerParams& other )
115 {
116 realParams_.reset();
117 setApproach( other.approach() );
118 return *this;
119 }
120
121 void setApproach(SatCurveMultiplexerApproach newApproach)
122 {
123 assert(realParams_ == 0);
124 approach_ = newApproach;
125
126 switch (approach()) {
127 case SatCurveMultiplexerApproach::LET:
128 realParams_ = ParamPointerType(new LETParams, Deleter< LETParams > () );
129 break;
130
131 case SatCurveMultiplexerApproach::PiecewiseLinear:
132 realParams_ = ParamPointerType(new PLParams, Deleter< PLParams > () );
133 break;
134 }
135 }
136
137 SatCurveMultiplexerApproach approach() const
138 { return approach_; }
139
140 // get the parameter object for the LET curve
141 template <SatCurveMultiplexerApproach approachV>
142 typename std::enable_if<approachV == SatCurveMultiplexerApproach::LET, LETParams>::type&
143 getRealParams()
144 {
145 assert(approach() == approachV);
146 return this->template castTo<LETParams>();
147 }
148
149 template <SatCurveMultiplexerApproach approachV>
150 typename std::enable_if<approachV == SatCurveMultiplexerApproach::LET, const LETParams>::type&
151 getRealParams() const
152 {
153 assert(approach() == approachV);
154 return this->template castTo<LETParams>();
155 }
156
157 // get the parameter object for the PL curve
158 template <SatCurveMultiplexerApproach approachV>
159 typename std::enable_if<approachV == SatCurveMultiplexerApproach::PiecewiseLinear, PLParams>::type&
160 getRealParams()
161 {
162 assert(approach() == approachV);
163 return this->template castTo<PLParams>();
164 }
165
166 template <SatCurveMultiplexerApproach approachV>
167 typename std::enable_if<approachV == SatCurveMultiplexerApproach::PiecewiseLinear, const PLParams>::type&
168 getRealParams() const
169 {
170 assert(approach() == approachV);
171 return this->template castTo<PLParams>();
172 }
173
174 template<class Serializer>
175 void serializeOp(Serializer& serializer)
176 {
177 switch (approach()) {
178 case SatCurveMultiplexerApproach::LET:
179 serializer(castTo<LETParams>());
180 break;
181
182 case SatCurveMultiplexerApproach::PiecewiseLinear:
183 serializer(castTo<PLParams>());
184 break;
185 }
186 }
187
188 Scalar SnTrapped([[maybe_unused]] bool maximumTrapping) const
189 {
190 throw std::logic_error("SatCurveMultiplexerParams::SnTrapped() not implemented");
191 }
192
193 Scalar SnStranded([[maybe_unused]] Scalar sg, [[maybe_unused]] Scalar krg) const
194 {
195 throw std::logic_error("SatCurveMultiplexerParams::SnStranded() not implemented");
196 }
197
198 Scalar SwTrapped() const
199 {
200 throw std::logic_error("SatCurveMultiplexerParams::SwTrapped() not implemented");
201 }
202
203
204private:
205 template <class ParamT>
206 ParamT& castTo()
207 {
208 return *(static_cast<ParamT *> (realParams_.operator->()));
209 }
210
211 template <class ParamT>
212 const ParamT& castTo() const
213 {
214 return *(static_cast<const ParamT *> (realParams_.operator->()));
215 }
216
217 SatCurveMultiplexerApproach approach_{SatCurveMultiplexerApproach::PiecewiseLinear};
218 ParamPointerType realParams_;
219};
220
221} // namespace Opm
222
223#endif // OPM_SAT_CURVE_MULTIPLEXER_PARAMS_HPP
Default implementation for asserting finalization of parameter objects.
Specification of the material parameters for a two-phase material law which uses a table and piecewis...
Implementation of a tabulated, piecewise linear capillary pressure law.
OPM_HOST_DEVICE EnsureFinalized()
The default constructor.
Definition EnsureFinalized.hpp:58
Implementation of a tabulated, piecewise linear capillary pressure law.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:53
PiecewiseLinearTwoPhaseMaterialParams< Traits > Params
Definition PiecewiseLinearTwoPhaseMaterial.hpp:61
Specification of the material parameters for the saturation function multiplexer.
Definition SatCurveMultiplexerParams.hpp:75
SatCurveMultiplexerParams()
The multiplexer constructor.
Definition SatCurveMultiplexerParams.hpp:104
Implementation of the LET curve saturation functions.
Definition TwoPhaseLETCurves.hpp:50
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
Definition SatCurveMultiplexerParams.hpp:58
Definition SatCurveMultiplexerParams.hpp:52