Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / thermophysicalModels / specie / thermo / hPolynomial / hPolynomialThermoI.H
blobb41de9b877f04ddb5caeae7ebc2f9e4bc25b5c3c
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2008-2011 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software: you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by
13     the Free Software Foundation, either version 3 of the License, or
14     (at your option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "hPolynomialThermo.H"
28 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
30 template<class EquationOfState, int PolySize>
31 inline Foam::hPolynomialThermo<EquationOfState, PolySize>::hPolynomialThermo
33     const EquationOfState& pt,
34     const scalar Hf,
35     const scalar Sf,
36     const Polynomial<PolySize>& CpCoeffs,
37     const typename Polynomial<PolySize>::intPolyType& hCoeffs,
38     const Polynomial<PolySize>& sCoeffs
41     EquationOfState(pt),
42     Hf_(Hf),
43     Sf_(Sf),
44     CpCoeffs_(CpCoeffs),
45     hCoeffs_(hCoeffs),
46     sCoeffs_(sCoeffs)
50 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
52 template<class EquationOfState, int PolySize>
53 inline Foam::hPolynomialThermo<EquationOfState, PolySize>::hPolynomialThermo
55     const hPolynomialThermo& pt
58     EquationOfState(pt),
59     Hf_(pt.Hf_),
60     Sf_(pt.Sf_),
61     CpCoeffs_(pt.CpCoeffs_),
62     hCoeffs_(pt.hCoeffs_),
63     sCoeffs_(pt.sCoeffs_)
67 template<class EquationOfState, int PolySize>
68 inline Foam::hPolynomialThermo<EquationOfState, PolySize>::hPolynomialThermo
70     const word& name,
71     const hPolynomialThermo& pt
74     EquationOfState(name, pt),
75     Hf_(pt.Hf_),
76     Sf_(pt.Sf_),
77     CpCoeffs_(pt.CpCoeffs_),
78     hCoeffs_(pt.hCoeffs_),
79     sCoeffs_(pt.sCoeffs_)
83 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
85 template<class EquationOfState, int PolySize>
86 inline Foam::scalar Foam::hPolynomialThermo<EquationOfState, PolySize>::cp
88     const scalar T
89 ) const
91     return CpCoeffs_.value(T);
95 template<class EquationOfState, int PolySize>
96 inline Foam::scalar Foam::hPolynomialThermo<EquationOfState, PolySize>::h
98     const scalar T
99 ) const
101     return hCoeffs_.value(T);
105 template<class EquationOfState, int PolySize>
106 inline Foam::scalar Foam::hPolynomialThermo<EquationOfState, PolySize>::hs
108     const scalar T
109 ) const
111     return h(T) - hc();
115 template<class EquationOfState, int PolySize>
116 inline Foam::scalar Foam::hPolynomialThermo<EquationOfState, PolySize>::hc()
117 const
119     return Hf_;
123 template<class EquationOfState, int PolySize>
124 inline Foam::scalar Foam::hPolynomialThermo<EquationOfState, PolySize>::s
126     const scalar T
127 ) const
129     return sCoeffs_.value(T);
133 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
135 template<class EquationOfState, int PolySize>
136 inline Foam::hPolynomialThermo<EquationOfState, PolySize>&
137 Foam::hPolynomialThermo<EquationOfState, PolySize>::operator=
139     const hPolynomialThermo<EquationOfState, PolySize>& pt
142     EquationOfState::operator=(pt);
144     Hf_ = pt.Hf_;
145     Sf_ = pt.Sf_;
146     CpCoeffs_ = pt.CpCoeffs_;
147     hCoeffs_ = pt.hCoeffs_;
148     sCoeffs_ = pt.sCoeffs_;
150     return *this;
154 template<class EquationOfState, int PolySize>
155 inline void Foam::hPolynomialThermo<EquationOfState, PolySize>::operator+=
157     const hPolynomialThermo<EquationOfState, PolySize>& pt
160     scalar molr1 = this->nMoles();
162     EquationOfState::operator+=(pt);
164     molr1 /= this->nMoles();
165     scalar molr2 = pt.nMoles()/this->nMoles();
167     Hf_ = molr1*Hf_ + molr2*pt.Hf_;
168     Sf_ = molr1*Sf_ + molr2*pt.Sf_;
169     CpCoeffs_ = molr1*CpCoeffs_ + molr2*pt.CpCoeffs_;
170     hCoeffs_ = molr1*hCoeffs_ + molr2*pt.hCoeffs_;
171     sCoeffs_ = molr1*sCoeffs_ + molr2*pt.sCoeffs_;
175 template<class EquationOfState, int PolySize>
176 inline void Foam::hPolynomialThermo<EquationOfState, PolySize>::operator-=
178     const hPolynomialThermo<EquationOfState, PolySize>& pt
181     scalar molr1 = this->nMoles();
183     EquationOfState::operator-=(pt);
185     molr1 /= this->nMoles();
186     scalar molr2 = pt.nMoles()/this->nMoles();
188     Hf_ = molr1*Hf_ - molr2*pt.Hf_;
189     Sf_ = molr1*Sf_ - molr2*pt.Sf_;
190     CpCoeffs_ = molr1*CpCoeffs_ - molr2*pt.CpCoeffs_;
191     hCoeffs_ = molr1*hCoeffs_ - molr2*pt.hCoeffs_;
192     sCoeffs_ = molr1*sCoeffs_ - molr2*pt.sCoeffs_;
196 template<class EquationOfState, int PolySize>
197 inline void Foam::hPolynomialThermo<EquationOfState, PolySize>::operator*=
199     const scalar s
202     EquationOfState::operator*=(s);
206 // * * * * * * * * * * * * * * * Friend Operators  * * * * * * * * * * * * * //
208 template<class EquationOfState, int PolySize>
209 inline Foam::hPolynomialThermo<EquationOfState, PolySize> Foam::operator+
211     const hPolynomialThermo<EquationOfState, PolySize>& pt1,
212     const hPolynomialThermo<EquationOfState, PolySize>& pt2
215     EquationOfState eofs = pt1;
216     eofs += pt2;
218     scalar molr1 = pt1.nMoles()/eofs.nMoles();
219     scalar molr2 = pt2.nMoles()/eofs.nMoles();
221     return hPolynomialThermo<EquationOfState, PolySize>
222     (
223         eofs,
224         molr1*pt1.Hf_ + molr2*pt2.Hf_,
225         molr1*pt1.Sf_ + molr2*pt2.Sf_,
226         molr1*pt1.CpCoeffs_ + molr2*pt2.CpCoeffs_,
227         molr1*pt1.hCoeffs_ + molr2*pt2.hCoeffs_,
228         molr1*pt1.sCoeffs_ + molr2*pt2.sCoeffs_
229     );
233 template<class EquationOfState, int PolySize>
234 inline Foam::hPolynomialThermo<EquationOfState, PolySize> Foam::operator-
236     const hPolynomialThermo<EquationOfState, PolySize>& pt1,
237     const hPolynomialThermo<EquationOfState, PolySize>& pt2
240     EquationOfState eofs = pt1;
241     eofs -= pt2;
243     scalar molr1 = pt1.nMoles()/eofs.nMoles();
244     scalar molr2 = pt2.nMoles()/eofs.nMoles();
246     return hPolynomialThermo<EquationOfState, PolySize>
247     (
248         eofs,
249         molr1*pt1.Hf_ - molr2*pt2.Hf_,
250         molr1*pt1.Sf_ - molr2*pt2.Sf_,
251         molr1*pt1.CpCoeffs_ - molr2*pt2.CpCoeffs_,
252         molr1*pt1.hCoeffs_ - molr2*pt2.hCoeffs_,
253         molr1*pt1.sCoeffs_ - molr2*pt2.sCoeffs_
254     );
258 template<class EquationOfState, int PolySize>
259 inline Foam::hPolynomialThermo<EquationOfState, PolySize> Foam::operator*
261     const scalar s,
262     const hPolynomialThermo<EquationOfState, PolySize>& pt
265     return hPolynomialThermo<EquationOfState, PolySize>
266     (
267         s*static_cast<const EquationOfState&>(pt),
268         pt.Hf_,
269         pt.Sf_,
270         pt.CpCoeffs_,
271         pt.hCoeffs_,
272         pt.sCoeffs_
273     );
277 template<class EquationOfState, int PolySize>
278 inline Foam::hPolynomialThermo<EquationOfState, PolySize> Foam::operator==
280     const hPolynomialThermo<EquationOfState, PolySize>& pt1,
281     const hPolynomialThermo<EquationOfState, PolySize>& pt2
284     return pt2 - pt1;
288 // ************************************************************************* //