30 template<
class EquationOfState,
int PolySize>
33 const EquationOfState& pt,
36 const Polynomial<PolySize>& CvCoeffs,
37 const typename Polynomial<PolySize>::intPolyType& eCoeffs,
38 const Polynomial<PolySize>& sCoeffs
50 template<
class EquationOfState,
int PolySize>
57 EquationOfState(
name, pt),
60 CvCoeffs_(pt.CvCoeffs_),
61 eCoeffs_(pt.eCoeffs_),
68 template<
class EquationOfState,
int PolySize>
78 template<
class EquationOfState,
int PolySize>
89 template<
class EquationOfState,
int PolySize>
100 template<
class EquationOfState,
int PolySize>
107 return es(
p,
T) + hf();
111 template<
class EquationOfState,
int PolySize>
119 template<
class EquationOfState,
int PolySize>
126 return sCoeffs_.value(
T) + EquationOfState::sv(
p,
T);
130 template<
class EquationOfState,
int PolySize>
142 template<
class EquationOfState,
int PolySize>
150 return CvCoeffs_.derivative(
T);
156 template<
class EquationOfState,
int PolySize>
162 scalar Y1 = this->
Y();
164 EquationOfState::operator+=(pt);
166 if (
mag(this->
Y()) > small)
169 const scalar Y2 = pt.Y()/this->
Y();
171 hf_ = Y1*hf_ + Y2*pt.hf_;
172 sf_ = Y1*sf_ + Y2*pt.sf_;
173 CvCoeffs_ = Y1*CvCoeffs_ + Y2*pt.CvCoeffs_;
174 eCoeffs_ = Y1*eCoeffs_ + Y2*pt.eCoeffs_;
175 sCoeffs_ = Y1*sCoeffs_ + Y2*pt.sCoeffs_;
180 template<
class EquationOfState,
int PolySize>
186 EquationOfState::operator*=(
s);
192 template<
class EquationOfState,
int PolySize>
199 EquationOfState eofs = pt1;
202 if (
mag(eofs.Y()) < small)
215 const scalar Y1 = pt1.Y()/eofs.Y();
216 const scalar Y2 = pt2.Y()/eofs.Y();
218 return ePolynomialThermo<EquationOfState, PolySize>
221 Y1*pt1.hf_ + Y2*pt2.hf_,
222 Y1*pt1.sf_ + Y2*pt2.sf_,
223 Y1*pt1.CvCoeffs_ + Y2*pt2.CvCoeffs_,
224 Y1*pt1.eCoeffs_ + Y2*pt2.eCoeffs_,
225 Y1*pt1.sCoeffs_ + Y2*pt2.sCoeffs_
231 template<
class EquationOfState,
int PolySize>
235 const ePolynomialThermo<EquationOfState, PolySize>& pt
238 return ePolynomialThermo<EquationOfState, PolySize>
240 s*
static_cast<const EquationOfState&
>(pt),
250 template<
class EquationOfState,
int PolySize>
253 const ePolynomialThermo<EquationOfState, PolySize>& pt1,
254 const ePolynomialThermo<EquationOfState, PolySize>& pt2
259 static_cast<const EquationOfState&
>(pt1)
260 ==
static_cast<const EquationOfState&
>(pt2)
263 const scalar Y1 = pt1.Y()/eofs.Y();
264 const scalar Y2 = pt2.Y()/eofs.Y();
266 return ePolynomialThermo<EquationOfState, PolySize>
269 Y2*pt2.hf_ - Y1*pt1.hf_,
270 Y2*pt2.sf_ - Y1*pt1.sf_,
271 Y2*pt2.CvCoeffs_ - Y1*pt1.CvCoeffs_,
272 Y2*pt2.eCoeffs_ - Y1*pt1.eCoeffs_,
273 Y2*pt2.sCoeffs_ - Y1*pt1.sCoeffs_
scalar es(const scalar p, const scalar T) const
scalar Cv(const scalar p, const scalar T) const
Internal energy based thermodynamics package using a polynomial function of temperature for the const...
scalar gStd(const scalar T) const
Gibbs free energy of the mixture in the standard state [J/kg].
scalar s(const scalar p, const scalar T) const
Entropy [J/kg/K].
scalar dCpdT(const scalar p, const scalar T) const
Temperature derivative of heat capacity at constant pressure.
scalar ea(const scalar p, const scalar T) const
Absolute internal energy [J/kg].
scalar Cv(const scalar p, const scalar T) const
Heat capacity at constant volume [J/kg/K].
scalar hf() const
Enthalpy of formation [J/kg].
scalar limit(const scalar) const
Limit the temperature to be in the range Tlow_ to Thigh_.
scalar es(const scalar p, const scalar T) const
Sensible internal energy [J/kg].
A class for handling words, derived from string.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
volScalarField sf(fieldObject, mesh)
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.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const dimensionedScalar e
Elementary charge.
const dimensionedScalar Pstd
Standard pressure.
word name(const bool)
Return a word representation of a bool.
dimensioned< scalar > mag(const dimensioned< Type > &)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
PtrList< volScalarField > & Y