Canoe
Comprehensive Atmosphere N' Ocean Engine
|
#include <thermodynamics.hpp>
Public Types | |
enum | { Size = 1 + NVAPOR + NCLOUD } |
enum class | Method { ReversibleAdiabat = 0 , PseudoAdiabat = 1 , DryAdiabat = 2 , Isothermal = 3 , NeutralStability = 4 } |
using | SVPFunc1Container = std::vector< std::vector< SatVaporPresFunc1 > > |
using | SVPFunc2Container = std::map< IndexPair, SatVaporPresFunc2 > |
Public Member Functions | |
~Thermodynamics () | |
Real | GetRd () const |
Real | GetGammad (AirParcel const &var) const |
adiaibatic index of dry air [1] More... | |
Real | GetGammadRef () const |
reference adiabatic index of dry air [1] More... | |
Real | GetMu (int n) const |
Real | GetInvMu (int n) const |
Real | GetMuRatio (int n) const |
Real | GetInvMuRatio (int n) const |
Real | GetCvRatioMass (int n) const |
Real | GetCvRatioMole (int n) const |
int | GetCloudIndex (int i, int j) const |
Real | GetCvMass (AirParcel const &qfrac, int n) const |
Real | GetCvMassRef (int n) const |
Reference specific heat capacity [J/(kg K)] at constant volume. More... | |
Real | GetCvMole (AirParcel const &qfrac, int n) const |
Real | GetCpRatioMass (int n) const |
Real | GetCpRatioMole (int n) const |
Real | GetCpMass (AirParcel const &qfrac, int n) const |
Real | GetCpMassRef (int n) const |
Real | GetCpMole (AirParcel const &qfrac, int n) const |
Real | GetChi (AirParcel const &qfrac) const |
Real | GetLatentEnergyMass (int n, Real temp=0.) const |
Temperature dependent specific latent energy [J/kg] of condensates at constant volume \(L_{ij}(T) = L_{ij}^r - (c_{ij} - c_{p,i})\times(T - T^r)\). More... | |
Real | GetLatentEnergyMole (int n, Real temp=0.) const |
Temperature dependent specific latent energy [J/mol] of condensates at constant volume. More... | |
Real | GetLatentHeatMole (int i, std::vector< Real > const &rates, Real temp) const |
Real | GetLatentHeatMass (int i, std::vector< Real > const &rates, Real temp) const |
RealArrayX | TryEquilibriumTP_VaporCloud (AirParcel const &qfrac, int ivapor, Real cv_hat=0., bool misty=false) const |
Calculate the equilibrium mole transfer by cloud reaction vapor -> cloud. More... | |
RealArray3 | TryEquilibriumTP_VaporVaporCloud (AirParcel const &air, IndexPair ij, Real cv_hat=0., bool misty=false) const |
Calculate the equilibrium mole transfer by cloud reaction vapor1 + vapor2 -> cloud. More... | |
void | Extrapolate (AirParcel *qfrac, Real dzORdlnp, Method method, Real grav=0., Real userp=0.) const |
void | EquilibrateTP (AirParcel *qfrac) const |
void | SaturationAdjustment (AirColumn &ac) const |
Real | RovRd (AirParcel const &qfrac) const |
Real | PotentialTemp (MeshBlock *pmb, Real p0, int k, int j, int i) const |
Calculate potential temperature from primitive variable. More... | |
Real | EquivalentPotentialTemp (MeshBlock *pmb, Real p0, int v, int k, int j, int i) const |
Calculate equivalent potential temperature from primitive variable. More... | |
Real | GetTemp (MeshBlock const *pmb, int k, int j, int i) const |
Calculate temperature from primitive variable. More... | |
Real | GetMu (MeshBlock const *pmb, int k, int j, int i) const |
Mean molecular weight. More... | |
Real | RovRd (MeshBlock const *pmb, int k, int j, int i) const |
Inverse of the mean molecular weight (no cloud) More... | |
Real | GetCvMass (MeshBlock const *pmb, int k, int j, int i) const |
Specific heat capacity [J/(kg K)] of the air parcel at constant volume. More... | |
Real | GetCpMass (MeshBlock const *pmb, int k, int j, int i) const |
Specific heat capacity [J/(kg K)] of the air parcel at constant pressure. More... | |
Real | GetChi (MeshBlock const *pmb, int k, int j, int i) const |
Adiabatic index. More... | |
Real | GetGamma (MeshBlock const *pmb, int k, int j, int i) const |
Polytropic index. More... | |
Real | GetPres (MeshBlock const *pmb, int k, int j, int i) const |
Real | GetEnthalpyMass (MeshBlock const *pmb, int k, int j, int i) const |
specific enthalpy [J/kg] of the air parcel More... | |
Real | MoistStaticEnergy (MeshBlock const *pmb, Real gz, int k, int j, int i) const |
Moist static energy. More... | |
Real | RelativeHumidity (MeshBlock const *pmb, int n, int k, int j, int i) const |
Relative humidity. More... | |
Real | RelativeHumidity (AirParcel const &qfrac, int n) const |
Static Public Member Functions | |
static Thermodynamics const * | GetInstance () |
Return a pointer to the one and only instance of Thermodynamics. More... | |
static Thermodynamics const * | InitFromAthenaInput (ParameterInput *pin) |
static void | Destroy () |
Destroy the one and only instance of Thermodynamics. More... | |
Static Public Attributes | |
static constexpr Real | RefTemp = 300. |
static constexpr Real | RefPres = 1.e5 |
Protected Types | |
enum | { NPHASE = NPHASE_LEGACY } |
Protected Member Functions | |
Thermodynamics () | |
void | updateTPConservingU (AirParcel *qfrac, Real rmole, Real umole) const |
update T/P More... | |
Real | getInternalEnergyMole (AirParcel const &qfrac) const |
Real | getDensityMole (AirParcel const &qfrac) const |
void | setTotalEquivalentVapor (AirParcel *qfrac) const |
Real | calDlnTDlnP (AirParcel const &qfrac, Real latent[]) const |
Calculate moist adiabatic temperature gradient. More... | |
void | rk4IntegrateLnp (AirParcel *qfrac, Real dlnp, Method method, Real adlnTdlnP) const |
void | rk4IntegrateZ (AirParcel *qfrac, Real dlnp, Method method, Real grav, Real adlnTdlnP) const |
void | enrollVaporFunctions () |
custom functions for vapor (to be overridden by user mods) More... | |
Static Protected Member Functions | |
static Thermodynamics * | fromLegacyInput (ParameterInput *pin) |
static Thermodynamics * | fromYAMLInput (YAML::Node node) |
Protected Attributes | |
Real | Rd_ |
ideal gas constant of dry air in J/kg More... | |
Real | gammad_ref_ |
reference polytropic index of dry air More... | |
std::array< Real, Size > | mu_ratio_ |
ratio of mean molecular weights More... | |
std::array< Real, Size > | inv_mu_ratio_ |
inverse ratio of mean molecular weights More... | |
std::array< Real, Size > | mu_ |
molecular weight [kg/mol] More... | |
std::array< Real, Size > | inv_mu_ |
inverse molecular weight [mol/kg] More... | |
std::array< Real, Size > | cp_ratio_mole_ |
ratio of specific heat capacities [J/mol] at constant pressure More... | |
std::array< Real, Size > | cp_ratio_mass_ |
ratio of specific heat capacities [J/mol] at constant pressure More... | |
std::array< Real, Size > | cv_ratio_mole_ |
ratio of specific heat capacities [J/mol] at constant volume More... | |
std::array< Real, Size > | cv_ratio_mass_ |
ratio of specific heat capacities [J/mol] at constant volume More... | |
std::array< Real, Size > | latent_energy_mole_ |
latent energy at constant volume [J/mol] More... | |
std::array< Real, Size > | latent_energy_mass_ |
latent energy at constant volume [J/kg] More... | |
std::array< Real, Size > | beta_ |
Dimensionless latent heat. More... | |
std::array< Real, Size > | delta_ |
Dimensionless differences in specific heat capacity. More... | |
std::array< Real, 1+NVAPOR > | t3_ |
triple point temperature [K] More... | |
std::array< Real, 1+NVAPOR > | p3_ |
triple point pressure [pa] More... | |
SVPFunc1Container | svp_func1_ |
saturation vapor pressure function: Vapor -> Cloud More... | |
std::vector< IndexSet > | cloud_index_set_ |
cloud index set More... | |
SVPFunc2Container | svp_func2_ |
saturation vapor pressure function: Vapor + Vapor -> Cloud More... | |
std::map< IndexPair, ReactionInfo > | cloud_reaction_map_ |
reaction information map More... | |
int | sa_max_iter_ = 5 |
maximum saturation adjustment iterations More... | |
Real | sa_relax_ = 0.8 |
saturation adjustment relaxation parameter More... | |
Real | sa_ftol_ = 0.01 |
saturation adjustment temperature tolerance More... | |
Static Protected Attributes | |
static Thermodynamics * | mythermo_ = nullptr |
pointer to the single Thermodynamics instance More... | |
Definition at line 71 of file thermodynamics.hpp.
using Thermodynamics::SVPFunc1Container = std::vector<std::vector<SatVaporPresFunc1> > |
Definition at line 82 of file thermodynamics.hpp.
using Thermodynamics::SVPFunc2Container = std::map<IndexPair, SatVaporPresFunc2> |
Definition at line 83 of file thermodynamics.hpp.
|
protected |
Enumerator | |
---|---|
NPHASE |
Definition at line 73 of file thermodynamics.hpp.
anonymous enum |
Enumerator | |
---|---|
Size |
Definition at line 85 of file thermodynamics.hpp.
|
strong |
Enumerator | |
---|---|
ReversibleAdiabat | |
PseudoAdiabat | |
DryAdiabat | |
Isothermal | |
NeutralStability |
Definition at line 87 of file thermodynamics.hpp.
|
inlineprotected |
Constructor for class sets up the initial conditions Protected ctor access thru static member function Instance
Definition at line 77 of file thermodynamics.hpp.
Thermodynamics::~Thermodynamics | ( | ) |
|
staticprotected |
Definition at line 58 of file thermodynamics.cpp.
|
staticprotected |
Definition at line 52 of file thermodynamics.cpp.
|
static |
Return a pointer to the one and only instance of Thermodynamics.
Definition at line 100 of file thermodynamics.cpp.
|
static |
Definition at line 111 of file thermodynamics.cpp.
|
static |
Destroy the one and only instance of Thermodynamics.
Definition at line 237 of file thermodynamics.cpp.
|
inline |
Ideal gas constant of dry air in [J/(kg K)]
Definition at line 111 of file thermodynamics.hpp.
Real Thermodynamics::GetGammad | ( | AirParcel const & | var | ) | const |
adiaibatic index of dry air [1]
|
inline |
reference adiabatic index of dry air [1]
Definition at line 117 of file thermodynamics.hpp.
|
inline |
Molecular weight [kg/mol]
[in] | n | the index of the thermodynamic species |
Definition at line 122 of file thermodynamics.hpp.
|
inline |
Inverse molecular weight [mol/kg]
[in] | n | the index of the thermodynamic species |
Definition at line 127 of file thermodynamics.hpp.
|
inline |
Ratio of molecular weights [1]
[in] | n | the index of the thermodynamic species |
Definition at line 132 of file thermodynamics.hpp.
|
inline |
Ratio of molecular weights [1]
[in] | n | the index of the thermodynamic species |
Definition at line 137 of file thermodynamics.hpp.
|
inline |
Ratio of specific heat capacity [J/(kg K)] at constant volume [1]
[in] | n | the index of the thermodynamic species |
Definition at line 142 of file thermodynamics.hpp.
|
inline |
Ratio of specific heat capacity [J/(mol K)] at constant volume [1]
[in] | n | the index of the thermodynamic species |
Definition at line 147 of file thermodynamics.hpp.
|
inline |
[in] | i | the index of the vapor |
[in] | j | the sequential index of the cloud |
Definition at line 152 of file thermodynamics.hpp.
|
inline |
Specific heat capacity [J/(kg K)] at constant volume \(c_{v,d} = \frac{R_d}{\gamma_d - 1}\)
\(c_{v,i} = \frac{c_{v,i}}{c_{v,d}}\times c_{v,d}\)
Definition at line 158 of file thermodynamics.hpp.
|
inline |
Reference specific heat capacity [J/(kg K)] at constant volume.
Definition at line 164 of file thermodynamics.hpp.
|
inline |
Specific heat capacity [J/(mol K)] of the air parcel at constant volume
Definition at line 171 of file thermodynamics.hpp.
|
inline |
Ratio of specific heat capacity [J/(kg K)] at constant pressure
Definition at line 177 of file thermodynamics.hpp.
|
inline |
Ratio of specific heat capacity [J/(mol K)] at constant pressure
Definition at line 181 of file thermodynamics.hpp.
|
inline |
Specific heat capacity [J/(kg K)] at constant pressure \(c_{p,d} = \frac{\gamma_d}{\gamma_d - 1}R_d\)
\(c_{p,i} = \frac{c_{p,i}}{c_{p,d}}\times c_{p,d}\)
Definition at line 187 of file thermodynamics.hpp.
|
inline |
Definition at line 193 of file thermodynamics.hpp.
|
inline |
Specific heat capacity [J/(mol K)] of the air parcel at constant pressure
Definition at line 200 of file thermodynamics.hpp.
Real Thermodynamics::GetChi | ( | AirParcel const & | qfrac | ) | const |
|
inline |
Temperature dependent specific latent energy [J/kg] of condensates at constant volume \(L_{ij}(T) = L_{ij}^r - (c_{ij} - c_{p,i})\times(T - T^r)\).
\(= L_{ij}^r - \delta_{ij}R_i(T - T^r)\)
Definition at line 214 of file thermodynamics.hpp.
|
inline |
Temperature dependent specific latent energy [J/mol] of condensates at constant volume.
Definition at line 222 of file thermodynamics.hpp.
Real Thermodynamics::GetLatentHeatMole | ( | int | i, |
std::vector< Real > const & | rates, | ||
Real | temp | ||
) | const |
Definition at line 441 of file thermodynamics.cpp.
|
inline |
Definition at line 229 of file thermodynamics.hpp.
RealArrayX Thermodynamics::TryEquilibriumTP_VaporCloud | ( | AirParcel const & | qfrac, |
int | ivapor, | ||
Real | cv_hat = 0. , |
||
bool | misty = false |
||
) | const |
Calculate the equilibrium mole transfer by cloud reaction vapor -> cloud.
[in] | qfrac | mole fraction representation of air parcel |
[in] | ivapor | the index of the vapor |
[in] | cv_hat | \(cv_hat\) molar heat capacity |
[in] | misty | if true, there is an infinite supple of cloud |
Definition at line 16 of file try_equilibrium_tp_vapor_cloud.cpp.
RealArray3 Thermodynamics::TryEquilibriumTP_VaporVaporCloud | ( | AirParcel const & | air, |
IndexPair | ij, | ||
Real | cv_hat = 0. , |
||
bool | misty = false |
||
) | const |
Calculate the equilibrium mole transfer by cloud reaction vapor1 + vapor2 -> cloud.
[in] | air | mole fraction representation of air parcel |
[in] | ij | the index pair of the vapor1 and vapor2 |
[in] | cv_hat | \(\har{cv}\) molar heat capacity |
[in] | misty | if true, there is an infinite supple of cloud |
Definition at line 27 of file try_equilibrium_tp_vapor_vapor_cloud.cpp.
void Thermodynamics::Extrapolate | ( | AirParcel * | qfrac, |
Real | dzORdlnp, | ||
Method | method, | ||
Real | grav = 0. , |
||
Real | userp = 0. |
||
) | const |
Construct an 1d atmosphere
[in,out] | qfrac | mole fraction representation of air parcel |
[in] | dzORdlnp | vertical grid spacing |
[in] | method | choose from [reversible, pseudo, dry, isothermal] |
[in] | grav | gravitational acceleration |
[in] | userp | user parameter to adjust the temperature gradient |
Definition at line 417 of file thermodynamics.cpp.
void Thermodynamics::EquilibrateTP | ( | AirParcel * | qfrac | ) | const |
Thermodnamic equilibrium at current TP
[in,out] | qfrac | mole fraction representation of air parcel |
Definition at line 7 of file equilibrate_tp.cpp.
void Thermodynamics::SaturationAdjustment | ( | AirColumn & | ac | ) | const |
Adjust to the maximum saturation state conserving internal energy
[in,out] | ac | mole fraction representation of a collection of air parcels |
Definition at line 20 of file saturation_adjustment.cpp.
Real Thermodynamics::RovRd | ( | AirParcel const & | qfrac | ) | const |
Inverse of the mean molecular weight (with cloud)
[in] | qfrac | mole fraction representation of air parcel |
Definition at line 322 of file thermodynamics.cpp.
|
inline |
Calculate potential temperature from primitive variable.
\(\theta = T(\frac{p_0}{p})^{\chi}\)
Definition at line 284 of file thermodynamics.hpp.
Real Thermodynamics::EquivalentPotentialTemp | ( | MeshBlock * | pmb, |
Real | p0, | ||
int | v, | ||
int | k, | ||
int | j, | ||
int | i | ||
) | const |
Calculate equivalent potential temperature from primitive variable.
\(\theta_e = T(\frac{p}{p_d})^{Rd/(cpd + cl r_t} \exp(\frac{L_v q_v}{c_p T})\)
Definition at line 455 of file thermodynamics.cpp.
|
inline |
Calculate temperature from primitive variable.
\(T = p/(\rho R) = p/(\rho \frac{R}{R_d} Rd)\)
Definition at line 301 of file thermodynamics.hpp.
Real Thermodynamics::GetMu | ( | MeshBlock const * | pmb, |
int | k, | ||
int | j, | ||
int | i | ||
) | const |
Mean molecular weight.
\(mu = \mu_d (1 + \sum_i q_i (\epsilon_i - 1))\)
Definition at line 401 of file thermodynamics.cpp.
|
inline |
Inverse of the mean molecular weight (no cloud)
\( \frac{R}{R_d} = \frac{\mu_d}{\mu}\)
Definition at line 317 of file thermodynamics.hpp.
Real Thermodynamics::GetCvMass | ( | MeshBlock const * | pmb, |
int | k, | ||
int | j, | ||
int | i | ||
) | const |
Specific heat capacity [J/(kg K)] of the air parcel at constant volume.
\(c_v = c_{v,d}*(1 + \sum_i (q_i*(\hat{c}_{v,i} - 1.))) = \gamma_d/(\gamma_d - 1.)*R_d*T*(1 + \sum_i (q_i*(\hat{c}_{v,i} - 1.)))\)
Definition at line 377 of file thermodynamics.cpp.
Real Thermodynamics::GetCpMass | ( | MeshBlock const * | pmb, |
int | k, | ||
int | j, | ||
int | i | ||
) | const |
Specific heat capacity [J/(kg K)] of the air parcel at constant pressure.
\(c_p = c_{p,d}*(1 + \sum_i (q_i*(\hat{c}_{p,i} - 1.))) = \gamma_d/(\gamma_d - 1.)*R_d*T*(1 + \sum_i (q_i*(\hat{c}_{p,i} - 1.)))\)
Definition at line 364 of file thermodynamics.cpp.
Real Thermodynamics::GetChi | ( | MeshBlock const * | pmb, |
int | k, | ||
int | j, | ||
int | i | ||
) | const |
Adiabatic index.
\(\chi = \frac{R}{c_p}\)
Definition at line 267 of file thermodynamics.cpp.
Real Thermodynamics::GetGamma | ( | MeshBlock const * | pmb, |
int | k, | ||
int | j, | ||
int | i | ||
) | const |
Polytropic index.
\(\gamma = \frac{c_p}{c_v}\)
Definition at line 301 of file thermodynamics.cpp.
Real Thermodynamics::GetPres | ( | MeshBlock const * | pmb, |
int | k, | ||
int | j, | ||
int | i | ||
) | const |
Real Thermodynamics::GetEnthalpyMass | ( | MeshBlock const * | pmb, |
int | k, | ||
int | j, | ||
int | i | ||
) | const |
specific enthalpy [J/kg] of the air parcel
\(h = c_{p,d}*T*(1 + \sum_i (q_i*(\hat{c}_{p,i} - 1.)))\) \( = \gamma_d/(\gamma_d - 1.)*R_d*T*(1 + \sum_i (q_i*(\hat{c}_{pi} - 1.)))\)
Definition at line 389 of file thermodynamics.cpp.
Real Thermodynamics::MoistStaticEnergy | ( | MeshBlock const * | pmb, |
Real | gz, | ||
int | k, | ||
int | j, | ||
int | i | ||
) | const |
Moist static energy.
\(h_s = c_{pd}T + gz + L_vq_v + L_s\sum_i q_i\)
Definition at line 339 of file thermodynamics.cpp.
Real Thermodynamics::RelativeHumidity | ( | MeshBlock const * | pmb, |
int | n, | ||
int | k, | ||
int | j, | ||
int | i | ||
) | const |
Relative humidity.
\(H_i = \frac{e_i}{e_i^s}\)
Definition at line 410 of file thermodynamics.cpp.
|
inline |
Definition at line 380 of file thermodynamics.hpp.
|
protected |
update T/P
Definition at line 505 of file thermodynamics.cpp.
|
protected |
Definition at line 532 of file thermodynamics.cpp.
|
protected |
Definition at line 553 of file thermodynamics.cpp.
|
protected |
Definition at line 560 of file thermodynamics.cpp.
|
protected |
Calculate moist adiabatic temperature gradient.
\(\Gamma_m = (\frac{d\ln T}{d\ln P})_m\)
Definition at line 7 of file cal_dlnT_dlnP.cpp.
|
protected |
Definition at line 8 of file rk4_integrate_lnp.cpp.
|
protected |
Definition at line 9 of file rk4_integrate_z.cpp.
|
protected |
custom functions for vapor (to be overridden by user mods)
Definition at line 104 of file giants_enroll_vapor_functions_v1.hpp.
|
staticconstexpr |
Definition at line 95 of file thermodynamics.hpp.
|
staticconstexpr |
Definition at line 96 of file thermodynamics.hpp.
|
protected |
ideal gas constant of dry air in J/kg
Definition at line 411 of file thermodynamics.hpp.
|
protected |
reference polytropic index of dry air
Definition at line 414 of file thermodynamics.hpp.
|
protected |
ratio of mean molecular weights
Definition at line 417 of file thermodynamics.hpp.
|
protected |
inverse ratio of mean molecular weights
Definition at line 420 of file thermodynamics.hpp.
|
protected |
molecular weight [kg/mol]
Definition at line 423 of file thermodynamics.hpp.
|
protected |
inverse molecular weight [mol/kg]
Definition at line 426 of file thermodynamics.hpp.
|
protected |
ratio of specific heat capacities [J/mol] at constant pressure
Definition at line 429 of file thermodynamics.hpp.
|
protected |
ratio of specific heat capacities [J/mol] at constant pressure
Definition at line 432 of file thermodynamics.hpp.
|
protected |
ratio of specific heat capacities [J/mol] at constant volume
Definition at line 435 of file thermodynamics.hpp.
|
protected |
ratio of specific heat capacities [J/mol] at constant volume
Definition at line 438 of file thermodynamics.hpp.
|
protected |
latent energy at constant volume [J/mol]
Definition at line 441 of file thermodynamics.hpp.
|
protected |
latent energy at constant volume [J/kg]
\(L_i^r+\Delta c_{ij}T^r == \beta_{ij}\frac{R_d}{\epsilon_i}T^r\)
Definition at line 446 of file thermodynamics.hpp.
|
protected |
Dimensionless latent heat.
\(\beta_{ij} = \frac{\Delta U_{ij}}{R_i T^r}\) \(\Delta U_{ij}\) is the difference in internal energy between the vapor \(i\) and the condensate \(j\) \(R_i=\hat{R}/m_i=R_d/\epsilon_i\) \(T^r\) is the triple point temperature \(\beta_{ij} == \frac{L_{ij}^r + \Delta c_{ij}T^r}{R_i T^r}\)
Definition at line 455 of file thermodynamics.hpp.
|
protected |
Dimensionless differences in specific heat capacity.
\(\delta_{ij} = \frac{c_{ij} - c_{p,i}}{R_i}\) \(c_{ij} - c_{p,j}\) is the difference in specific heat capacity at constant pressure. \(c_{ij}\) is the heat capacity of condensate \(j\) of vapor \(i\)
Definition at line 463 of file thermodynamics.hpp.
|
protected |
triple point temperature [K]
Definition at line 466 of file thermodynamics.hpp.
|
protected |
triple point pressure [pa]
Definition at line 469 of file thermodynamics.hpp.
|
protected |
saturation vapor pressure function: Vapor -> Cloud
Definition at line 472 of file thermodynamics.hpp.
|
protected |
cloud index set
Definition at line 475 of file thermodynamics.hpp.
|
protected |
saturation vapor pressure function: Vapor + Vapor -> Cloud
Definition at line 478 of file thermodynamics.hpp.
|
protected |
reaction information map
Definition at line 481 of file thermodynamics.hpp.
|
staticprotected |
pointer to the single Thermodynamics instance
Definition at line 484 of file thermodynamics.hpp.
|
protected |
maximum saturation adjustment iterations
Definition at line 487 of file thermodynamics.hpp.
|
protected |
saturation adjustment relaxation parameter
Definition at line 490 of file thermodynamics.hpp.
|
protected |
saturation adjustment temperature tolerance
Definition at line 493 of file thermodynamics.hpp.