[05910c]: src / thermophysicalModels / specie / equationOfState / aungierRedlichKwong / aungierRedlichKwongI.H  Maximize  Restore  History

Download this file

372 lines (301 with data), 10.4 kB

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright held by original author
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software; you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation; either version 2 of the License, or (at your
option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM; if not, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
Author
Christian Lucas
Institut f��r Thermodynamik
Technische Universit��t Braunschweig
Germany
\*---------------------------------------------------------------------------*/
#include "aungierRedlichKwong.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
// Construct from components
inline aungierRedlichKwong::aungierRedlichKwong
(
const specie& sp,
scalar pcrit,
scalar Tcrit,
scalar azentricFactor,
scalar rhocrit
)
:
specie(sp),
pcrit_(pcrit),
Tcrit_(Tcrit),
azentricFactor_(azentricFactor),
rhocrit_(rhocrit),
a_(0.42747*pow(this->RR,2)*pow(Tcrit_,2)/pcrit_),
b_(0.08664*this->RR*Tcrit_/pcrit_),
c_(this->RR*Tcrit_/(pcrit_+(a_/(this->W()/rhocrit_*(this->W()/rhocrit_+b_))))+b_-this->W()/rhocrit_),
n_(0.4986+1.2735*azentricFactor_+0.4754*pow(azentricFactor_,2)),
// Starting GUESS for the density by ideal gas law
rhostd_(this->rho(Pstd,Tstd,Pstd/(Tstd*this->R())))
{}
// Construct as named copy
inline aungierRedlichKwong::aungierRedlichKwong(const word& name, const aungierRedlichKwong& pg)
:
specie(name, pg),
pcrit_(pg.pcrit_),
Tcrit_(pg.Tcrit_),
azentricFactor_(pg.azentricFactor_),
rhocrit_(pg.rhocrit_),
a_(pg.a_),
b_(pg.b_),
c_(pg.c_),
n_(pg.n_),
rhostd_(pg.rhostd_)
{}
// Construct and return a clone
inline autoPtr<aungierRedlichKwong> aungierRedlichKwong::clone() const
{
return autoPtr<aungierRedlichKwong>(new aungierRedlichKwong(*this));
}
// Selector from Istream
inline autoPtr<aungierRedlichKwong> aungierRedlichKwong::New(Istream& is)
{
return autoPtr<aungierRedlichKwong>(new aungierRedlichKwong(is));
}
// * * * * * * * * * * * * * Member Functions * * * * * * * * * * * //
inline scalar aungierRedlichKwong::rhostd()const
{
return rhostd_;
}
//returns the pressure for a given density and temperature
inline scalar aungierRedlichKwong::p(const scalar rho,const scalar T) const
{
scalar Vm = this->W()/rho;
return this->RR*T/(Vm-b_+c_)
-a_*pow(T,-n_)/(pow(Tcrit_,-n_)*Vm*(Vm+b_));
}
//Real deviative dp/dv at constant temperature
//(molar values)
inline scalar aungierRedlichKwong::dpdv(const scalar rho, const scalar T) const
{
scalar Vm = this->W()/rho;
return a_*pow(T,-n_)*pow(Tcrit_,n_)*(b_+2*Vm)/(pow(Vm,2)*pow((b_+Vm),2))
-this->RR*T/(pow((b_-Vm-c_),2));
}
//Real deviative dp/dT at constant molar volume
//(molar values)
inline scalar aungierRedlichKwong::dpdT(const scalar rho, const scalar T) const
{
scalar Vm = this->W()/rho;
return a_*n_*pow(T,-n_-1)*pow(Tcrit_,n_)/(Vm*(Vm+b_))
-this->RR/(b_-Vm-c_);
}
//Real deviative dv/dT at constant pressure
//using implicit differentiation
//(molar values)
inline scalar aungierRedlichKwong::dvdT(const scalar rho,const scalar T) const
{
return (-1)*this->dpdT(rho,T)/this->dpdv(rho,T);
}
//Real deviative dv/dp at constant temperature
//(molar values)
inline scalar aungierRedlichKwong::dvdp(const scalar rho,const scalar T) const
{
return 1/this->dpdv(rho,T);
}
//needed to calculate the internal energy
//(molar values)
inline scalar aungierRedlichKwong::integral_p_dv
(
const scalar rho,
const scalar T
) const
{
scalar Vm = this->W()/rho;
return -pow((T/Tcrit_),-n_)*(a_*log(Vm)/(b_)
-a_*log(b_+Vm)/b_)+this->RR*T*log(Vm-b_+c_);
}
//needed to calculate the entropy
//(molar values)
inline scalar aungierRedlichKwong::integral_dpdT_dv
(
const scalar rho,
const scalar T
) const
{
scalar Vm = this->W()/rho;
return pow(T,-n_-1)*pow(Tcrit_,n_)*(a_*n_*log(Vm)/(b_)
-a_*n_*log(b_+Vm)/(b_))+ this->RR*log(-b_+c_+Vm);
}
//(molar values)
inline scalar aungierRedlichKwong::d2pdT2(const scalar rho,const scalar T) const
{
scalar Vm = this->W()/rho;
return -a_*n_*pow(T,-n_-2)*pow(Tcrit_,n_)*(n_+1)/(Vm*(b_+Vm));
}
//(molar values)
inline scalar aungierRedlichKwong::d2pdv2(const scalar rho,const scalar T) const
{
scalar Vm = this->W()/rho;
return -2*a_*pow(T,-n_)*pow(Tcrit_,n_)*(pow(b_,2)
+3*b_*Vm+3*pow(Vm,2))/(pow(Vm,3)*pow((b_+Vm),3))
-2*this->RR*T/(pow((b_-Vm-c_),3));
}
//(molar values)
//using second order implicit differentiation
inline scalar aungierRedlichKwong::d2vdT2(const scalar rho, const scalar T) const
{
return
-(
pow(this->dpdT(rho,T),2)*this->d2pdv2(rho,T)
+ pow(this->dpdv(rho,T),2)*this->d2pdT2(rho,T)
- 2*this->dpdv(rho,T)*this->dpdT(rho,T)*this->d2pdvdT(rho,T)
)
/(pow(this->dpdv(rho,T),3));
}
//(molar values)
inline scalar aungierRedlichKwong::d2pdvdT(const scalar rho, const scalar T) const
{
scalar Vm = this->W()/rho;
return -a_*n_*pow(T,-n_-1)*pow(Tcrit_,n_)*(b_+2*Vm)/(pow(Vm,2)*pow((b_+Vm),2))
-this->RR/(pow((b_-c_-Vm),2));
}
// the result of this intergal is needed for the nasa based cp polynomial
//(molar values)
inline scalar aungierRedlichKwong::integral_d2pdT2_dv
(
const scalar rho,
const scalar T
) const
{
scalar Vm = this->W()/rho;
return pow(T,-n_-2)*pow(Tcrit_,n_)*(a_*n_*(n_+1)*log(b_+Vm)
/b_-a_*n_*(1+n_)*log(Vm)/b_);
}
//Isobar expansion Coefficent beta = 1/v (dv/dt) at constant p
//(molar values)
inline scalar aungierRedlichKwong::isobarExpCoef(const scalar rho,const scalar T) const
{
return this->dvdT(rho, T)*rho/this->W();
}
//isothemal compressiblity kappa (not Thermal conductivity)
//(molar values)
inline scalar aungierRedlichKwong::isothermalCompressiblity(const scalar rho,const scalar T) const
{
return this->isobarExpCoef(rho, T)/this->dpdT(rho, T);
//also possible : return -this->dvdp(rho,T)*rho/this->W();
}
//- Return density [kg/m^3]
inline scalar aungierRedlichKwong::rho(
const scalar p,
const scalar T,
const scalar rho0
) const
{
scalar molarVolumePrevIteration;
scalar molarVolume;
label iter=0;
label maxIter_=400;
scalar tol_=1e-8;
scalar rho1=rhoMax_;
scalar rho2=rhoMin_;
molarVolume=this->W()/rho0;
do
{
molarVolumePrevIteration= molarVolume;
label i=0;
do
{
molarVolume=molarVolumePrevIteration
-(
(this->p((this->W()/molarVolumePrevIteration),T) - p)
/(this->dpdv((this->W()/molarVolumePrevIteration),T))
)/pow(2,i);
i++;
if (i>8)
{
//CL: using bisection methode as backup,
//CL: solution must be between rho=0.001 to rho=1500;
//CL: if not, change rhoMax_ and rhoMin_
for(i=0; i<200; i++)
{
scalar f1 = this->p(rho1,T) - p;
scalar f2 = this->p(rho2,T) - p;
scalar rho3 = (rho1 + rho2)/2;
scalar f3 = this->p(rho3,T) - p;
if ((f2 < 0 && f3 > 0) || (f2 > 0 && f3 < 0))
{
rho1=rho3;
}
else if ((f1 < 0 && f3 > 0)||(f1 > 0 && f3 < 0))
{
rho2=rho3;
}
else
{
rho2=(rho2 + rho3)/2;
}
if(mag(f3) < p*tol_)
{
molarVolume=this->W()/rho3;
molarVolumePrevIteration=this->W()/rho3;
break;
}
else
{
molarVolumePrevIteration=this->W()/rho3;
}
}
}
}
while
(
mag(this->p((this->W()/molarVolume),T) - p)
> mag(this->p((this->W()/molarVolumePrevIteration),T) - p)
);
if (iter++ > maxIter_)
{
FatalErrorIn
(
"inline scalar redlichKwong::rho(const scalar p, const scalar T, const scalar rho0) const "
) << "Maximum number of iterations exceeded"
<< abort(FatalError);
}
}
while(mag(molarVolumePrevIteration-molarVolume) > tol_*(this->W()/rho0));
return this->W()/molarVolume;
}
//- Return density [kg/m^3]
inline scalar aungierRedlichKwong::rho(const scalar p,const scalar T) const
{
//using perfect gas equation as starting point
return rho(p,T,p/(this->R()*T));
}
//- Return compressibility drho/dp at T=constant [s^2/m^2]
inline scalar aungierRedlichKwong::psi(const scalar rho, const scalar T) const
{
return -this->dvdp(rho,T)*pow(rho,2)/this->W();
}
//- Return compression factor []
inline scalar aungierRedlichKwong::Z( const scalar p, const scalar T,const scalar rho0) const
{
return p/(this->R()*T*this->rho(p,T,rho0));
}
// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// ************************************************************************* //