OGS
ProcessLib::TH2M::ConstitutiveRelations::PhaseTransition Struct Reference

Detailed Description

Definition at line 21 of file PhaseTransition.h.

#include <PhaseTransition.h>

Inheritance diagram for ProcessLib::TH2M::ConstitutiveRelations::PhaseTransition:
[legend]
Collaboration diagram for ProcessLib::TH2M::ConstitutiveRelations::PhaseTransition:
[legend]

Public Member Functions

 PhaseTransition (std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
 
void eval (SpaceTimeData const &x_t, MediaData const &media_data, GasPressureData const &p_GR, CapillaryPressureData const &p_cap, TemperatureData const &T_data, PureLiquidDensityData const &rho_W_LR, EnthalpyData &enthalpy_data, MassMoleFractionsData &mass_mole_fractions_data, FluidDensityData &fluid_density_data, VapourPartialPressureData &vapour_pressure_data, ConstituentDensityData &constituent_density_data, PhaseTransitionData &cv) const override
 
- Public Member Functions inherited from ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionModel
 PhaseTransitionModel (std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
 
virtual ~PhaseTransitionModel ()=default
 

Private Attributes

int const n_components_gas_
 
int const gas_phase_vapour_component_index_
 
int const gas_phase_dry_air_component_index_
 
int const liquid_phase_solute_component_index_
 
int const liquid_phase_solvent_component_index_
 

Constructor & Destructor Documentation

◆ PhaseTransition()

ProcessLib::TH2M::ConstitutiveRelations::PhaseTransition::PhaseTransition ( std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const & media)
explicit

Definition at line 61 of file PhaseTransition.cpp.

63 : PhaseTransitionModel(media),
65 // Identifies and initialises the indices of the various components via
66 // their properties
69 // Dry air component is complement of vapour component index
71 // The solute is the component of the liquid phase which has been assigned
72 // the property `henry_coefficient'.
74 media, "AqueousLiquid",
76 // The solvent is again the bitwise complement of the solute component
79{
80 DBUG("Create PhaseTransition constitutive model.");
81
82 if (n_components_gas_ != 2)
83 {
85 "The current implementation of PhaseTransitionModelEvaporation "
86 "requires the specification of exactly two components in the "
87 "gas "
88 "phase.");
89 }
90
91 // It is always the first (begin) medium that holds fluid phases.
92 auto const medium = media.begin()->second;
93 auto const& gas_phase = medium->phase("Gas");
94 auto const& liquid_phase = medium->phase("AqueousLiquid");
95
96 // check for minimum requirement definitions in media object
97 std::array const required_vapour_component_properties = {
102
103 std::array const required_dry_air_component_properties = {
106
107 std::array const required_solute_component_properties = {
111
112 std::array const required_solvent_component_properties = {
114
115 std::array const required_gas_properties = {MaterialPropertyLib::density};
116
117 std::array const required_liquid_properties = {
119
121 gas_phase.component(gas_phase_vapour_component_index_),
122 required_vapour_component_properties);
124 gas_phase.component(gas_phase_dry_air_component_index_),
125 required_dry_air_component_properties);
127 liquid_phase.component(liquid_phase_solute_component_index_),
128 required_solute_component_properties);
130 liquid_phase.component(liquid_phase_solvent_component_index_),
131 required_solvent_component_properties);
132 checkRequiredProperties(gas_phase, required_gas_properties);
133 checkRequiredProperties(liquid_phase, required_liquid_properties);
134}
#define OGS_FATAL(...)
Definition Error.h:26
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
void checkRequiredProperties(Component const &c, std::span< PropertyType const > const required_properties)
Definition Component.cpp:60
int findComponentIndex(std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media, std::string const &phase_name, MaterialPropertyLib::PropertyType property_type)
int numberOfComponents(std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media, std::string const &phase_name)
PhaseTransitionModel(std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)

References DBUG(), MaterialPropertyLib::density, MaterialPropertyLib::diffusion, gas_phase_dry_air_component_index_, gas_phase_vapour_component_index_, MaterialPropertyLib::henry_coefficient, liquid_phase_solute_component_index_, liquid_phase_solvent_component_index_, MaterialPropertyLib::molar_mass, n_components_gas_, OGS_FATAL, MaterialPropertyLib::specific_heat_capacity, MaterialPropertyLib::specific_latent_heat, and MaterialPropertyLib::vapour_pressure.

Member Function Documentation

◆ eval()

void ProcessLib::TH2M::ConstitutiveRelations::PhaseTransition::eval ( SpaceTimeData const & x_t,
MediaData const & media_data,
GasPressureData const & p_GR,
CapillaryPressureData const & p_cap,
TemperatureData const & T_data,
PureLiquidDensityData const & rho_W_LR,
EnthalpyData & enthalpy_data,
MassMoleFractionsData & mass_mole_fractions_data,
FluidDensityData & fluid_density_data,
VapourPartialPressureData & vapour_pressure_data,
ConstituentDensityData & constituent_density_data,
PhaseTransitionData & cv ) const
overridevirtual

Implements ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionModel.

Definition at line 136 of file PhaseTransition.cpp.

148{
150
151 // primary variables
152 auto const pGR = p_GR();
153 auto const pCap = p_cap();
154 auto const T = T_data.T;
155 variables.gas_phase_pressure = pGR;
156 variables.capillary_pressure = pCap;
157 variables.temperature = T;
158
159 auto const& liquid_phase = media_data.liquid;
160 auto const& gas_phase = media_data.gas;
161
163
164 // identify vapour and air components for convenient access
165 auto const& vapour_component =
166 gas_phase.component(gas_phase_vapour_component_index_);
167 auto const& dry_air_component =
168 gas_phase.component(gas_phase_dry_air_component_index_);
169
170 // specific latent heat (of evaporation)
171 const auto dh_evap =
172 vapour_component
174 .template value<double>(variables, x_t.x, x_t.t, x_t.dt);
175
176 // molar mass of evaporating component
177 auto const M_W =
178 vapour_component.property(MaterialPropertyLib::PropertyType::molar_mass)
179 .template value<double>(variables, x_t.x, x_t.t, x_t.dt);
180
181 // provide evaporation enthalpy and molar mass of the evaporating
182 // component in the variable array
183 variables.enthalpy_of_evaporation = dh_evap;
184 variables.molar_mass = M_W;
185
186 // vapour pressure over flat interface
187 const auto p_vap_flat =
188 vapour_component
190 .template value<double>(variables, x_t.x, x_t.t, x_t.dt);
191
192 const auto dp_vap_flat_dT =
193 vapour_component
195 .template dValue<double>(variables,
197 x_t.x, x_t.t, x_t.dt);
198
199 // molar mass of dry air component
200 auto const M_C =
201 dry_air_component
203 .template value<double>(variables, x_t.x, x_t.t, x_t.dt);
204
205 // Water pressure is passed to the VariableArray in order to calculate
206 // the water density.
207 auto const pLR = pGR - pCap;
208 variables.liquid_phase_pressure = pLR;
209
210 // Kelvin-Laplace correction for menisci
211 const double K =
212 pCap > 0. ? std::exp(-pCap * M_W / rho_W_LR() / R / T) : 1.;
213 const double dK_dT =
214 pCap > 0. ? pCap * M_W / rho_W_LR() / R / T / T * K : 0;
215 const double dK_dpCap =
216 pCap > 0. ? -M_W / rho_W_LR() / R / T * K
217 : 0.; // rho_W_LR is treated as a constant here. However, the
218 // resulting errors are very small and can be ignored.
219
220 // vapour pressure inside porespace (== water partial pressure in gas
221 // phase)
222 vapour_pressure_data.pWGR = p_vap_flat * K;
223 auto const dpWGR_dT = dp_vap_flat_dT * K + p_vap_flat * dK_dT;
224 auto const dpWGR_dpCap = p_vap_flat * dK_dpCap;
225
226 // gas phase molar fractions
227 auto const xnWG_min =
228 1.e-12; // Magic number; prevents the mass fraction of a component
229 // from ever becoming zero (which would cause the partial
230 // density to disappear and thus one of the Laplace terms
231 // of the mass balance on the diagonal of the local element
232 // matrix to be zero). The value is simply made up, seems
233 // reasonable.
234 double const xnWG =
235 std::clamp(vapour_pressure_data.pWGR / pGR, xnWG_min, 1. - xnWG_min);
236 mass_mole_fractions_data.xnCG = 1. - xnWG;
237
238 // gas phase molar fraction derivatives
239 auto const dxnWG_dpGR = -vapour_pressure_data.pWGR / pGR / pGR;
240 auto const dxnWG_dpCap = dpWGR_dpCap / pGR;
241 auto const dxnWG_dT = dpWGR_dT / pGR;
242
243 // molar mass of the gas phase as a mixture of 'air' and vapour
244 auto const MG = mass_mole_fractions_data.xnCG * M_C + xnWG * M_W;
245 variables.molar_mass = MG;
246
247 // gas phase mass fractions
248 double const xmWG = xnWG * M_W / MG;
249 mass_mole_fractions_data.xmCG = 1. - xmWG;
250
251 auto const dxn_dxm_conversion = M_W * M_C / MG / MG;
252 // gas phase mass fraction derivatives
253 cv.dxmWG_dpGR = dxnWG_dpGR * dxn_dxm_conversion;
254 cv.dxmWG_dpCap = dxnWG_dpCap * dxn_dxm_conversion;
255 cv.dxmWG_dT = dxnWG_dT * dxn_dxm_conversion;
256
257 // density of overall gas phase
258 fluid_density_data.rho_GR =
260 .template value<double>(variables, x_t.x, x_t.t, x_t.dt);
261
262 // derivatives of average molar mass of the gas phase
263 auto const dMG = M_W - M_C;
264 auto const dMG_dpGR = dxnWG_dpGR * dMG;
265
266 variables.molar_mass_derivative = dMG_dpGR;
267
268 // Derivatives of the density of the (composite gas phase) and the
269 // partial densities of its components. The density of the mixture can
270 // be obtained via the property 'IdealGasLawBinaryMixture', for this
271 // purpose the derivatives of the mean molar mass are passed in each
272 // case via the variable_array.
273 cv.drho_GR_dp_GR =
275 .template dValue<double>(
277 x_t.x, x_t.t, x_t.dt);
278
279 auto const dMG_dpCap = dxnWG_dpCap * dMG;
280 variables.molar_mass_derivative = dMG_dpCap;
281 cv.drho_GR_dp_cap =
283 .template dValue<double>(
285 x_t.x, x_t.t, x_t.dt);
286
287 auto const dMG_dT = dxnWG_dT * dMG;
288 variables.molar_mass_derivative = dMG_dT;
289 cv.drho_GR_dT =
291 .template dValue<double>(variables,
293 x_t.x, x_t.t, x_t.dt);
294
295 // The derivatives of the partial densities of the gas phase are
296 // hard-coded (they should remain so, as they are a fundamental part of
297 // this evaporation model). By outsourcing the derivatives of the phase
298 // density to the MPL, a constant phase density can also be assumed, the
299 // derivatives of the partial densities are then unaffected and the
300 // model is still consistent.
301 constituent_density_data.rho_C_GR =
302 mass_mole_fractions_data.xmCG * fluid_density_data.rho_GR;
303 constituent_density_data.rho_W_GR = xmWG * fluid_density_data.rho_GR;
304
305 // 'Air'-component partial density derivatives
306 cv.drho_C_GR_dp_GR = mass_mole_fractions_data.xmCG * cv.drho_GR_dp_GR -
307 cv.dxmWG_dpGR * fluid_density_data.rho_GR;
308 cv.drho_C_GR_dp_cap = mass_mole_fractions_data.xmCG * cv.drho_GR_dp_cap -
309 cv.dxmWG_dpCap * fluid_density_data.rho_GR;
310 cv.drho_C_GR_dT = mass_mole_fractions_data.xmCG * cv.drho_GR_dT -
311 cv.dxmWG_dT * fluid_density_data.rho_GR;
312
313 // Vapour-component partial density derivatives
314 cv.drho_W_GR_dp_GR =
315 xmWG * cv.drho_GR_dp_GR + cv.dxmWG_dpGR * fluid_density_data.rho_GR;
316 cv.drho_W_GR_dp_cap =
317 xmWG * cv.drho_GR_dp_cap + cv.dxmWG_dpCap * fluid_density_data.rho_GR;
318 cv.drho_W_GR_dT =
319 xmWG * cv.drho_GR_dT + cv.dxmWG_dT * fluid_density_data.rho_GR;
320
321 // specific heat capacities of dry air and vapour
322 auto const cpCG =
323 dry_air_component
325 .template value<double>(variables, x_t.x, x_t.t, x_t.dt);
326 auto const cpWG =
327 vapour_component
329 .template value<double>(variables, x_t.x, x_t.t, x_t.dt);
330
331 // specific enthalpy of dry air and vapour components
332 cv.hCG = cpCG * T;
333 cv.hWG = cpWG * T + dh_evap;
334
335 // specific enthalpy of gas phase
336 enthalpy_data.h_G = mass_mole_fractions_data.xmCG * cv.hCG + xmWG * cv.hWG;
337 cv.dh_G_dT = 0;
338
339 // specific inner energies of gas phase
340 cv.uG = enthalpy_data.h_G - pGR / fluid_density_data.rho_GR;
341 cv.du_G_dT = 0;
342 cv.du_G_dp_GR = 0;
343
344 // diffusion
345 auto const tortuosity =
346 media_data.medium
348 .template value<double>(variables, x_t.x, x_t.t, x_t.dt);
349
350 auto const D_W_G_m =
351 vapour_component.property(MaterialPropertyLib::PropertyType::diffusion)
352 .template value<double>(variables, x_t.x, x_t.t, x_t.dt);
353 cv.diffusion_coefficient_vapour =
354 tortuosity * D_W_G_m; // Note here that D_W_G = D_C_G.
355
356 variables.molar_fraction = mass_mole_fractions_data.xnCG;
357
358 // Dissolution part -- Liquid phase properties
359 // -------------------------------------------
360
361 // Reference to the gas component dissolved in the liquid phase
362 auto const& solute_component =
363 liquid_phase.component(liquid_phase_solute_component_index_);
364
365 // The amount of dissolved gas is described by Henry's law. If no
366 // dissolution is intended, the user must define the Henry coefficient
367 // as 'constant 0'.
368 //
369 // Henry-Coefficient and derivatives
370 auto const H =
371 solute_component
373 .template value<double>(variables, x_t.x, x_t.t, x_t.dt);
374
375 auto const dH_dT =
376 solute_component
378 .template dValue<double>(variables,
380 x_t.x, x_t.t, x_t.dt);
381 auto const dH_dpGR =
382 solute_component
384 .template dValue<double>(
386 x_t.x, x_t.t, x_t.dt);
387
388 // Concentration of the dissolved gas as amount of substance of the
389 // mixture component C related to the total volume of the liquid phase.
390 auto const cCL = H * mass_mole_fractions_data.xnCG * pGR;
391 // Fortunately for the developer, the signs of the derivatives of the
392 // composition of binary mixtures are often opposed.
393 auto const dxnCG_dpGR = -dxnWG_dpGR;
394 auto const dxnCG_dT = -dxnWG_dT;
395
396 auto const dcCL_dpGR =
397 (dH_dpGR * mass_mole_fractions_data.xnCG + H * dxnCG_dpGR) * pGR +
398 H * mass_mole_fractions_data.xnCG;
399 auto const dcCL_dT =
400 pGR * (dH_dT * mass_mole_fractions_data.xnCG + H * dxnCG_dT);
401
402 variables.concentration = cCL;
403 // Liquid density including dissolved gas components. Attention! This
404 // only works if the concentration of the C-component is taken into
405 // account in the selected equation of state, e.g. via a (multi)-linear
406 // equation of state. Should a constant density be assumed, no
407 // dissolution will take place, since the composition of the water phase
408 // is determined via the mass fractions (as the ratio of the densities
409 // of solvent and solution). NB! This problem did not occur with the gas
410 // phase because the composition there was determined via the molar
411 // fractions (ratio of the partial pressures).
412 fluid_density_data.rho_LR =
413 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
414 .template value<double>(variables, x_t.x, x_t.t, x_t.dt);
415 variables.density = fluid_density_data.rho_LR;
416
417 // Gas component partial density in liquid phase
418 constituent_density_data.rho_C_LR = fluid_density_data.rho_LR - rho_W_LR();
419
420 // liquid phase composition (mass fraction)
421 mass_mole_fractions_data.xmWL =
422 std::clamp(rho_W_LR() / fluid_density_data.rho_LR, 0., 1.);
423 auto const xmCL = 1. - mass_mole_fractions_data.xmWL;
424
425 // Attention! Usually a multi-linear equation of state is used to
426 // determine the density of the solution. This requires independent
427 // variables, but in this case the concentration is not independent, but
428 // is determined via the equilibrium state. The exact derivation of the
429 // density according to e.g. the pressure pLR is thus:
430 //
431 // drho_d_pGR = rho_ref * (beta_pLR + betaC * dC_dpLR)
432 //
433 // instead of
434 //
435 // d_rho_d_pGR =
436 // liquid_phase.property(MaterialPropertyLib::PropertyType::density)
437 // .template dValue<double>(
438 // variables,
439 // MaterialPropertyLib::Variable::liquid_phase_pressure,
440 // x_t.x, x_t.t, x_t.dt);
441
442 auto const rho_ref_betaP =
443 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
444 .template dValue<double>(
446 x_t.x, x_t.t, x_t.dt);
447
448 auto const rho_ref_betaT =
449 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
450 .template dValue<double>(variables,
452 x_t.x, x_t.t, x_t.dt);
453
454 auto const rho_ref_betaC =
455 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
456 .template dValue<double>(
458 x_t.t, x_t.dt);
459
460 // liquid phase density derivatives
461 auto const drhoLR_dpGR = rho_ref_betaP + rho_ref_betaC * dcCL_dpGR;
462 auto const drhoLR_dpCap = -rho_ref_betaP;
463 cv.drho_LR_dT = rho_ref_betaT + rho_ref_betaC * dcCL_dT;
464
465 // solvent partial density derivatives
466 auto const drhoWLR_dpGR = rho_ref_betaP;
467 auto const drhoWLR_dpCap = -rho_ref_betaP;
468 auto const drhoWLR_dT = rho_ref_betaT;
469
470 // liquid phase mass fraction derivatives
471 cv.dxmWL_dpGR =
472 1. / fluid_density_data.rho_LR *
473 (drhoWLR_dpGR - mass_mole_fractions_data.xmWL * drhoLR_dpGR);
474 cv.dxmWL_dpCap =
475 1. / fluid_density_data.rho_LR *
476 (drhoWLR_dpCap - mass_mole_fractions_data.xmWL * drhoLR_dpCap);
477 cv.dxmWL_dT = 1. / fluid_density_data.rho_LR *
478 (drhoWLR_dT - mass_mole_fractions_data.xmWL * cv.drho_LR_dT);
479
480 // liquid phase molar fractions and derivatives
481 mass_mole_fractions_data.xnWL =
482 mass_mole_fractions_data.xmWL * M_C /
483 (mass_mole_fractions_data.xmWL * M_C + xmCL * M_W);
484
485 // Reference to the pure liquid component
486 auto const& solvent_component =
487 liquid_phase.component(liquid_phase_solvent_component_index_);
488
489 // specific heat capacities of liquid phase components
490 auto const cpCL =
491 solute_component
493 .template value<double>(variables, x_t.x, x_t.t, x_t.dt);
494 auto const cpWL =
495 solvent_component
497 .template value<double>(variables, x_t.x, x_t.t, x_t.dt);
498
499 // specific heat of solution
500 const auto dh_sol =
501 solute_component
503 .template value<double>(variables, x_t.x, x_t.t, x_t.dt);
504
505 // specific enthalpy of liquid phase and its components
506 double const hCL = cpCL * T + dh_sol;
507 double const hWL = cpWL * T;
508 enthalpy_data.h_L = xmCL * hCL + mass_mole_fractions_data.xmWL * hWL;
509 cv.dh_L_dT = 0;
510
511 // specific inner energies of liquid phase
512 cv.uL = enthalpy_data.h_L;
513 cv.du_L_dT = 0;
514
515 // diffusion
516 auto const D_C_L_m =
517 solute_component.property(MaterialPropertyLib::PropertyType::diffusion)
518 .template value<double>(variables, x_t.x, x_t.t, x_t.dt);
519 cv.diffusion_coefficient_solute =
520 tortuosity * D_C_L_m; // Note here that D_C_L = D_W_L.
521
522 // Some default initializations.
523 cv.drho_LR_dp_LR = 0;
524 cv.drho_W_LR_dp_GR = 0.;
525 cv.drho_W_LR_dT = 0.;
526 cv.drho_W_LR_dp_LR = 0.;
527}

References MaterialPropertyLib::capillary_pressure, MaterialPropertyLib::VariableArray::capillary_pressure, MaterialPropertyLib::Phase::component(), MaterialPropertyLib::concentration, MaterialPropertyLib::VariableArray::concentration, MaterialPropertyLib::density, MaterialPropertyLib::VariableArray::density, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::dh_G_dT, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::dh_L_dT, MaterialPropertyLib::diffusion, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::diffusion_coefficient_solute, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::diffusion_coefficient_vapour, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::drho_C_GR_dp_cap, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::drho_C_GR_dp_GR, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::drho_C_GR_dT, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::drho_GR_dp_cap, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::drho_GR_dp_GR, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::drho_GR_dT, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::drho_LR_dp_LR, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::drho_LR_dT, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::drho_W_GR_dp_cap, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::drho_W_GR_dp_GR, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::drho_W_GR_dT, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::drho_W_LR_dp_GR, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::drho_W_LR_dp_LR, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::drho_W_LR_dT, ProcessLib::ConstitutiveRelations::SpaceTimeData::dt, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::du_G_dp_GR, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::du_G_dT, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::du_L_dT, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::dxmWG_dpCap, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::dxmWG_dpGR, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::dxmWG_dT, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::dxmWL_dpCap, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::dxmWL_dpGR, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::dxmWL_dT, MaterialPropertyLib::VariableArray::enthalpy_of_evaporation, ProcessLib::TH2M::ConstitutiveRelations::MediaData::gas, gas_phase_dry_air_component_index_, MaterialPropertyLib::gas_phase_pressure, MaterialPropertyLib::VariableArray::gas_phase_pressure, gas_phase_vapour_component_index_, ProcessLib::TH2M::ConstitutiveRelations::EnthalpyData::h_G, ProcessLib::TH2M::ConstitutiveRelations::EnthalpyData::h_L, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::hCG, MaterialPropertyLib::henry_coefficient, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::hWG, MaterialLib::PhysicalConstant::IdealGasConstant, ProcessLib::TH2M::ConstitutiveRelations::MediaData::liquid, MaterialPropertyLib::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_phase_pressure, liquid_phase_solute_component_index_, liquid_phase_solvent_component_index_, ProcessLib::TH2M::ConstitutiveRelations::MediaData::medium, MaterialPropertyLib::VariableArray::molar_fraction, MaterialPropertyLib::molar_mass, MaterialPropertyLib::VariableArray::molar_mass, MaterialPropertyLib::VariableArray::molar_mass_derivative, MaterialPropertyLib::Component::property(), MaterialPropertyLib::Medium::property(), ProcessLib::TH2M::ConstitutiveRelations::VapourPartialPressureData::pWGR, ProcessLib::TH2M::ConstitutiveRelations::ConstituentDensityData::rho_C_GR, ProcessLib::TH2M::ConstitutiveRelations::ConstituentDensityData::rho_C_LR, ProcessLib::TH2M::ConstitutiveRelations::FluidDensityData::rho_GR, ProcessLib::TH2M::ConstitutiveRelations::FluidDensityData::rho_LR, ProcessLib::TH2M::ConstitutiveRelations::ConstituentDensityData::rho_W_GR, MaterialPropertyLib::specific_heat_capacity, MaterialPropertyLib::specific_latent_heat, ProcessLib::ConstitutiveRelations::SpaceTimeData::t, ProcessLib::TH2M::ConstitutiveRelations::TemperatureData::T, MaterialPropertyLib::temperature, MaterialPropertyLib::VariableArray::temperature, MaterialPropertyLib::tortuosity, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::uG, ProcessLib::TH2M::ConstitutiveRelations::PhaseTransitionData::uL, MaterialPropertyLib::vapour_pressure, ProcessLib::ConstitutiveRelations::SpaceTimeData::x, ProcessLib::TH2M::ConstitutiveRelations::MassMoleFractionsData::xmCG, ProcessLib::TH2M::ConstitutiveRelations::MassMoleFractionsData::xmWL, ProcessLib::TH2M::ConstitutiveRelations::MassMoleFractionsData::xnCG, and ProcessLib::TH2M::ConstitutiveRelations::MassMoleFractionsData::xnWL.

Member Data Documentation

◆ gas_phase_dry_air_component_index_

int const ProcessLib::TH2M::ConstitutiveRelations::PhaseTransition::gas_phase_dry_air_component_index_
private

Definition at line 41 of file PhaseTransition.h.

Referenced by PhaseTransition(), and eval().

◆ gas_phase_vapour_component_index_

int const ProcessLib::TH2M::ConstitutiveRelations::PhaseTransition::gas_phase_vapour_component_index_
private

Definition at line 40 of file PhaseTransition.h.

Referenced by PhaseTransition(), and eval().

◆ liquid_phase_solute_component_index_

int const ProcessLib::TH2M::ConstitutiveRelations::PhaseTransition::liquid_phase_solute_component_index_
private

Definition at line 42 of file PhaseTransition.h.

Referenced by PhaseTransition(), and eval().

◆ liquid_phase_solvent_component_index_

int const ProcessLib::TH2M::ConstitutiveRelations::PhaseTransition::liquid_phase_solvent_component_index_
private

Definition at line 43 of file PhaseTransition.h.

Referenced by PhaseTransition(), and eval().

◆ n_components_gas_

int const ProcessLib::TH2M::ConstitutiveRelations::PhaseTransition::n_components_gas_
private

Definition at line 39 of file PhaseTransition.h.

Referenced by PhaseTransition().


The documentation for this struct was generated from the following files: