Go to the documentation of this file.
30 template<
class Thermo,
template<
class>
class Type>
40 template<
class Thermo,
template<
class>
class Type>
54 scalar Ttol = T0*tol_;
62 (Test - ((this->*F)(
p, Test) -
f)/(this->*dFdT)(
p, Test));
64 if (iter++ > maxIter_)
67 <<
"Maximum number of iterations exceeded"
71 }
while (
mag(Tnew - Test) > Ttol);
79 template<
class Thermo,
template<
class>
class Type>
92 template<
class Thermo,
template<
class>
class Type>
96 return Type<thermo<Thermo, Type> >
::name();
100 template<
class Thermo,
template<
class>
class Type>
104 return Type<thermo<Thermo, Type> >
::he(*
this,
p,
T);
108 template<
class Thermo,
template<
class>
class Type>
112 return this->
cp(p,
T) - this->cpMcv(
p,
T);
116 template<
class Thermo,
template<
class>
class Type>
120 return Type<thermo<Thermo, Type> >::cpv(*
this,
p,
T);
124 template<
class Thermo,
template<
class>
class Type>
128 scalar
cp = this->
cp(p,
T);
129 return cp/(
cp - this->cpMcv(
p,
T));
133 template<
class Thermo,
template<
class>
class Type>
141 return Type<thermo<Thermo, Type> >::cpBycpv(*
this,
p,
T);
145 template<
class Thermo,
template<
class>
class Type>
149 return this->hs(
p,
T) -
p*this->W()/this->
rho(p,
T);
153 template<
class Thermo,
template<
class>
class Type>
157 return this->ha(
p,
T) -
p*this->W()/this->
rho(p,
T);
161 template<
class Thermo,
template<
class>
class Type>
165 return this->ha(
p,
T) -
T*this->
s(p,
T);
169 template<
class Thermo,
template<
class>
class Type>
173 return this->ea(
p,
T) -
T*this->
s(p,
T);
177 template<
class Thermo,
template<
class>
class Type>
181 return this->cpv(
p,
T)/this->W();
185 template<
class Thermo,
template<
class>
class Type>
189 return this->
cp(p,
T)/this->W();
193 template<
class Thermo,
template<
class>
class Type>
197 return this->cv(
p,
T)/this->W();
201 template<
class Thermo,
template<
class>
class Type>
205 return Type<thermo<Thermo, Type> >::HE(*
this,
p,
T);
209 template<
class Thermo,
template<
class>
class Type>
213 return this->hs(
p,
T)/this->W();
217 template<
class Thermo,
template<
class>
class Type>
221 return this->hc()/this->W();
225 template<
class Thermo,
template<
class>
class Type>
229 return this->ha(
p,
T)/this->W();
233 template<
class Thermo,
template<
class>
class Type>
237 return this->
s(p,
T)/this->W();
241 template<
class Thermo,
template<
class>
class Type>
245 return this->
e(p,
T)/this->W();
249 template<
class Thermo,
template<
class>
class Type>
253 return this->es(
p,
T)/this->W();
256 template<
class Thermo,
template<
class>
class Type>
260 return this->ea(
p,
T)/this->W();
264 template<
class Thermo,
template<
class>
class Type>
268 return this->
g(p,
T)/this->W();
272 template<
class Thermo,
template<
class>
class Type>
276 return this->a(
p,
T)/this->W();
280 template<
class Thermo,
template<
class>
class Type>
284 scalar arg = -this->nMoles()*this->
g(
Pstd, T)/(
RR*
T);
297 template<
class Thermo,
template<
class>
class Type>
305 template<
class Thermo,
template<
class>
class Type>
309 if (
equal(this->nMoles(), SMALL))
320 template<
class Thermo,
template<
class>
class Type>
327 if (
equal(this->nMoles(), SMALL))
338 template<
class Thermo,
template<
class>
class Type>
346 if (
equal(this->nMoles(), SMALL))
357 template<
class Thermo,
template<
class>
class Type>
365 return Type<thermo<Thermo, Type> >::THE(*
this,
he,
p, T0);
369 template<
class Thermo,
template<
class>
class Type>
389 template<
class Thermo,
template<
class>
class Type>
409 template<
class Thermo,
template<
class>
class Type>
429 template<
class Thermo,
template<
class>
class Type>
451 template<
class Thermo,
template<
class>
class Type>
457 Thermo::operator+=(st);
461 template<
class Thermo,
template<
class>
class Type>
467 Thermo::operator-=(st);
471 template<
class Thermo,
template<
class>
class Type>
474 Thermo::operator*=(
s);
480 template<
class Thermo,
template<
class>
class Type>
489 static_cast<const Thermo&
>(st1) +
static_cast<const Thermo&
>(st2)
494 template<
class Thermo,
template<
class>
class Type>
503 static_cast<const Thermo&
>(st1) -
static_cast<const Thermo&
>(st2)
508 template<
class Thermo,
template<
class>
class Type>
517 s*
static_cast<const Thermo&
>(st)
522 template<
class Thermo,
template<
class>
class Type>
const scalar RR
Universal gas constant (default in [J/(kmol K)])
scalar ea(const scalar p, const scalar T) const
Absolute internal energy [J/kmol].
A class for handling words, derived from string.
scalar S(const scalar p, const scalar T) const
Entropy [J/(kg K)].
scalar Cpv(const scalar p, const scalar T) const
Heat capacity at constant pressure/volume [J/(kg K)].
scalar a(const scalar p, const scalar T) const
Helmholtz free energy [J/kmol].
const dimensionedScalar Pstd
Standard pressure.
scalar Ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kg].
scalar he(const scalar p, const scalar T) const
Enthalpy/Internal energy [J/kmol].
scalar E(const scalar p, const scalar T) const
Internal energy [J/kg].
scalar Hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
scalar Hc() const
Chemical enthalpy [J/kg].
thermo(const Thermo &sp)
Construct from components.
const dimensionedVector & g
dimensionedScalar exp(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
scalar THs(const scalar Hs, const scalar p, const scalar T0) const
Temperature from sensible enthalpy given an initial T0.
scalar Kp(const scalar p, const scalar T) const
Equilibrium constant [] i.t.o. partial pressures.
const dimensionedScalar F
Faraday constant: default SI units: [C/mol].
void operator*=(const scalar)
CGAL::Exact_predicates_exact_constructions_kernel K
scalar TEa(const scalar E, const scalar p, const scalar T0) const
Temperature from absolute internal energy.
static word heName()
Name of Enthalpy/Internal energy.
scalar Cv(const scalar p, const scalar T) const
Heat capacity at constant volume [J/(kg K)].
scalar THE(const scalar H, const scalar p, const scalar T0) const
Temperature from enthalpy or internal energy.
scalar K(const scalar p, const scalar T) const
Equilibrium constant [] i.t.o fugacities.
scalar cpv(const scalar p, const scalar T) const
Heat capacity at constant pressure/volume [J/(kmol K)].
const volScalarField & cp
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void limit(const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su, const scalar psiMax, const scalar psiMin, const bool returnCorr)
scalar es(const scalar p, const scalar T) const
Sensible internal energy [J/kmol].
errorManip< error > abort(error &err)
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
scalar G(const scalar p, const scalar T) const
Gibbs free energy [J/kg].
scalar Kn(const scalar p, const scalar T, const scalar n) const
Equilibrium constant [] i.t.o. number of moles.
scalar Cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kg K)].
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
scalar A(const scalar p, const scalar T) const
Helmholtz free energy [J/kg].
scalar TEs(const scalar E, const scalar p, const scalar T0) const
Temperature from sensible internal energy.
scalar Kx(const scalar p, const scalar T) const
Equilibrium constant [] i.t.o. mole-fractions.
const dimensionedScalar e
Elementary charge.
scalar T(scalar f, scalar p, scalar T0, scalar(thermo::*F)(const scalar, const scalar) const, scalar(thermo::*dFdT)(const scalar, const scalar) const, scalar(thermo::*limit)(const scalar) const) const
Return the temperature corresponding to the value of the.
scalar cpBycpv(const scalar p, const scalar T) const
Ratio of heat capacity at constant pressure to that at.
scalar Es(const scalar p, const scalar T) const
Sensible internal energy [J/kg].
scalar gamma(const scalar p, const scalar T) const
Gamma = cp/cv [].
scalar g(const scalar p, const scalar T) const
Gibbs free energy [J/kmol].
scalar HE(const scalar p, const scalar T) const
Enthalpy/Internal energy [J/kg].
scalar Kc(const scalar p, const scalar T) const
Equilibrium constant i.t.o. molar concentration.
scalar THa(const scalar H, const scalar p, const scalar T0) const
Temperature from absolute enthalpy.
bool equal(const T &s1, const T &s2)
scalar cv(const scalar p, const scalar T) const
Heat capacity at constant volume [J/(kmol K)].
word name(const complex &)
Return a string representation of a complex.
scalar Ea(const scalar p, const scalar T) const
Absolute internal energy [J/kg].