-
Notifications
You must be signed in to change notification settings - Fork 0
/
solver_eos.f90
459 lines (440 loc) · 15.4 KB
/
solver_eos.f90
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
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
MODULE solver_eos
!
!
USE def_constants, ONLY: pr
USE properties
!USE peng_robinson
USE interp_table
!USE var_const, ONLY: guessP, guessp_BC, i_bc
!
IMPLICIT NONE
PRIVATE F_Brent
PUBLIC eos_1d,BrentRoots2
CONTAINS
!===============================================================================
SUBROUTINE eos_1d(MODE, out_1, out_2, resnorm, Niter,&
& exitflag, GGG, X0, in_1, out3)
!===============================================================================
! 1D inverse the internal energy to find the temperature
!------------------------------------------------------------------------------
!
! Input :
! -------
! MODE = 2 Peng-robinson, 4 span-wagner, 5 (P,T)--> v and e
! GGG = e
! in_1 = v
! X0 = initial guess of T
!
!
! Output :
! -------
!
! out_1 = final T
! out_2 = idem
! out3 = idem
!
!-------------------------------------------------------------------
! M. De Lorenzo 03/2016
!-------------------------------------------------------------------
!
!
IMPLICIT NONE
!
INTEGER :: MAX_ITER
INTEGER, INTENT(IN) :: MODE
INTEGER, INTENT(OUT) :: exitflag, Niter
REAL(8), INTENT(IN) :: GGG, X0, in_1
REAL(8), INTENT(OUT) :: out_1, out_2, resnorm, out3
REAL(8) :: TOL_X, TOL_FUN, ALPHA, MIN_LAMBDA, MAX_LAMBDA, Jstar,&
& lambda, slope, fffold, lambda_min, fff, fff2, A_l, aa, bb, Dx_J,&
& discriminant, lambda2, lambda_OLD, XXX, F, delta, XXX_forw, J0,&
& XXX_back, F_forw, F_back, dF, TYPX, dx_star, dx, XXX_old, g, J
REAL(8), DIMENSION (2,2) :: B_l
REAL(8), DIMENSION (2,1) :: C_l, aabb
!
!
XXX = X0 ! Initial value
!
! Options
!
TOL_X = 1d-20 ! relative max step tolerance
TOL_FUN = 1d-12 ! function tolerance
MAX_ITER = 500 ! max number of iterations
ALPHA = 1d-4 ! criteria for decrease
MIN_LAMBDA = 1d-1 ! min lambda
MAX_LAMBDA = 5d-1
TYPX = abs(XXX) ! x scaling value, remove zeros
!
! ---- Function F definition
! F = GGG_input - GGG
!
call F_New_Rap1D(MODE, F, out_2, XXX, GGG, in_1, out3)
!
!
! --------- Jacobian estimation-----------------
!
Dx_J = 1d-6 ! finite difference delta
delta = Dx_J
XXX_forw = XXX + delta
XXX_back = XXX - delta
call F_New_Rap1D(MODE, F_forw, out_2, XXX_forw, GGG, in_1, out3)
call F_New_Rap1D(MODE, F_back, out_2, XXX_back, GGG, in_1, out3)
dF = F_forw - F_back ! delta F
J = dF/Dx_J/2d0 ! derivatives dF/d_n
J0 = 1d0 / TYPX ! Jacobian scaling factor
Jstar = J/J0; ! scale Jacobian
resnorm = abs(F) ! Norm of residues
exitflag = 1
!
!
!--- SOLVER -------------------
!
Niter = 0
lambda = 1d0 ! backtracking
lambda_OLD = lambda
!
!
!
1 DO WHILE ( ((resnorm>TOL_FUN) .OR. (lambda<1d0)) .AND. &
& (exitflag>=0) .AND. (Niter<=MAX_ITER))
Niter = Niter + 1
!
!--------- Newton-Raphson solver
!
IF (lambda .eq. 1.0d0) THEN
dx_star = -F/Jstar ! calculate Newton step
dx = dx_star*TYPX
g = F*Jstar
slope = g*dx_star
fffold = F*F
XXX_old = XXX
lambda_min = TOL_X/(ABS(dx)/ABS(XXX_old))
ENDIF
!
!--------- Check about proximity of XXX and XXX_OLD
!
IF (lambda < lambda_min) THEN
exitflag = 2;
!print*,'XXX is too close to XXX_OLD'
!EXIT ! OUT NewRap
ENDIF
!
!--------- Eventually, backtracking of New-Rap step
!
XXX = XXX_old + dx*lambda
call F_New_Rap1D(MODE, F, out_2, XXX, GGG, in_1, out3)
!
!--------- Jacobian estimation
!
XXX_forw = XXX + delta
XXX_back = XXX - delta
call F_New_Rap1D(MODE, F_forw, out_2, XXX_forw, GGG, in_1,&
& out3)
call F_New_Rap1D(MODE, F_back, out_2, XXX_back, GGG, in_1,&
& out3)
dF = F_forw - F_back ! delta F
J = dF/Dx_J/2d0 ! derivatives dF/d_n
Jstar = J/J0 ! scale Jacobian
fff = F*F ! f to be minimized
!
!--------- Optimization technique (minimizing fff = SUM(F*F) )
!
IF (fff > fffold + ALPHA*lambda*slope) THEN
IF (lambda .eq. 1d0) THEN
lambda = -slope/2d0/(fff-fffold-slope) ! calculate lambda
ELSE
A_l = 1d0/(lambda_OLD - lambda2)
B_l(1,1) = 1d0/lambda_OLD/lambda_OLD
B_l(1,2) = -1d0/lambda2/lambda2
B_l(2,1) = -lambda2/lambda_OLD/lambda_OLD
B_l(2,2) = lambda_OLD/lambda2/lambda2
C_l(1,1) = fff-fffold-lambda_OLD*slope
C_l(2,1) = fff2-fffold-lambda2*slope
!
aabb = A_l* MATMUL(B_l,C_l)
aa = aabb(1,1)
bb = aabb(2,1)
IF (aa .EQ. 0.0d0) THEN
lambda = -slope/2d0/bb;
ELSE
discriminant = bb**2d0 - 3d0*aa*slope;
IF (discriminant < 0d0) THEN
lambda = MAX_LAMBDA*lambda_OLD;
ELSEIF (bb <= 0d0) THEN
lambda = (-bb + sqrt(discriminant))/3d0/aa
ELSE
lambda = -slope/(bb + sqrt(discriminant))
ENDIF
ENDIF
lambda = min(lambda,MAX_LAMBDA*lambda_OLD) ! minimum step length
ENDIF
ELSE
lambda = 1d0; ! fraction of Newton step
ENDIF
IF (lambda < 1d0) THEN
lambda2 = lambda_OLD;
fff2 = fff; ! save 2nd most previous value
lambda = max(lambda,MIN_LAMBDA*lambda_OLD) ! minimum step length
GO TO 1
ENDIF
!
resnorm = abs(F)
!
END DO
!
out_1 = XXX
!
END SUBROUTINE eos_1d
!
!
!===============================================================================
SUBROUTINE F_New_Rap1D(MODE, F, out_2, XXX, GGG, in_1, out_3)
!===============================================================================
IMPLICIT NONE
INTEGER, INTENT(IN) :: MODE
REAL(8), INTENT(IN) :: in_1, XXX, GGG
REAL(8), INTENT(OUT) :: F, out_2,out_3
!
!LOCAL
INTEGER :: i,j
REAL(8) :: temp, e, v,sound,a_out,res,p_guess
REAL(8) :: delta_p,press,e_l,e_v,V_v,V_l,qual
! (v,e) --> T: v=in_1, e=GGG, T=XXX
! v = in_1
! temp = XXX
IF (MODE .EQ. 4) THEN
v = in_1
temp = XXX
!
CALL inter_energy(temp, v, e)
F = (e - GGG) * 1d-3
out_2 = 0.0d0
out_3 = 0.0d0
!
ELSEIF (MODE .EQ. 2) THEN
!
v = in_1
temp = XXX
! CALL interenergy_pr(temp, v, e)
F = (e - GGG) * 1d-3
out_2 = 0.0d0
out_3 = 0.0d0
!
ELSEIF (MODE .EQ. 5) THEN
! (T,p) --> rho,e, temp=in_1,
temp = in_1
v = XXX
CALL pressure(temp,v,press)
F = (press - GGG) * 10d-6
CALL inter_energy(temp,v,e)
out_2 = e
ELSEIF (MODE .EQ. 3) THEN
! (v,p) --> e for outletBC using LOOKUP-table
v = in_1
e = XXX
!
CALL CO2BLLT_EQUI(press,temp,sound,&
& qual, a_out, res, e,v,i)
F = (press - GGG) * 1d-6
! guessp_BC(i_bc) = press
!
out_2 = 0.0d0
out_3 = 0.0d0
!
ELSEIF (MODE .EQ. 6) THEN
!(e,p)--> v for outlet BC
e = in_1
v = XXX
!
! guessp_BC(i_bc) = (guessp_BC(i_bc)*2.0+guessP(100,i_bc))/3.0
! print*,'guessp_BC', guessp_BC(i_bc)
! CALL CO2BLLT_EQUI(press,temp,sound,&
! & qual, a_out, res, e,v,guessp_BC(i_bc))
F = (press - GGG) * 1d-6
! print*, 'press_outlet',press
! guessp_BC(i_bc) = press
!
out_2 = 0.0d0
out_3 = 0.0d0
ELSE
print*,'MODE of EOS unknown'
STOP
ENDIF
!
END SUBROUTINE F_New_Rap1D
!!
!!
!=======================================================================================
!
SUBROUTINE BrentRoots2(MODE, out_1, out_2, residue, Niter, GGG, lower, upper, in_1, in_2)
!
!=======================================================================================
!
IMPLICIT NONE
INTEGER, INTENT(IN) :: MODE
INTEGER, INTENT(OUT) :: Niter
REAL(8), INTENT(IN) :: lower, upper, GGG, in_1, in_2
REAL(8), INTENT(OUT) :: out_1, out_2, residue
INTEGER :: done
REAL(8), PARAMETER :: FPP = 1d-16 ! floating-point precision
REAL(8), PARAMETER :: MAX_ITER = 2d2 ! maximum allowed number of iterations
REAL(8), PARAMETER :: nearzero = 1d-16
REAL(8) :: TOL, AA, BB, CC, DD, EE, FA, FB, FC, xm, PP, SS, QQ, RR,&
& Tol1, Root, Minimum
TOL = 1d-12
done = 0
Niter = 0
!
AA = lower
BB = upper
CALL F_Brent(MODE, FA, out_2, AA, GGG, in_1, in_2)
IF (abs(FA) .LT. TOL) THEN
out_1 = AA
RETURN
ENDIF
!
CALL F_Brent(MODE, FB, out_2, BB, GGG, in_1, in_2)
IF (abs(FB) .LT. TOL) THEN
out_1 = BB
RETURN
ENDIF
!
! Check if root is bracketed
IF (((FA .GT. 0d0) .AND. (FB .GT. 0d0)) .OR. ((FA .LT. 0d0) .AND.&
& (FB .LT. 0d0))) THEN
print*, 'ERROR: Root is not bracketed in MODE', MODE
STOP
ELSE
FC = FB
DO WHILE ((done .EQ. 0) .AND. (Niter .LT. MAX_ITER))
Niter = Niter + 1
!
IF (((FC .GT. 0d0) .AND. (FB .GT. 0d0)) .OR.&
& ((FC .LT. 0d0) .AND. (FB .LT. 0d0))) THEN
! rename AA, BB, CC and adjust bounding interval DD
CC = AA
FC = FA
DD = BB - AA
EE = DD
ENDIF
IF (abs(FC) .LT. abs(FB)) THEN
AA = BB
BB = CC
CC = AA
FA = FB
FB = FC
FC = FA
ENDIF
Tol1 = 2d0 * FPP * abs(BB) + 0.5d0 * TOL ! convergence check
xm = 0.5d0 * (CC - BB)
IF((abs(xm) .LE.Tol1) .OR. (abs(FA) .LT.nearzero)) THEN
! a root has been found
Root = BB
done = 1
CALL F_Brent(MODE, residue, out_2,BB,GGG,in_1,in_2)
ELSE
IF ((abs(EE) .GE. Tol1) .AND.&
& (abs(FA) .GT. abs(FB))) THEN
SS = FB / FA ! attempt inverse quadratic function
IF (abs(AA - CC) .LT. nearzero) THEN
PP = 2d0 * xm * SS
QQ = 1d0 - SS
ELSE
QQ = FA / FC
RR = FB / FC
PP = SS * (2d0 * xm * QQ * (QQ - RR) -&
& (BB- AA) * (RR - 1d0))
QQ = (QQ - 1d0) * (RR - 1d0) * (SS - 1d0)
ENDIF
IF (PP .GT. nearzero) THEN
QQ = -QQ
ENDIF
PP = abs(PP)
! Evaluation of the minimum
IF ((3d0 * xm * QQ - abs(Tol1 * QQ)) .LT.&
& (abs(EE * QQ))) THEN
Minimum = 3d0 * xm * QQ - abs(Tol1 * QQ)
ELSE
Minimum = abs(EE * QQ)
ENDIF
!
IF ((2d0 * PP) .LT. Minimum) THEN
! accept interpolation
EE = DD
DD = PP / QQ
ELSE
! interpolation failed, use bisection
DD = xm
EE = DD
ENDIF
ELSE
! bounds decreasing too slowly, use bisection
DD = xm
EE = DD
ENDIF
AA = BB ! move last best guess to AA
FA = FB
IF (abs(DD) .GT. Tol1) THEN
! evaluate new trial root
BB = BB + DD
ELSE
IF (xm .GT. 0d0) THEN
BB = BB + abs(Tol1)
ELSE
BB = BB - abs(Tol1)
ENDIF
ENDIF
CALL F_Brent(MODE, FB, out_2, BB, GGG, in_1, in_2)
IF (Niter .GT. MAX_ITER) THEN
print*, 'Max N iteration exceeded in Brent'
STOP
ENDIF
ENDIF
END DO
ENDIF
!
out_1 = Root
!
END SUBROUTINE BrentRoots2
!
!
!===============================================================================
SUBROUTINE F_Brent(MODE, F, out_2, XXX, GGG, in_1, in_2)
!===============================================================================
IMPLICIT NONE
INTEGER, INTENT(IN) :: MODE
REAL(8), INTENT(IN) :: XXX, GGG, in_1, in_2
REAL(8), INTENT(OUT) :: F, out_2
REAL(8) :: temp, v_v, v_l, e_v, e_l,e,press,v,qual
REAL(8) :: a_out,sound,res
INTEGER :: flag
IF (MODE .EQ. 1) THEN
!
! In two-phase region, (e,v) --> (psat,x,Tsat)
press = XXX
v = in_1
e = GGG
CALL satprop(5, press, temp, v_v, v_l, e_v, e_l)
! print*, 'iter process',v,v_l, v_v
out_2 = temp
qual = (v - v_l)/(v_v-v_l)
F = (e - qual*e_v - (1_pr - qual)*e_l) * 1e-3_pr
!
!
ELSEIF (MODE .EQ. 2) THEN
! (v,p) --> e for outletBC using LOOKUP-table
e = XXX
v = in_1
!
CALL CO2BLLT_EQUI(press,temp,sound,&
& qual, a_out, res,e,v,flag)
! print*, 'pressure ', press
F = (GGG - press) * 1d-6
out_2 = 0.0d0
ELSE
STOP 'MODE of BrentRoots unknown'
ENDIF
!
END SUBROUTINE F_Brent
!
END MODULE solver_eos