33template<
class Thermo,
template<
class>
class Type>
43template<
class Thermo,
template<
class>
class Type>
44inline Foam::scalar Foam::species::thermo<Thermo, Type>::T
58 <<
"Negative initial temperature T0: " <<
T0
64 scalar Ttol =
T0*tol_;
72 (Test - ((this->*F)(
p, Test) -
f)/(this->*dFdT)(
p, Test));
74 if (iter++ > maxIter_)
77 <<
"Maximum number of iterations exceeded: " << maxIter_
78 <<
" when starting from T0:" <<
T0
79 <<
" old T:" << Test <<
" new T:" << Tnew
86 }
while (
mag(Tnew - Test) > Ttol);
94template<
class Thermo,
template<
class>
class Type>
107template<
class Thermo,
template<
class>
class Type>
111 return Type<thermo<Thermo, Type>>::energyName();
115template<
class Thermo,
template<
class>
class Type>
119 return Type<thermo<Thermo, Type>>
::Cpv(*
this,
p, T);
123template<
class Thermo,
template<
class>
class Type>
128 volatile const scalar
Cp = this->
Cp(p,
T);
130 const scalar
Cp = this->
Cp(p,
T);
133 return Cp/(
Cp - this->CpMCv(
p, T));
137template<
class Thermo,
template<
class>
class Type>
145 return Type<thermo<Thermo, Type>>
::CpByCpv(*
this,
p, T);
149template<
class Thermo,
template<
class>
class Type>
153 return Type<thermo<Thermo, Type>>
::HE(*
this,
p, T);
157template<
class Thermo,
template<
class>
class Type>
161 return this->
Ha(p, T) - T*this->S(
p, T);
165template<
class Thermo,
template<
class>
class Type>
169 return this->
Ea(p, T) - T*this->S(
p, T);
173template<
class Thermo,
template<
class>
class Type>
177 return this->
Cp(p, T)*this->W();
181template<
class Thermo,
template<
class>
class Type>
185 return this->
Ha(p, T)*this->W();
189template<
class Thermo,
template<
class>
class Type>
193 return this->
Hs(p, T)*this->W();
197template<
class Thermo,
template<
class>
class Type>
201 return this->Hc()*this->W();
205template<
class Thermo,
template<
class>
class Type>
209 return this->S(
p, T)*this->W();
213template<
class Thermo,
template<
class>
class Type>
217 return this->
HE(
p, T)*this->W();
221template<
class Thermo,
template<
class>
class Type>
225 return this->
Cv(p, T)*this->W();
229template<
class Thermo,
template<
class>
class Type>
233 return this->
Es(p, T)*this->W();
237template<
class Thermo,
template<
class>
class Type>
241 return this->
Ea(p, T)*this->W();
245template<
class Thermo,
template<
class>
class Type>
249 return this->
G(
p, T)*this->W();
253template<
class Thermo,
template<
class>
class Type>
257 return this->
A(p, T)*this->W();
261template<
class Thermo,
template<
class>
class Type>
265 scalar arg = -this->
Y()*this->G(
Pstd,
T)/(
RR*
T);
278template<
class Thermo,
template<
class>
class Type>
286template<
class Thermo,
template<
class>
class Type>
290 const scalar nm = this->
Y()/this->W();
292 if (
equal(nm, SMALL))
303template<
class Thermo,
template<
class>
class Type>
310 const scalar nm = this->
Y()/this->W();
312 if (
equal(nm, SMALL))
323template<
class Thermo,
template<
class>
class Type>
331 const scalar nm = this->
Y()/this->W();
333 if (
equal(nm, SMALL))
344template<
class Thermo,
template<
class>
class Type>
352 return Type<thermo<Thermo, Type>>
::THE(*
this,
he,
p,
T0);
356template<
class Thermo,
template<
class>
class Type>
376template<
class Thermo,
template<
class>
class Type>
396template<
class Thermo,
template<
class>
class Type>
416template<
class Thermo,
template<
class>
class Type>
436template<
class Thermo,
template<
class>
class Type>
444 const scalar dKcdTbyKc =
445 (this->S(
Pstd,
T) + this->Gstd(
T)/
T)*this->
Y()/(
RR*
T);
447 const scalar nm = this->
Y()/this->W();
448 if (
equal(nm, SMALL))
459template<
class Thermo,
template<
class>
class Type>
463 return this->dCpdT(
p, T)*this->W();
468template<
class Thermo,
template<
class>
class Type>
474 Thermo::operator+=(st);
478template<
class Thermo,
template<
class>
class Type>
481 Thermo::operator*=(
s);
487template<
class Thermo,
template<
class>
class Type>
496 static_cast<const Thermo&
>(st1) +
static_cast<const Thermo&
>(st2)
501template<
class Thermo,
template<
class>
class Type>
510 s*
static_cast<const Thermo&
>(st)
515template<
class Thermo,
template<
class>
class Type>
516inline Foam::species::thermo<Thermo, Type> Foam::species::operator==
524 static_cast<const Thermo&
>(st1) ==
static_cast<const Thermo&
>(st2)
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
CGAL::Exact_predicates_exact_constructions_kernel K
scalar Ha(const scalar p, const scalar T) const
scalar Hs(const scalar p, const scalar T) const
scalar Es(const scalar p, const scalar T) const
scalar Ea(const scalar p, const scalar T) const
scalar g(const scalar p, const scalar T) const
Gibbs free energy [J/kmol].
scalar Kx(const scalar p, const scalar T) const
Equilibrium constant [] i.t.o. mole-fractions.
scalar Kc(const scalar p, const scalar T) const
Equilibrium constant i.t.o. molar concentration.
scalar G(const scalar p, const scalar T) const
Gibbs free energy [J/kg].
scalar TEs(const scalar E, const scalar p, const scalar T0) const
Temperature from sensible internal energy.
scalar THE(const scalar H, const scalar p, const scalar T0) const
Temperature from enthalpy or internal energy.
scalar HE(const scalar p, const scalar T) const
Enthalpy/Internal energy [J/kg].
scalar s(const scalar p, const scalar T) const
Entropy [J/(kmol K)].
scalar a(const scalar p, const scalar T) const
Helmholtz free energy [J/kmol].
scalar TEa(const scalar E, const scalar p, const scalar T0) const
Temperature from absolute internal energy.
scalar hc() const
Chemical enthalpy [J/kmol].
static word heName()
Name of Enthalpy/Internal energy.
scalar he(const scalar p, const scalar T) const
Enthalpy/Internal energy [J/kmol].
scalar ea(const scalar p, const scalar T) const
Absolute internal energy [J/kmol].
scalar K(const scalar p, const scalar T) const
Equilibrium constant [] i.t.o fugacities.
scalar cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kmol K)].
scalar cv(const scalar p, const scalar T) const
Heat capacity at constant volume [J/(kmol K)].
thermo(const Thermo &sp)
Construct from components.
scalar CpByCpv(const scalar p, const scalar T) const
Ratio of heat capacity at constant pressure to that at.
scalar Cpv(const scalar p, const scalar T) const
Heat capacity at constant pressure/volume [J/(kg K)].
scalar ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kmol].
scalar A(const scalar p, const scalar T) const
Helmholtz free energy [J/kg].
scalar gamma(const scalar p, const scalar T) const
Gamma = Cp/Cv [].
scalar es(const scalar p, const scalar T) const
Sensible internal energy [J/kmol].
scalar THa(const scalar H, const scalar p, const scalar T0) const
Temperature from absolute enthalpy.
scalar dcpdT(const scalar p, const scalar T) const
Derivative of cp w.r.t. temperature.
scalar Kn(const scalar p, const scalar T, const scalar n) const
Equilibrium constant [] i.t.o. number of moles.
scalar dKcdTbyKc(const scalar p, const scalar T) const
Derivative of B (acooding to Niemeyer et al.) w.r.t. temperature.
scalar THs(const scalar Hs, const scalar p, const scalar T0) const
Temperature from sensible enthalpy given an initial T0.
void operator*=(const scalar)
scalar hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kmol].
friend Ostream & operator(Ostream &, const thermo &)
scalar Kp(const scalar p, const scalar T) const
Equilibrium constant [] i.t.o. partial pressures.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for handling words, derived from Foam::string.
PtrList< volScalarField > & Y
const volScalarField & Cv
const volScalarField & Cp
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
volVectorField F(fluid.F())
const scalar RR
Universal gas constant: default in [J/(kmol K)].
const scalar Pstd
Standard pressure: default in [Pa].
dimensionedScalar exp(const dimensionedScalar &ds)
complex limit(const complex &c1, const complex &c2)
bool equal(const T &a, const T &b)
Compare two values for equality.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
psiReactionThermo & thermo