version 3.8
1padapter.hh
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_FLUIDSYTEMS_ONEP_ADAPTER_HH
13#define DUMUX_FLUIDSYTEMS_ONEP_ADAPTER_HH
14
15#include <cassert>
16
17#include <dune/common/exceptions.hh>
18
21
22namespace Dumux {
23namespace FluidSystems {
24
31template <class MPFluidSystem, int phase = 0>
33: public Base<typename MPFluidSystem::Scalar, OnePAdapter<MPFluidSystem, phase>>
34{
36
37 static_assert(phase < MPFluidSystem::numPhases, "Phase does not exist in multi-phase fluidsystem!");
38
39 struct AdapterPolicy
40 {
41 using FluidSystem = MPFluidSystem;
42
43 // the phase index is always zero, other phases than the chosen phase should never be called
44 static int phaseIdx(int mpFluidPhaseIdx)
45 { return 0; }
46
47 // the main component is currently excepted to have the same index as it's phase
48 // (see Fluidsystems::Base::getMainComponent for more information)
49 // so we swap the main component with the first component
50 // this mapping works in both ways since we are only swapping components
51 static constexpr int compIdx(int compIdx)
52 {
53 if (compIdx == 0)
54 return phase;
55 else if (compIdx == phase)
56 return 0;
57 else
58 return compIdx;
59 }
60 };
61
62 template<class FluidState>
63 static auto adaptFluidState(const FluidState& fluidState)
65
66public:
67 using Scalar = typename MPFluidSystem::Scalar;
69
71 using MultiPhaseFluidSystem = MPFluidSystem;
73 static constexpr int multiphaseFluidsystemPhaseIdx = phase;
74
76 static constexpr int numPhases = 1;
79 static constexpr int numComponents = MultiPhaseFluidSystem::isMiscible() ? MultiPhaseFluidSystem::numComponents : numPhases;
81 static constexpr int phase0Idx = 0;
82
84 static constexpr int compIdx(int multiPhaseFluidSystemCompIdx)
85 { return AdapterPolicy::compIdx(multiPhaseFluidSystemCompIdx); }
86
90 template<class ...Args>
91 static void init(Args&&... args)
92 { MultiPhaseFluidSystem::init(std::forward<Args>(args)...); }
93
94 /****************************************
95 * Fluid phase related static parameters
96 ****************************************/
102 static std::string phaseName(int phaseIdx = 0)
103 { return MultiPhaseFluidSystem::phaseName(phase); }
104
110 static std::string componentName(int compIdx)
111 { return MultiPhaseFluidSystem::componentName(AdapterPolicy::compIdx(compIdx)); }
112
116 static std::string name()
117 { return MultiPhaseFluidSystem::phaseName(phase); }
118
122 static constexpr bool isMiscible()
123 { return false; }
124
128 static constexpr bool isGas(int phaseIdx = 0)
129 { return MultiPhaseFluidSystem::isGas(phase); }
130
145 static constexpr bool isIdealMixture(int phaseIdx = 0)
146 { return MultiPhaseFluidSystem::isIdealMixture(phase); }
147
151 static constexpr bool isCompressible(int phaseIdx = 0)
152 { return MultiPhaseFluidSystem::isCompressible(phase); }
153
157 static constexpr bool viscosityIsConstant(int phaseIdx = 0)
158 { return MultiPhaseFluidSystem::viscosityIsConstant(phase); }
159
163 static constexpr bool isIdealGas(int phaseIdx = 0)
164 { return MultiPhaseFluidSystem::isIdealGas(phase); }
165
170 { return MultiPhaseFluidSystem::molarMass(AdapterPolicy::compIdx(compIdx)); }
171
174 template <class FluidState>
175 static Scalar density(const FluidState &fluidState, int phaseIdx = 0)
176 {
177 assert(phaseIdx == 0);
178 return MultiPhaseFluidSystem::density(adaptFluidState(fluidState), phase);
179 }
180
183 template <class FluidState>
184 static Scalar molarDensity(const FluidState &fluidState, int phaseIdx = 0)
185 {
186 assert(phaseIdx == 0);
187 return MultiPhaseFluidSystem::molarDensity(adaptFluidState(fluidState), phase);
188 }
189
192 template <class FluidState>
193 static Scalar enthalpy(const FluidState &fluidState, int phaseIdx = 0)
194 {
195 assert(phaseIdx == 0);
196 return MultiPhaseFluidSystem::enthalpy(adaptFluidState(fluidState), phase);
197 }
198
206 template <class FluidState>
207 static Scalar componentEnthalpy(const FluidState &fluidState,
208 int phaseIdx,
209 int compIdx)
210 {
211 assert(phaseIdx == 0);
212 return MultiPhaseFluidSystem::componentEnthalpy(adaptFluidState(fluidState), phase,
213 AdapterPolicy::compIdx(compIdx));
214 }
215
218 template <class FluidState>
219 static Scalar viscosity(const FluidState &fluidState, int phaseIdx = 0)
220 {
221 assert(phaseIdx == 0);
222 return MultiPhaseFluidSystem::viscosity(adaptFluidState(fluidState), phase);
223 }
224
227 template <class FluidState>
228 static Scalar fugacityCoefficient(const FluidState &fluidState,
229 int phaseIdx,
230 int compIdx)
231 {
232 assert(phaseIdx == 0);
233 return MultiPhaseFluidSystem::fugacityCoefficient(adaptFluidState(fluidState), phase,
234 AdapterPolicy::compIdx(compIdx));
235 }
236
239 template <class FluidState>
240 static Scalar diffusionCoefficient(const FluidState &fluidState,
241 int phaseIdx,
242 int compIdx)
243 {
244 assert(phaseIdx == 0);
245 return MultiPhaseFluidSystem::diffusionCoefficient(adaptFluidState(fluidState), phase,
246 AdapterPolicy::compIdx(compIdx));
247 }
248
251 template <class FluidState>
252 static Scalar binaryDiffusionCoefficient(const FluidState &fluidState,
253 int phaseIdx,
254 int compIIdx,
255 int compJIdx)
256 {
257 assert(phaseIdx == 0);
258 return MultiPhaseFluidSystem::binaryDiffusionCoefficient(adaptFluidState(fluidState), phase,
259 AdapterPolicy::compIdx(compIIdx),
260 AdapterPolicy::compIdx(compJIdx));
261 }
262
265 template <class FluidState>
266 static Scalar thermalConductivity(const FluidState &fluidState,
267 int phaseIdx = 0)
268 {
269 assert(phaseIdx == 0);
270 return MultiPhaseFluidSystem::thermalConductivity(adaptFluidState(fluidState), phase);
271 }
272
275 template <class FluidState>
276 static Scalar heatCapacity(const FluidState &fluidState,
277 int phaseIdx = 0)
278 {
279 assert(phaseIdx == 0);
280 return MultiPhaseFluidSystem::heatCapacity(adaptFluidState(fluidState), phase);
281 }
282};
283
284} // namespace FluidSystems
285} // namespace
286
287#endif
Adapter class for fluid states with different indices.
Adapter class for fluid states with different indices.
Definition: adapter.hh:32
Fluid system base class.
Definition: fluidsystems/base.hh:33
An adapter for multi-phase fluid systems to be used with (compositional) one-phase models.
Definition: 1padapter.hh:34
static Scalar componentEnthalpy(const FluidState &fluidState, int phaseIdx, int compIdx)
Returns the specific enthalpy of a component in a specific phase.
Definition: 1padapter.hh:207
static Scalar thermalConductivity(const FluidState &fluidState, int phaseIdx=0)
Thermal conductivity of a fluid phase .
Definition: 1padapter.hh:266
static Scalar viscosity(const FluidState &fluidState, int phaseIdx=0)
Calculate the dynamic viscosity of a fluid phase .
Definition: 1padapter.hh:219
static std::string phaseName(int phaseIdx=0)
Return the human readable name of a fluid phase.
Definition: 1padapter.hh:102
static void init(Args &&... args)
Initialize the fluid system's static parameters generically.
Definition: 1padapter.hh:91
static Scalar density(const FluidState &fluidState, int phaseIdx=0)
Calculate the density of a fluid phase.
Definition: 1padapter.hh:175
static constexpr int multiphaseFluidsystemPhaseIdx
the index of the phase we choose from the multi-phase fluid system
Definition: 1padapter.hh:73
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx=0)
Calculate the molar density of a fluid phase.
Definition: 1padapter.hh:184
static constexpr int numPhases
number of phases in the fluid system
Definition: 1padapter.hh:76
static constexpr int phase0Idx
number of components has to be the same as in the multi-phase fluid system as the composition needs t...
Definition: 1padapter.hh:81
static Scalar binaryDiffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIIdx, int compJIdx)
Given a phase's composition, temperature and pressure, return the binary diffusion coefficient for c...
Definition: 1padapter.hh:252
typename MPFluidSystem::Scalar Scalar
Definition: 1padapter.hh:67
static constexpr int compIdx(int multiPhaseFluidSystemCompIdx)
convert a component index of the multi-phase component index to the actual component index
Definition: 1padapter.hh:84
static Scalar enthalpy(const FluidState &fluidState, int phaseIdx=0)
Given a phase's composition, temperature, pressure and density, calculate its specific enthalpy .
Definition: 1padapter.hh:193
MPFluidSystem MultiPhaseFluidSystem
export the wrapped MultiPhaseFluidSystem type
Definition: 1padapter.hh:71
static Scalar molarMass(int compIdx)
The mass in of one mole of the component.
Definition: 1padapter.hh:169
static std::string componentName(int compIdx)
A human readable name for the component.
Definition: 1padapter.hh:110
static std::string name()
A human readable name for the component.
Definition: 1padapter.hh:116
static constexpr bool viscosityIsConstant(int phaseIdx=0)
Returns true if the fluid viscosity is constant.
Definition: 1padapter.hh:157
static constexpr bool isGas(int phaseIdx=0)
Returns whether the fluid is gaseous.
Definition: 1padapter.hh:128
static constexpr bool isMiscible()
There is only one phase, so no mass transfer between phases can occur.
Definition: 1padapter.hh:122
static Scalar fugacityCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Calculate the fugacity coefficient of an individual component in a fluid phase.
Definition: 1padapter.hh:228
static constexpr bool isIdealMixture(int phaseIdx=0)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: 1padapter.hh:145
static constexpr bool isIdealGas(int phaseIdx=0)
Returns true if the fluid is assumed to be an ideal gas.
Definition: 1padapter.hh:163
static constexpr int numComponents
Definition: 1padapter.hh:79
static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx=0)
Specific isobaric heat capacity of a fluid phase .
Definition: 1padapter.hh:276
static constexpr bool isCompressible(int phaseIdx=0)
Returns true if the fluid is assumed to be compressible.
Definition: 1padapter.hh:151
static Scalar diffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Calculate the binary molecular diffusion coefficient for a component in a fluid phase .
Definition: 1padapter.hh:240
The a parameter cache which does nothing.
Definition: nullparametercache.hh:22
Fluid system base class.
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:62
std::string molarDensity(int phaseIdx) noexcept
I/O name of molar density for multiphase systems.
Definition: name.hh:71
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:53
Definition: adapt.hh:17