1 |
/* |
2 |
ASCEND modelling environment |
3 |
Copyright (C) 2004-2010 John Pye |
4 |
|
5 |
This program is free software; you can redistribute it and/or |
6 |
modify it under the terms of the GNU General Public License |
7 |
as published by the Free Software Foundation; either version 2 |
8 |
of the License, or (at your option) any later version. |
9 |
|
10 |
This program is distributed in the hope that it will be useful, |
11 |
but WITHOUT ANY WARRANTY; without even the implied warranty of |
12 |
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
13 |
GNU General Public License for more details. |
14 |
|
15 |
You should have received a copy of the GNU General Public License |
16 |
along with this program. If not, see <http://www.gnu.org/licenses/>. |
17 |
*/ |
18 |
|
19 |
#include "solve_ph.h" |
20 |
#include "fprops.h" |
21 |
#include "sat.h" |
22 |
#include "derivs.h" |
23 |
#include "rundata.h" |
24 |
|
25 |
#include <stdio.h> |
26 |
#include <math.h> |
27 |
#include <stdlib.h> |
28 |
|
29 |
//#define SOLVE_PH_DEBUG |
30 |
#define SOLVE_PH_ERRORS |
31 |
|
32 |
#define FPE_DEBUG |
33 |
#define ASSERT_DEBUG |
34 |
|
35 |
#ifdef ASSERT_DEBUG |
36 |
# include <assert.h> |
37 |
#else |
38 |
# define assert(ARGS...) |
39 |
#endif |
40 |
|
41 |
#include <setjmp.h> |
42 |
#include <signal.h> |
43 |
|
44 |
#ifdef FPE_DEBUG |
45 |
#define _GNU_SOURCE |
46 |
#include <fenv.h> |
47 |
int feenableexcept(int excepts); |
48 |
int fedisableexcept(int excepts); |
49 |
int fegetexcept(void); |
50 |
#endif |
51 |
|
52 |
#ifdef SOLVE_PH_DEBUG |
53 |
# include "color.h" |
54 |
# define MSG(FMT, ...) \ |
55 |
color_on(stderr,ASC_FG_BRIGHTRED);\ |
56 |
fprintf(stderr,"%s:%d: ",__FILE__,__LINE__);\ |
57 |
color_on(stderr,ASC_FG_BRIGHTBLUE);\ |
58 |
fprintf(stderr,"%s: ",__func__);\ |
59 |
color_off(stderr);\ |
60 |
fprintf(stderr,FMT "\n",##__VA_ARGS__) |
61 |
#else |
62 |
# define MSG(ARGS...) ((void)0) |
63 |
#endif |
64 |
|
65 |
#ifdef SOLVE_PH_ERRORS |
66 |
# include "color.h" |
67 |
# define ERRMSG(STR,...) \ |
68 |
color_on(stderr,ASC_FG_BRIGHTRED);\ |
69 |
fprintf(stderr,"ERROR:");\ |
70 |
color_off(stderr);\ |
71 |
fprintf(stderr," %s:%d:" STR "\n", __func__, __LINE__ ,##__VA_ARGS__) |
72 |
#else |
73 |
# define ERRMSG(ARGS...) ((void)0) |
74 |
#endif |
75 |
|
76 |
int fprops_region_ph(double p, double h, const PureFluid *fluid, FpropsError *err){ |
77 |
double Tsat, rhof, rhog; |
78 |
double p_c = fluid->data->p_c; |
79 |
|
80 |
if(p >= p_c)return FPROPS_NON; |
81 |
|
82 |
fprops_sat_p(p, &Tsat, &rhof, &rhog, fluid, err); |
83 |
if(*err){ |
84 |
*err = FPROPS_SAT_CVGC_ERROR; |
85 |
return FPROPS_ERROR; |
86 |
} |
87 |
|
88 |
double hf = fluid->h_fn(Tsat, rhof, fluid->data, err); |
89 |
if(h <= hf)return FPROPS_NON; |
90 |
|
91 |
double hg = fluid->h_fn(Tsat,rhog, fluid->data, err); |
92 |
if(h >= hg)return FPROPS_NON; |
93 |
|
94 |
return FPROPS_SAT; |
95 |
} |
96 |
|
97 |
typedef void SignalHandler(int); |
98 |
|
99 |
#ifdef FPE_DEBUG |
100 |
jmp_buf mark; |
101 |
void fprops_fpe(int sig){ |
102 |
MSG("Catching signal %d!",sig); |
103 |
feclearexcept(FE_DIVBYZERO|FE_INVALID|FE_OVERFLOW); |
104 |
longjmp(mark, -1); |
105 |
} |
106 |
#endif |
107 |
|
108 |
void fprops_solve_ph(double p, double h, double *T, double *rho, int use_guess |
109 |
, const PureFluid *fluid, FpropsError *err |
110 |
){ |
111 |
double Tsat, rhof, rhog, hf, hg; |
112 |
double T1, rho1; |
113 |
int subcrit_pressure = 0; |
114 |
int liquid_iteration = 0; |
115 |
double rhof_t; |
116 |
double p_c = fluid->data->p_c; |
117 |
|
118 |
MSG("Solving for p=%f bar, h=%f kJ/kgK (EOS type %d)",p/1e5,h/1e3,fluid->type); |
119 |
|
120 |
#ifdef FPE_DEBUG |
121 |
feenableexcept(FE_DIVBYZERO | FE_INVALID | FE_OVERFLOW); |
122 |
SignalHandler *old = signal(SIGFPE,&fprops_fpe); |
123 |
(void)old;/* not used for anything at this stage */ |
124 |
int jmpret = setjmp(mark); |
125 |
if(jmpret==0){ |
126 |
#endif |
127 |
MSG("p_c = %f bar",p_c/1e5); |
128 |
if(p < p_c){ |
129 |
MSG("Calculate saturation Tsat(p < p_c) with p = %f",p); |
130 |
fprops_sat_p(p, &Tsat, &rhof, &rhog, fluid, err); |
131 |
if(*err){ |
132 |
ERRMSG("Unable to solve saturation state"); |
133 |
*err = FPROPS_SAT_CVGC_ERROR; |
134 |
return; |
135 |
} |
136 |
hf = fluid->h_fn(Tsat, rhof, fluid->data, err); |
137 |
hg = fluid->h_fn(Tsat, rhog, fluid->data, err); |
138 |
MSG("at p = %f bar, T_sat = %f, rhof = %f, hf = %f kJ/kg, hg = %f",p/1e5,Tsat,rhof, hf/1e3,hg/1e3); |
139 |
|
140 |
if(hf < h && h < hg){ |
141 |
MSG("SATURATION REGION"); |
142 |
/* saturation region... easy */ |
143 |
double x = (h - hf)/(hg - hf); |
144 |
*rho = 1./(x/rhog + (1.-x)/rhof); |
145 |
*T = Tsat; |
146 |
return; |
147 |
} |
148 |
|
149 |
subcrit_pressure = 1; |
150 |
if(h < hf){ |
151 |
liquid_iteration = 1; |
152 |
if(!use_guess){ |
153 |
*T = Tsat; |
154 |
*rho = rhof; |
155 |
MSG("LIQUID GUESS: T = %f, rho = %f",*T, *rho); |
156 |
} |
157 |
}else if(!use_guess){ |
158 |
*T = 1.1 * Tsat; |
159 |
*rho = rhog * 0.5; |
160 |
MSG("GAS GUESS: T = %f, rho = %f",*T, *rho); |
161 |
} |
162 |
}else{ /* p >= p_c */ |
163 |
/* FIXME: still some problems here at very high pressures */ |
164 |
if(!use_guess){ |
165 |
/* FIXME we should cache/precalculate hc, store in rundata. */ |
166 |
double hc = fluid->h_fn(fluid->data->T_c, fluid->data->rho_c, fluid->data, err); |
167 |
assert(!isnan(hc)); |
168 |
MSG("hc = %f kJ/kgK",hc/1e3); |
169 |
if(h < 0.8*hc){ /* FIXME should make use of h_ref here in some way? Otherwise arbitrary... */ |
170 |
#if 0 |
171 |
double Tsat1, psat1, rhof1, rhog1; |
172 |
int res = fprops_sat_p(p, &Tsat1, &rhof1, &rhog1, fluid); |
173 |
MSG("Unable to solve saturation at p = %f",p); |
174 |
*T = Tsat1; |
175 |
*rho = rhof1; |
176 |
#else |
177 |
MSG("h < 0.9 hc... using saturation Tsat(hf) for starting guess"); |
178 |
double Tsat1, psat1, rhof1, rhog1; |
179 |
fprops_sat_hf(h, &Tsat1, &psat1, &rhof1, &rhog1, fluid, err); |
180 |
if(*err){ |
181 |
MSG("Unable to solve Tsat(hf)"); |
182 |
/* accuracy of the estimate of Tsat(hf) doesn't matter |
183 |
very much so we can ignore this error. */ |
184 |
|
185 |
/* try initialising to T_t, rho_f(T_t) */ |
186 |
Tsat1 = fluid->data->T_t; |
187 |
MSG("T_t = %f",Tsat1); |
188 |
fprops_sat_T(Tsat1, &psat1, &rhof1, &rhog1, fluid, err); |
189 |
if(*err)MSG("Unable to solve rhof(Tt)"); |
190 |
} |
191 |
*T = Tsat1; |
192 |
*rho = rhof1; |
193 |
#endif |
194 |
}else{ |
195 |
*T = fluid->data->T_c * 1.01; |
196 |
*rho = fluid->data->rho_c * 1.05; |
197 |
} |
198 |
MSG("SUPERCRITICAL GUESS: T = %f, rho = %f", *T, *rho); |
199 |
} |
200 |
} |
201 |
|
202 |
MSG("STARTING NON-SAT ITERATION"); |
203 |
//*rho = 976.82687191126922; |
204 |
//*T = 344.80371310850518; |
205 |
|
206 |
T1 = *T; |
207 |
rho1 = *rho; |
208 |
assert(!isnan(T1)); |
209 |
assert(!isnan(rho1)); |
210 |
|
211 |
if(liquid_iteration){ |
212 |
double pt,rhogt; |
213 |
fprops_triple_point(&pt, &rhof_t, &rhogt, fluid, err); |
214 |
if(*err){ |
215 |
ERRMSG("Unable to solve triple point liquid density."); |
216 |
*err = FPROPS_SAT_CVGC_ERROR; |
217 |
return; |
218 |
} |
219 |
} |
220 |
|
221 |
/* try our own home-baked newton iteration */ |
222 |
int i = 0; |
223 |
*err = FPROPS_NO_ERROR; |
224 |
double delta_T = 0; |
225 |
double delta_rho = 0; |
226 |
MSG("STARTING ITERATION"); |
227 |
MSG("rhof_t = %f",rhof_t); |
228 |
while(i++ < 200){ |
229 |
double p1 = fluid->p_fn(T1,rho1, fluid->data, err); |
230 |
if(*err){ |
231 |
MSG("Got an error ('%s') in p_fn calculation",fprops_error(*err)); |
232 |
} |
233 |
double h1 = fluid->h_fn(T1,rho1, fluid->data, err); |
234 |
assert(!isnan(h1)); |
235 |
|
236 |
MSG(" %d: T = %f, rho = %f\tp = %f bar, h = %f kJ/kg", i, T1, rho1, p1/1e5, h1/1e3); |
237 |
//MSG(" p error = %f bar",(p1 - p)/1e5); |
238 |
//MSG(" h error = %f kJ/kg",(h1 - h)/1e3); |
239 |
|
240 |
if(p1 < 0){ |
241 |
MSG("p1 < 0, reducing dT, drho"); |
242 |
T1 -= (delta_T *= 0.5); |
243 |
rho1 -= (delta_rho *= 0.5); |
244 |
continue; |
245 |
} |
246 |
|
247 |
if(fabs(p1 - p) < 1e-4 && fabs(h1 - h) < 1e-8){ |
248 |
MSG("Converged to T = %f, rho = %f, in homebaked Newton solver", T1, rho1); |
249 |
*T = T1; |
250 |
*rho = rho1; |
251 |
return; |
252 |
} |
253 |
/* calculate step, we're solving log(p1) in this code... */ |
254 |
double f = log(p1) - log(p); |
255 |
double g = h1 - h; |
256 |
assert(!isnan(f)); |
257 |
assert(!isnan(g)); |
258 |
assert(rho1 != 0); |
259 |
assert(p1 != 0); |
260 |
assert(!isnan(fprops_non_dZdT_v('p', T1,rho1, fluid,err))); |
261 |
double f_T = 1./p1 * fprops_non_dZdT_v('p', T1,rho1, fluid,err); |
262 |
double f_rho = -1./p1/SQ(rho1) * fprops_non_dZdv_T('p', T1, rho1, fluid,err); |
263 |
double g_T = fprops_non_dZdT_v('h', T1,rho1, fluid,err); |
264 |
double g_rho = -1./SQ(rho1) * fprops_non_dZdv_T('h', T1, rho1, fluid,err); |
265 |
assert(!isnan(f_T)); |
266 |
|
267 |
if(isnan(f_rho)){ |
268 |
ERRMSG(" rho1 = %f, T1 = %f",rho1, T1); |
269 |
} |
270 |
assert(!isnan(f_rho)); |
271 |
|
272 |
assert(!isnan(g_T)); |
273 |
assert(!isnan(g_rho)); |
274 |
|
275 |
double det = g_rho * f_T - f_rho * g_T; |
276 |
assert(det!=0); |
277 |
assert(!isnan(det)); |
278 |
MSG(" df/dT = %e\t\tdf/drho = %e",f_T, f_rho); |
279 |
MSG(" dg/dT = %e\t\tdg/drho = %e",g_T, g_rho); |
280 |
|
281 |
delta_T = -1./det * (g_rho * f - f_rho * g); |
282 |
delta_rho = -1./det * (f_T * g - g_T * f); |
283 |
assert(!isnan(delta_T)); |
284 |
assert(!isnan(delta_rho)); |
285 |
MSG(" dT = %f", delta_T); |
286 |
MSG(" drho = %f", delta_rho); |
287 |
|
288 |
if(subcrit_pressure){ |
289 |
if(h > hg){ |
290 |
/* vapour */ |
291 |
int i = 0; |
292 |
while(rho1 + delta_rho > rhog && i++ < 20){ |
293 |
delta_rho *= 0.5; |
294 |
} |
295 |
}else{ |
296 |
/* liquid */ |
297 |
while(rho1 + delta_rho < rhof && i++ < 20){ |
298 |
delta_rho *= 0.5; |
299 |
} |
300 |
if(T1 + delta_T < fluid->data->T_t) delta_T = 0.5 * (T1 + fluid->data->T_t); |
301 |
} |
302 |
}else{ |
303 |
#if 0 |
304 |
/* supercritical... stay above critical temperature of density > rho_crit */ |
305 |
if(rho1 + delta_rho > fluid->data->rho_c && T1 + delta_T < fluid->data->T_c){ |
306 |
delta_T = 0.5 *(T1 + fluid->data->T_c); |
307 |
} |
308 |
#endif |
309 |
} |
310 |
/* don't go too dense */ |
311 |
#if 1 |
312 |
if(liquid_iteration){ |
313 |
MSG("rho1 = %f, delta_rho = %f, rho1+delta_rho = %f", rho1, delta_rho, rho1+delta_rho); |
314 |
if(rho1 + delta_rho > rhof_t){ |
315 |
MSG("Limit rho to be less than rhof_t"); |
316 |
delta_rho = rhof_t - rho1; |
317 |
} |
318 |
} |
319 |
#endif |
320 |
|
321 |
/* don't go too hot */ |
322 |
if(T1 + delta_T > 5000) delta_T = 5000 - T1; |
323 |
|
324 |
/* avoid huge step */ |
325 |
while(fabs(delta_T / T1) > 0.6){ |
326 |
MSG("Reduce delta_T"); |
327 |
delta_T *= 0.5; |
328 |
} |
329 |
|
330 |
/* NOTE if step limit is less than 0.5, we're getting errors */ |
331 |
while(fabs(delta_rho / rho1) > 0.7){ |
332 |
MSG("Reduce delta_rho"); |
333 |
delta_rho *= 0.5; |
334 |
} |
335 |
|
336 |
T1 = T1 + delta_T; |
337 |
if(T1 < fluid->data->T_t)T1 = fluid->data->T_t; |
338 |
rho1 = rho1 + delta_rho; |
339 |
} |
340 |
#ifdef FPE_DEBUG |
341 |
}else{ |
342 |
/* an FPE occurred */ |
343 |
MSG("An FPE occurred"); |
344 |
abort(); |
345 |
} |
346 |
#endif |
347 |
|
348 |
*T = T1; |
349 |
*rho = rho1; |
350 |
ERRMSG("Iteration failed for '%s' with p = %.12e, h = %.12e",fluid->name, p,h); |
351 |
*err = FPROPS_NUMERIC_ERROR; |
352 |
return; |
353 |
|
354 |
#if 0 |
355 |
int res = fprops_nonsolver('p','h',p,h,T,rho,fluid); |
356 |
ERRMSG("Iteration failed in nonsolver"); |
357 |
return res; |
358 |
#endif |
359 |
} |
360 |
|