Actual source code: shashi.F90
petsc-3.8.3 2017-12-09
3: program main
4: #include <petsc/finclude/petscs.h>
5: use petsc
6: implicit none
8: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: ! Variable declarations
10: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
11: !
12: ! Variables:
13: ! snes - nonlinear solver
14: ! ksp - linear solver
15: ! pc - preconditioner context
16: ! ksp - Krylov subspace method context
17: ! x, r - solution, residual vectors
18: ! J - Jacobian matrix
19: ! its - iterations for convergence
20: !
21: SNES snes
22: PC pc
23: KSP ksp
24: Vec x,r,lb,ub
25: Mat J
26: SNESLineSearch linesearch
27: PetscErrorCode ierr
28: PetscInt its,i2,i20
29: PetscMPIInt size,rank
30: PetscScalar pfive
31: PetscReal tol
32: PetscBool setls
33: PetscScalar,pointer :: xx(:)
34: PetscScalar zero,big
35: SNESLineSearch ls
36:
37: ! Note: Any user-defined Fortran routines (such as FormJacobian)
38: ! MUST be declared as external.
40: external FormFunction, FormJacobian
41: external ShashiPostCheck
43: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
44: ! Macro definitions
45: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
46: !
47: ! Macros to make clearer the process of setting values in vectors and
48: ! getting values from vectors. These vectors are used in the routines
49: ! FormFunction() and FormJacobian().
50: ! - The element lx_a(ib) is element ib in the vector x
51: !
52: #define lx_a(ib) lx_v(lx_i + (ib))
53: #define lf_a(ib) lf_v(lf_i + (ib))
54: !
55: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
56: ! Beginning of program
57: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
59: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
60: if (ierr .ne. 0) then
61: print*,'Unable to initialize PETSc'
62: stop
63: endif
64: call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
65: call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
66: if (size .ne. 1) then
67: SETERRA(PETSC_COMM_WORLD,1,'requires one process')
68: endif
70: big = 2.88
71: big = PETSC_INFINITY
72: zero = 0.0
73: i2 = 26
74: ! - - - - - - - - - -- - - - - - - - - - - - - - - - - - - - - - - - - -
75: ! Create nonlinear solver context
76: ! - - - - - - - - - -- - - - - - - - - - - - - - - - - - - - - - - - - -
78: call SNESCreate(PETSC_COMM_WORLD,snes,ierr)
80: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
81: ! Create matrix and vector data structures; set corresponding routines
82: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84: ! Create vectors for solution and nonlinear function
86: call VecCreateSeq(PETSC_COMM_SELF,i2,x,ierr)
87: call VecDuplicate(x,r,ierr)
90: ! Create Jacobian matrix data structure
92: call MatCreateDense(PETSC_COMM_SELF,26,26,26,26, &
93: & PETSC_NULL_SCALAR,J,ierr)
95: ! Set function evaluation routine and vector
97: call SNESSetFunction(snes,r,FormFunction,0,ierr)
99: ! Set Jacobian matrix data structure and Jacobian evaluation routine
101: call SNESSetJacobian(snes,J,J,FormJacobian,0,ierr)
103: call VecDuplicate(x,lb,ierr)
104: call VecDuplicate(x,ub,ierr)
105: call VecSet(lb,zero,ierr)
106: ! call VecGetArrayF90(lb,xx,ierr)
107: ! call ShashiLowerBound(xx)
108: ! call VecRestoreArrayF90(lb,xx,ierr)
109: call VecSet(ub,big,ierr)
110: ! call SNESVISetVariableBounds(snes,lb,ub,ierr)
112: call SNESGetLineSearch(snes,ls,ierr)
113: call SNESLineSearchSetPostCheck(ls,ShashiPostCheck, &
114: , 0,ierr)
115: call SNESSetType(snes,SNESVINEWTONRSLS,ierr)
117: call SNESSetFromOptions(snes,ierr)
119: ! set initial guess
121: call VecGetArrayF90(x,xx,ierr)
122: call ShashiInitialGuess(xx)
123: call VecRestoreArrayF90(x,xx,ierr)
125: call SNESSolve(snes,PETSC_NULL_VEC,x,ierr)
127: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128: ! Free work space. All PETSc objects should be destroyed when they
129: ! are no longer needed.
130: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131: call VecDestroy(lb,ierr)
132: call VecDestroy(ub,ierr)
133: call VecDestroy(x,ierr)
134: call VecDestroy(r,ierr)
135: call MatDestroy(J,ierr)
136: call SNESDestroy(snes,ierr)
137: call PetscFinalize(ierr)
138: end
139: !
140: ! ------------------------------------------------------------------------
141: !
142: ! FormFunction - Evaluates nonlinear function, F(x).
143: !
144: ! Input Parameters:
145: ! snes - the SNES context
146: ! x - input vector
147: ! dummy - optional user-defined context (not used here)
148: !
149: ! Output Parameter:
150: ! f - function vector
151: !
152: subroutine FormFunction(snes,x,f,dummy,ierr)
153: implicit none
155: #include <petsc/finclude/petscsys.h>
156: #include <petsc/finclude/petscvec.h>
157: #include <petsc/finclude/petscsnes.h>
159: SNES snes
160: Vec x,f
161: PetscErrorCode ierr
162: integer dummy(*)
164: ! Declarations for use with local arrays
166: PetscScalar lx_v(2),lf_v(2)
167: PetscOffset lx_i,lf_i
169: ! Get pointers to vector data.
170: ! - For default PETSc vectors, VecGetArray() returns a pointer to
171: ! the data array. Otherwise, the routine is implementation dependent.
172: ! - You MUST call VecRestoreArray() when you no longer need access to
173: ! the array.
174: ! - Note that the Fortran interface to VecGetArray() differs from the
175: ! C version. See the Fortran chapter of the users manual for details.
177: call VecGetArrayRead(x,lx_v,lx_i,ierr)
178: call VecGetArray(f,lf_v,lf_i,ierr)
179: call ShashiFormFunction(lx_a(1),lf_a(1))
180: call VecRestoreArrayRead(x,lx_v,lx_i,ierr)
181: call VecRestoreArray(f,lf_v,lf_i,ierr)
183: return
184: end
186: ! ---------------------------------------------------------------------
187: !
188: ! FormJacobian - Evaluates Jacobian matrix.
189: !
190: ! Input Parameters:
191: ! snes - the SNES context
192: ! x - input vector
193: ! dummy - optional user-defined context (not used here)
194: !
195: ! Output Parameters:
196: ! A - Jacobian matrix
197: ! B - optionally different preconditioning matrix
198: ! flag - flag indicating matrix structure
199: !
200: subroutine FormJacobian(snes,X,jac,B,dummy,ierr)
201: implicit none
203: #include <petsc/finclude/petscsys.h>
204: #include <petsc/finclude/petscvec.h>
205: #include <petsc/finclude/petscmat.h>
206: #include <petsc/finclude/petscpc.h>
207: #include <petsc/finclude/petscsnes.h>
208: #include <petsc/finclude/petscmat.h90>
209:
210: SNES snes
211: Vec X
212: Mat jac,B
213: PetscScalar A(4)
214: PetscErrorCode ierr
215: PetscInt idx(2),i2
216: integer dummy(*)
217:
218: ! Declarations for use with local arrays
220: PetscScalar lx_v(1),lf_v(1)
221: PetscOffset lx_i,lf_i
223: ! Get pointer to vector data
225: call VecGetArrayRead(x,lx_v,lx_i,ierr)
226: call MatDenseGetArray(B,lf_v,lf_i,ierr)
227: call ShashiFormJacobian(lx_a(1),lf_a(1))
228: call MatDenseRestoreArray(B,lf_v,lf_i,ierr)
229: call VecRestoreArrayRead(x,lx_v,lx_i,ierr)
230:
231: ! Assemble matrix
233: call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
234: call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
235:
236: return
237: end
239: subroutine ShashiLowerBound(an_r)
240: ! implicit PetscScalar (a-h,o-z)
241: implicit none
242: PetscScalar an_r(26)
243: PetscInt i
245: do i=2,26
246: an_r(i) = 1000.0/6.023D+23
247: enddo
248: return
249: end
251: subroutine ShashiInitialGuess(an_r)
252: ! implicit PetscScalar (a-h,o-z)
253: implicit none
254: PetscScalar an_c_additive
255: PetscScalar an_h_additive
256: PetscScalar an_o_additive
257: PetscScalar atom_c_init
258: PetscScalar atom_h_init
259: PetscScalar atom_n_init
260: PetscScalar atom_o_init
261: PetscScalar h_init
262: PetscScalar p_init
263: PetscInt nfuel
264: PetscScalar temp,pt
265: PetscScalar an_r(26),k_eq(23),f_eq(26)
266: PetscScalar d_eq(26,26),H_molar(26)
267: PetscInt an_h(1),an_c(1)
268: PetscScalar part_p(26)
269: PetscInt i_cc,i_hh,i_h2o
270: PetscInt i_pwr2,i_co2_h2o
271:
272: pt = 0.1
273: atom_c_init =6.7408177364816552D-022
274: atom_h_init = 2.0
275: atom_o_init = 1.0
276: atom_n_init = 3.76
277: nfuel = 1
278: an_c(1) = 1
279: an_h(1) = 4
280: an_c_additive = 2
281: an_h_additive = 6
282: an_o_additive = 1
283: h_init = 128799.7267952987
284: p_init = 0.1
285: temp = 1500
288: an_r( 1) = 1.66000D-24
289: an_r( 2) = 1.66030D-22
290: an_r( 3) = 5.00000D-01
291: an_r( 4) = 1.66030D-22
292: an_r( 5) = 1.66030D-22
293: an_r( 6) = 1.88000D+00
294: an_r( 7) = 1.66030D-22
295: an_r( 8) = 1.66030D-22
296: an_r( 9) = 1.66030D-22
297: an_r(10) = 1.66030D-22
298: an_r(11) = 1.66030D-22
299: an_r(12) = 1.66030D-22
300: an_r(13) = 1.66030D-22
301: an_r(14) = 1.00000D+00
302: an_r(15) = 1.66030D-22
303: an_r(16) = 1.66030D-22
304: an_r(17) = 1.66000D-24
305: an_r(18) = 1.66030D-24
306: an_r(19) = 1.66030D-24
307: an_r(20) = 1.66030D-24
308: an_r(21) = 1.66030D-24
309: an_r(22) = 1.66030D-24
310: an_r(23) = 1.66030D-24
311: an_r(24) = 1.66030D-24
312: an_r(25) = 1.66030D-24
313: an_r(26) = 1.66030D-24
315: an_r = 0
316: an_r( 3) = .5
317: an_r( 6) = 1.88000
318: an_r(14) = 1.
320:
321: #if defined(solution)
322: an_r( 1 ) = 3.802208D-33
323: an_r( 2 ) = 1.298287D-29
324: an_r( 3 ) = 2.533067D-04
325: an_r( 4 ) = 6.865078D-22
326: an_r( 5 ) = 9.993125D-01
327: an_r( 6 ) = 1.879964D+00
328: an_r( 7 ) = 4.449489D-13
329: an_r( 8 ) = 3.428687D-07
330: an_r( 9 ) = 7.105138D-05
331: an_r(10 ) = 1.094368D-04
332: an_r(11 ) = 2.362305D-06
333: an_r(12 ) = 1.107145D-09
334: an_r(13 ) = 1.276162D-24
335: an_r(14 ) = 6.315538D-04
336: an_r(15 ) = 2.356540D-09
337: an_r(16 ) = 2.048248D-09
338: an_r(17 ) = 1.966187D-22
339: an_r(18 ) = 7.856497D-29
340: an_r(19 ) = 1.987840D-36
341: an_r(20 ) = 8.182441D-22
342: an_r(21 ) = 2.684880D-16
343: an_r(22 ) = 2.680473D-16
344: an_r(23 ) = 6.594967D-18
345: an_r(24 ) = 2.509714D-21
346: an_r(25 ) = 3.096459D-21
347: an_r(26 ) = 6.149551D-18
348: #endif
349:
350: return
351: end
355: subroutine ShashiFormFunction(an_r,f_eq)
356: ! implicit PetscScalar (a-h,o-z)
357: implicit none
358: PetscScalar an_c_additive
359: PetscScalar an_h_additive
360: PetscScalar an_o_additive
361: PetscScalar atom_c_init
362: PetscScalar atom_h_init
363: PetscScalar atom_n_init
364: PetscScalar atom_o_init
365: PetscScalar h_init
366: PetscScalar p_init
367: PetscInt nfuel
368: PetscScalar temp,pt
369: PetscScalar an_r(26),k_eq(23),f_eq(26)
370: PetscScalar d_eq(26,26),H_molar(26)
371: PetscInt an_h(1),an_c(1)
372: PetscScalar part_p(26),idiff
373: PetscInt i_cc,i_hh,i_h2o
374: PetscInt i_pwr2,i_co2_h2o
375: PetscScalar an_t,sum_h,pt_cubed,pt_five
376: PetscScalar pt_four,pt_val1,pt_val2
377: PetscScalar a_io2
378: PetscInt i,ip
379: pt = 0.1
380: atom_c_init =6.7408177364816552D-022
381: atom_h_init = 2.0
382: atom_o_init = 1.0
383: atom_n_init = 3.76
384: nfuel = 1
385: an_c(1) = 1
386: an_h(1) = 4
387: an_c_additive = 2
388: an_h_additive = 6
389: an_o_additive = 1
390: h_init = 128799.7267952987
391: p_init = 0.1
392: temp = 1500
393:
394: k_eq( 1) = 1.75149D-05
395: k_eq( 2) = 4.01405D-06
396: k_eq( 3) = 6.04663D-14
397: k_eq( 4) = 2.73612D-01
398: k_eq( 5) = 3.25592D-03
399: k_eq( 6) = 5.33568D+05
400: k_eq( 7) = 2.07479D+05
401: k_eq( 8) = 1.11841D-02
402: k_eq( 9) = 1.72684D-03
403: k_eq(10) = 1.98588D-07
404: k_eq(11) = 7.23600D+27
405: k_eq(12) = 5.73926D+49
406: k_eq(13) = 1.00000D+00
407: k_eq(14) = 1.64493D+16
408: k_eq(15) = 2.73837D-29
409: k_eq(16) = 3.27419D+50
410: k_eq(17) = 1.72447D-23
411: k_eq(18) = 4.24657D-06
412: k_eq(19) = 1.16065D-14
413: k_eq(20) = 3.28020D+25
414: k_eq(21) = 1.06291D+00
415: k_eq(22) = 9.11507D+02
416: k_eq(23) = 6.02837D+03
418: H_molar( 1) = 3.26044D+03
419: H_molar( 2) = -8.00407D+04
420: H_molar( 3) = 4.05873D+04
421: H_molar( 4) = -3.31849D+05
422: H_molar( 5) = -1.93654D+05
423: H_molar( 6) = 3.84035D+04
424: H_molar( 7) = 4.97589D+05
425: H_molar( 8) = 2.74483D+05
426: H_molar( 9) = 1.30022D+05
427: H_molar(10) = 7.58429D+04
428: H_molar(11) = 2.42948D+05
429: H_molar(12) = 1.44588D+05
430: H_molar(13) = -7.16891D+04
431: H_molar(14) = 3.63075D+04
432: H_molar(15) = 9.23880D+04
433: H_molar(16) = 6.50477D+04
434: H_molar(17) = 3.04310D+05
435: H_molar(18) = 7.41707D+05
436: H_molar(19) = 6.32767D+05
437: H_molar(20) = 8.90624D+05
438: H_molar(21) = 2.49805D+04
439: H_molar(22) = 6.43473D+05
440: H_molar(23) = 1.02861D+06
441: H_molar(24) = -6.07503D+03
442: H_molar(25) = 1.27020D+05
443: H_molar(26) = -1.07011D+05
444: !=============
445: an_t = 0
446: sum_h = 0
448: do i = 1,26
449: an_t = an_t + an_r(i)
450: enddo
452: f_eq(1) = atom_h_init
453: , - (an_h(1)*an_r(1) + an_h_additive*an_r(2)
454: , + 2*an_r(5) + an_r(10) + an_r(11) + 2*an_r(14)
455: , + an_r(16)+ 2*an_r(17) + an_r(19)
456: , +an_r(20) + 3*an_r(22)+an_r(26))
459: f_eq(2) = atom_o_init
460: , - (an_o_additive*an_r(2) + 2*an_r(3)
461: , + 2*an_r(4) + an_r(5)
462: , + an_r(8) + an_r(9) + an_r(10) + an_r(12) + an_r(13)
463: , + 2*an_r(15) + 2*an_r(16)+ an_r(20) + an_r(22)
464: , + an_r(23) + 2*an_r(24) + 1*an_r(25)+an_r(26))
467: f_eq(3) = an_r(2)-1.0d-150
469: f_eq(4) = atom_c_init
470: , - (an_c(1)*an_r(1) + an_c_additive * an_r(2)
471: , + an_r(4) + an_r(13)+ 2*an_r(17) + an_r(18)
472: , + an_r(19) + an_r(20))
476:
477: do ip = 1,26
478: part_p(ip) = (an_r(ip)/an_t) * pt
479: enddo
481: i_cc = an_c(1)
482: i_hh = an_h(1)
483: a_io2 = i_cc + i_hh/4.0
484: i_h2o = i_hh/2
485: idiff = (i_cc + i_h2o) - (a_io2 + 1)
487: f_eq(5) = k_eq(11)*an_r(1)*an_r(3)**a_io2
488: , - (an_r(4)**i_cc)*(an_r(5)**i_h2o)*((pt/an_t)**idiff)
489: ! write(6,*)f_eq(5),an_r(1), an_r(3), an_r(4),an_r(5),' II'
490: ! stop
491: f_eq(6) = atom_n_init
492: , - (2*an_r(6) + an_r(7) + an_r(9) + 2*an_r(12)
493: , + an_r(15)
494: , + an_r(23))
499: f_eq( 7) = part_p(11)
500: , - (k_eq( 1) * sqrt(part_p(14)+1d-23))
501: f_eq( 8) = part_p( 8)
502: , - (k_eq( 2) * sqrt(part_p( 3)+1d-23))
504: f_eq( 9) = part_p( 7)
505: , - (k_eq( 3) * sqrt(part_p( 6)+1d-23))
507: f_eq(10) = part_p(10)
508: , - (k_eq( 4) * sqrt(part_p( 3)+1d-23))
509: , * sqrt(part_p(14))
511: f_eq(11) = part_p( 9)
512: , - (k_eq( 5) * sqrt(part_p( 3)+1d-23))
513: , * sqrt(part_p( 6)+1d-23)
514: f_eq(12) = part_p( 5)
515: , - (k_eq( 6) * sqrt(part_p( 3)+1d-23))
516: , * (part_p(14))
519: f_eq(13) = part_p( 4)
520: , - (k_eq( 7) * sqrt(part_p(3)+1.0d-23))
521: , * (part_p(13))
523: f_eq(14) = part_p(15)
524: , - (k_eq( 8) * sqrt(part_p(3)+1.0d-50))
525: , * (part_p( 9))
526: f_eq(15) = part_p(16)
527: , - (k_eq( 9) * part_p( 3))
528: , * sqrt(part_p(14)+1d-23)
530: f_eq(16) = part_p(12)
531: , - (k_eq(10) * sqrt(part_p( 3)+1d-23))
532: , * (part_p( 6))
534: f_eq(17) = part_p(14)*part_p(18)**2
535: , - (k_eq(15) * part_p(17))
537: f_eq(18) = (part_p(13)**2)
538: , - (k_eq(16) * part_p(3)*part_p(18)**2)
539: print*,f_eq(18),part_p(3),part_p(18),part_p(13),k_eq(16)
541: f_eq(19) = part_p(19)*part_p(3) - k_eq(17)*part_p(13)*part_p(10)
543: f_eq(20) = part_p(21)*part_p(20) - k_eq(18)*part_p(19)*part_p(8)
545: f_eq(21) = part_p(21)*part_p(23) - k_eq(19)*part_p(7)*part_p(8)
548: f_eq(22) = part_p(5)*part_p(11) - k_eq(20)*part_p(21)*part_p(22)
551: f_eq(23) = part_p(24) - k_eq(21)*part_p(21)*part_p(3)
553: f_eq(24) = part_p(3)*part_p(25) - k_eq(22)*part_p(24)*part_p(8)
555: f_eq(25) = part_p(26) - k_eq(23)*part_p(21)*part_p(10)
557: f_eq(26) = -(an_r(20) + an_r(22) + an_r(23))
558: , +(an_r(21) + an_r(24) + an_r(25) + an_r(26))
560: do i = 1,26
561: write(44,*)i,f_eq(i)
562: enddo
563:
564: return
565: end
567: subroutine ShashiFormJacobian(an_r,d_eq)
568: ! implicit PetscScalar (a-h,o-z)
569: implicit none
570: PetscScalar an_c_additive
571: PetscScalar an_h_additive
572: PetscScalar an_o_additive
573: PetscScalar atom_c_init
574: PetscScalar atom_h_init
575: PetscScalar atom_n_init
576: PetscScalar atom_o_init
577: PetscScalar h_init
578: PetscScalar p_init
579: PetscInt nfuel
580: PetscScalar temp,pt
581: PetscScalar an_t,ai_o2,sum_h
582: PetscScalar an_tot1_d,an_tot1
583: PetscScalar an_tot2_d,an_tot2
584: PetscScalar const5,const3,const2
585: PetscScalar const_cube
586: PetscScalar const_five
587: PetscScalar const_four
588: PetscScalar const_six
589: PetscInt jj,jb,ii3,id,ib,ip,i
590: PetscScalar pt2,pt1
591: PetscScalar an_r(26),k_eq(23),f_eq(26)
592: PetscScalar d_eq(26,26),H_molar(26)
593: PetscInt an_h(1),an_c(1)
594: PetscScalar ai_pwr1,part_p(26),idiff
595: PetscInt i_cc,i_hh
596: PetscInt i_h2o,i_pwr2,i_co2_h2o
597: PetscScalar pt_cube,pt_five
598: PetscScalar pt_four
599: PetscScalar pt_val1,pt_val2,a_io2
600: PetscInt j
602: pt = 0.1
603: atom_c_init =6.7408177364816552D-022
604: atom_h_init = 2.0
605: atom_o_init = 1.0
606: atom_n_init = 3.76
607: nfuel = 1
608: an_c(1) = 1
609: an_h(1) = 4
610: an_c_additive = 2
611: an_h_additive = 6
612: an_o_additive = 1
613: h_init = 128799.7267952987
614: p_init = 0.1
615: temp = 1500
617: k_eq( 1) = 1.75149D-05
618: k_eq( 2) = 4.01405D-06
619: k_eq( 3) = 6.04663D-14
620: k_eq( 4) = 2.73612D-01
621: k_eq( 5) = 3.25592D-03
622: k_eq( 6) = 5.33568D+05
623: k_eq( 7) = 2.07479D+05
624: k_eq( 8) = 1.11841D-02
625: k_eq( 9) = 1.72684D-03
626: k_eq(10) = 1.98588D-07
627: k_eq(11) = 7.23600D+27
628: k_eq(12) = 5.73926D+49
629: k_eq(13) = 1.00000D+00
630: k_eq(14) = 1.64493D+16
631: k_eq(15) = 2.73837D-29
632: k_eq(16) = 3.27419D+50
633: k_eq(17) = 1.72447D-23
634: k_eq(18) = 4.24657D-06
635: k_eq(19) = 1.16065D-14
636: k_eq(20) = 3.28020D+25
637: k_eq(21) = 1.06291D+00
638: k_eq(22) = 9.11507D+02
639: k_eq(23) = 6.02837D+03
641: H_molar( 1) = 3.26044D+03
642: H_molar( 2) = -8.00407D+04
643: H_molar( 3) = 4.05873D+04
644: H_molar( 4) = -3.31849D+05
645: H_molar( 5) = -1.93654D+05
646: H_molar( 6) = 3.84035D+04
647: H_molar( 7) = 4.97589D+05
648: H_molar( 8) = 2.74483D+05
649: H_molar( 9) = 1.30022D+05
650: H_molar(10) = 7.58429D+04
651: H_molar(11) = 2.42948D+05
652: H_molar(12) = 1.44588D+05
653: H_molar(13) = -7.16891D+04
654: H_molar(14) = 3.63075D+04
655: H_molar(15) = 9.23880D+04
656: H_molar(16) = 6.50477D+04
657: H_molar(17) = 3.04310D+05
658: H_molar(18) = 7.41707D+05
659: H_molar(19) = 6.32767D+05
660: H_molar(20) = 8.90624D+05
661: H_molar(21) = 2.49805D+04
662: H_molar(22) = 6.43473D+05
663: H_molar(23) = 1.02861D+06
664: H_molar(24) = -6.07503D+03
665: H_molar(25) = 1.27020D+05
666: H_molar(26) = -1.07011D+05
667:
668: !
669: !=======
670: do jb = 1,26
671: do ib = 1,26
672: d_eq(ib,jb) = 0.0d0
673: end do
674: end do
676: an_t = 0.0
677: do id = 1,26
678: an_t = an_t + an_r(id)
679: enddo
681: !====
682: !====
683: d_eq(1,1) = -an_h(1)
684: d_eq(1,2) = -an_h_additive
685: d_eq(1,5) = -2
686: d_eq(1,10) = -1
687: d_eq(1,11) = -1
688: d_eq(1,14) = -2
689: d_eq(1,16) = -1
690: d_eq(1,17) = -2
691: d_eq(1,19) = -1
692: d_eq(1,20) = -1
693: d_eq(1,22) = -3
694: d_eq(1,26) = -1
696: d_eq(2,2) = -1*an_o_additive
697: d_eq(2,3) = -2
698: d_eq(2,4) = -2
699: d_eq(2,5) = -1
700: d_eq(2,8) = -1
701: d_eq(2,9) = -1
702: d_eq(2,10) = -1
703: d_eq(2,12) = -1
704: d_eq(2,13) = -1
705: d_eq(2,15) = -2
706: d_eq(2,16) = -2
707: d_eq(2,20) = -1
708: d_eq(2,22) = -1
709: d_eq(2,23) = -1
710: d_eq(2,24) = -2
711: d_eq(2,25) = -1
712: d_eq(2,26) = -1
716: d_eq(6,6) = -2
717: d_eq(6,7) = -1
718: d_eq(6,9) = -1
719: d_eq(6,12) = -2
720: d_eq(6,15) = -1
721: d_eq(6,23) = -1
725: d_eq(4,1) = -an_c(1)
726: d_eq(4,2) = -an_c_additive
727: d_eq(4,4) = -1
728: d_eq(4,13) = -1
729: d_eq(4,17) = -2
730: d_eq(4,18) = -1
731: d_eq(4,19) = -1
732: d_eq(4,20) = -1
735: !----------
736: const2 = an_t*an_t
737: const3 = (an_t)*sqrt(an_t)
738: const5 = an_t*const3
741: const_cube = an_t*an_t*an_t
742: const_four = const2*const2
743: const_five = const2*const_cube
744: const_six = const_cube*const_cube
745: pt_cube = pt*pt*pt
746: pt_four = pt_cube*pt
747: pt_five = pt_cube*pt*pt
749: i_cc = an_c(1)
750: i_hh = an_h(1)
751: ai_o2 = i_cc + float(i_hh)/4.0
752: i_co2_h2o = i_cc + i_hh/2
753: i_h2o = i_hh/2
754: ai_pwr1 = 1 + i_cc + float(i_hh)/4.0
755: i_pwr2 = i_cc + i_hh/2
756: idiff = (i_cc + i_h2o) - (ai_o2 + 1)
758: pt1 = pt**(ai_pwr1)
759: an_tot1 = an_t**(ai_pwr1)
760: pt_val1 = (pt/an_t)**(ai_pwr1)
761: an_tot1_d = an_tot1*an_t
763: pt2 = pt**(i_pwr2)
764: an_tot2 = an_t**(i_pwr2)
765: pt_val2 = (pt/an_t)**(i_pwr2)
766: an_tot2_d = an_tot2*an_t
769: d_eq(5,1) =
770: , -(an_r(4)**i_cc)*(an_r(5)**i_h2o)
771: , *((pt/an_t)**idiff) *(-idiff/an_t)
774: do jj = 2,26
775: d_eq(5,jj) = d_eq(5,1)
776: enddo
778: d_eq(5,1) = d_eq(5,1) + k_eq(11)* (an_r(3) **ai_o2)
780: d_eq(5,3) = d_eq(5,3) + k_eq(11)* (ai_o2*an_r(3)**(ai_o2-1))
781: , * an_r(1)
783: d_eq(5,4) = d_eq(5,4) - (i_cc*an_r(4)**(i_cc-1))*
784: , (an_r(5)**i_h2o)* ((pt/an_t)**idiff)
785: d_eq(5,5) = d_eq(5,5)
786: , - (i_h2o*(an_r(5)**(i_h2o-1)))
787: , * (an_r(4)**i_cc)* ((pt/an_t)**idiff)
791: d_eq(3,1) = -(an_r(4)**2)*(an_r(5)**3)*(pt/an_t)*(-1.0/an_t)
792: do jj = 2,26
793: d_eq(3,jj) = d_eq(3,1)
794: enddo
797: d_eq(3,2) = d_eq(3,2) + k_eq(12)* (an_r(3)**3)
801: d_eq(3,3) = d_eq(3,3) + k_eq(12)* (3*an_r(3)**2)*an_r(2)
803: d_eq(3,4) = d_eq(3,4) - 2*an_r(4)*(an_r(5)**3)*(pt/an_t)
806: d_eq(3,5) = d_eq(3,5) - 3*(an_r(5)**2)*(an_r(4)**2)*(pt/an_t)
807: ! , *(pt_five/const_five)
809: do ii3 = 1,26
810: d_eq(3,ii3) = 0.0d0
811: enddo
812: d_eq(3,2) = 1.0d0
816: d_eq(7,1) = pt*an_r(11)*(-1.0)/const2
817: , -k_eq(1)*sqrt(pt)*sqrt(an_r(14)+1d-50)*(-0.5/const3)
819: do jj = 2,26
820: d_eq(7,jj) = d_eq(7,1)
821: enddo
823: d_eq(7,11) = d_eq(7,11) + pt/an_t
824: d_eq(7,14) = d_eq(7,14)
825: , - k_eq(1)*sqrt(pt)*(0.5/(sqrt((an_r(14)+1d-50)*an_t)))
828: d_eq(8,1) = pt*an_r(8)*(-1.0)/const2
829: , -k_eq(2)*sqrt(pt)*sqrt(an_r(3)+1.0d-50)*(-0.5/const3)
831: do jj = 2,26
832: d_eq(8,jj) = d_eq(8,1)
833: enddo
835: d_eq(8,3) = d_eq(8,3)
836: , -k_eq(2)*sqrt(pt)*(0.5/(sqrt((an_r(3)+1.0d-50)*an_t)))
837: d_eq(8,8) = d_eq(8,8) + pt/an_t
840: d_eq(9,1) = pt*an_r(7)*(-1.0)/const2
841: , -k_eq(3)*sqrt(pt)*sqrt(an_r(6))*(-0.5/const3)
843: do jj = 2,26
844: d_eq(9,jj) = d_eq(9,1)
845: enddo
847: d_eq(9,7) = d_eq(9,7) + pt/an_t
848: d_eq(9,6) = d_eq(9,6)
849: , -k_eq(3)*sqrt(pt)*(0.5/(sqrt(an_r(6)*an_t)))
852: d_eq(10,1) = pt*an_r(10)*(-1.0)/const2
853: , -k_eq(4)*(pt)*sqrt((an_r(3)+1.0d-50)
854: , *an_r(14))*(-1.0/const2)
855: do jj = 2,26
856: d_eq(10,jj) = d_eq(10,1)
857: enddo
859: d_eq(10,3) = d_eq(10,3)
860: , -k_eq(4)*(pt)*sqrt(an_r(14))
861: , *(0.5/(sqrt(an_r(3)+1.0d-50)*an_t))
862: d_eq(10,10) = d_eq(10,10) + pt/an_t
863: d_eq(10,14) = d_eq(10,14)
864: , -k_eq(4)*(pt)*sqrt(an_r(3)+1.0d-50)
865: , *(0.5/(sqrt(an_r(14)+1.0d-50)*an_t))
867: d_eq(11,1) = pt*an_r(9)*(-1.0)/const2
868: , -k_eq(5)*(pt)*sqrt((an_r(3)+1.0d-50)*an_r(6))
869: , *(-1.0/const2)
871: do jj = 2,26
872: d_eq(11,jj) = d_eq(11,1)
873: enddo
875: d_eq(11,3) = d_eq(11,3)
876: , -k_eq(5)*(pt)*sqrt(an_r(6))*(0.5/
877: , (sqrt(an_r(3)+1.0d-50)*an_t))
878: d_eq(11,6) = d_eq(11,6)
879: , -k_eq(5)*(pt)*sqrt(an_r(3)+1.0d-50)
880: , *(0.5/(sqrt(an_r(6))*an_t))
881: d_eq(11,9) = d_eq(11,9) + pt/an_t
885: d_eq(12,1) = pt*an_r(5)*(-1.0)/const2
886: , -k_eq(6)*(pt**1.5)*sqrt(an_r(3)+1.0d-50)
887: , *(an_r(14))*(-1.5/const5)
890: do jj = 2,26
891: d_eq(12,jj) = d_eq(12,1)
892: enddo
894: d_eq(12,3) = d_eq(12,3)
895: , -k_eq(6)*(pt**1.5)*((an_r(14)+1.0d-50)/const3)
896: , *(0.5/sqrt(an_r(3)+1.0d-50))
898: d_eq(12,5) = d_eq(12,5) + pt/an_t
899: d_eq(12,14) = d_eq(12,14)
900: , -k_eq(6)*(pt**1.5)*(sqrt(an_r(3)+1.0d-50)/const3)
903: d_eq(13,1) = pt*an_r(4)*(-1.0)/const2
904: , -k_eq(7)*(pt**1.5)*sqrt(an_r(3)+1.0d-50)
905: , *(an_r(13))*(-1.5/const5)
907: do jj = 2,26
908: d_eq(13,jj) = d_eq(13,1)
909: enddo
911: d_eq(13,3) = d_eq(13,3)
912: , -k_eq(7)*(pt**1.5)*(an_r(13)/const3)
913: , *(0.5/sqrt(an_r(3)+1.0d-50))
915: d_eq(13,4) = d_eq(13,4) + pt/an_t
916: d_eq(13,13) = d_eq(13,13)
917: , -k_eq(7)*(pt**1.5)*(sqrt(an_r(3)+1.0d-50)/const3)
919:
922: d_eq(14,1) = pt*an_r(15)*(-1.0)/const2
923: , -k_eq(8)*(pt**1.5)*sqrt(an_r(3)+1.0d-50)
924: , *(an_r(9))*(-1.5/const5)
926: do jj = 2,26
927: d_eq(14,jj) = d_eq(14,1)
928: enddo
930: d_eq(14,3) = d_eq(14,3)
931: , -k_eq(8)*(pt**1.5)*(an_r(9)/const3)
932: , *(0.5/sqrt(an_r(3)+1.0d-50))
933: d_eq(14,9) = d_eq(14,9)
934: , -k_eq(8)*(pt**1.5)*(sqrt(an_r(3)+1.0d-50)/const3)
935: d_eq(14,15) = d_eq(14,15)+ pt/an_t
939: d_eq(15,1) = pt*an_r(16)*(-1.0)/const2
940: , -k_eq(9)*(pt**1.5)*sqrt(an_r(14)+1.0d-50)
941: , *(an_r(3))*(-1.5/const5)
943: do jj = 2,26
944: d_eq(15,jj) = d_eq(15,1)
945: enddo
947: d_eq(15,3) = d_eq(15,3)
948: , -k_eq(9)*(pt**1.5)*(sqrt(an_r(14)+1.0d-50)/const3)
949: d_eq(15,14) = d_eq(15,14)
950: , -k_eq(9)*(pt**1.5)*(an_r(3)/const3)
951: , *(0.5/sqrt(an_r(14)+1.0d-50))
952: d_eq(15,16) = d_eq(15,16) + pt/an_t
955: d_eq(16,1) = pt*an_r(12)*(-1.0)/const2
956: , -k_eq(10)*(pt**1.5)*sqrt(an_r(3)+1.0d-50)
957: , *(an_r(6))*(-1.5/const5)
959: do jj = 2,26
960: d_eq(16,jj) = d_eq(16,1)
961: enddo
963: d_eq(16,3) = d_eq(16,3)
964: , -k_eq(10)*(pt**1.5)*(an_r(6)/const3)
965: , *(0.5/sqrt(an_r(3)+1.0d-50))
967: d_eq(16,6) = d_eq(16,6)
968: , -k_eq(10)*(pt**1.5)*(sqrt(an_r(3)+1.0d-50)/const3)
969: d_eq(16,12) = d_eq(16,12) + pt/an_t
975: const_cube = an_t*an_t*an_t
976: const_four = const2*const2
979: d_eq(17,1) = an_r(14)*an_r(18)*an_r(18)*(pt**3)*(-3/const_four)
980: , - k_eq(15) * an_r(17)*pt * (-1/const2)
981: do jj = 2,26
982: d_eq(17,jj) = d_eq(17,1)
983: enddo
984: d_eq(17,14) = d_eq(17,14) + an_r(18)*an_r(18)*(pt**3)/const_cube
985: d_eq(17,17) = d_eq(17,17) - k_eq(15)*pt/an_t
986: d_eq(17,18) = d_eq(17,18) + 2*an_r(18)*an_r(14)
987: , *(pt**3)/const_cube
990: d_eq(18,1) = an_r(13)*an_r(13)*(pt**2)*(-2/const_cube)
991: , - k_eq(16) * an_r(3)*an_r(18)*an_r(18)
992: , * (pt*pt*pt) * (-3/const_four)
993: do jj = 2,26
994: d_eq(18,jj) = d_eq(18,1)
995: enddo
996: d_eq(18,3) = d_eq(18,3)
997: , - k_eq(16) *an_r(18)* an_r(18)*pt*pt*pt /const_cube
998: d_eq(18,13) = d_eq(18,13)
999: , + 2* an_r(13)*pt*pt /const2
1000: d_eq(18,18) = d_eq(18,18) -k_eq(16)*an_r(3)
1001: , * 2*an_r(18)*pt*pt*pt/const_cube
1005: !====for eq 19
1007: d_eq(19,1) = an_r(3)*an_r(19)*(pt**2)*(-2/const_cube)
1008: , - k_eq(17)*an_r(13)*an_r(10)*pt*pt * (-2/const_cube)
1009: do jj = 2,26
1010: d_eq(19,jj) = d_eq(19,1)
1011: enddo
1012: d_eq(19,13) = d_eq(19,13)
1013: , - k_eq(17) *an_r(10)*pt*pt /const2
1014: d_eq(19,10) = d_eq(19,10)
1015: , - k_eq(17) *an_r(13)*pt*pt /const2
1016: d_eq(19,3) = d_eq(19,3) + an_r(19)*pt*pt/const2
1017: d_eq(19,19) = d_eq(19,19) + an_r(3)*pt*pt/const2
1018: !====for eq 20
1020: d_eq(20,1) = an_r(21)*an_r(20)*(pt**2)*(-2/const_cube)
1021: , - k_eq(18) * an_r(19)*an_r(8)*pt*pt * (-2/const_cube)
1022: do jj = 2,26
1023: d_eq(20,jj) = d_eq(20,1)
1024: enddo
1025: d_eq(20,8) = d_eq(20,8)
1026: , - k_eq(18) *an_r(19)*pt*pt /const2
1027: d_eq(20,19) = d_eq(20,19)
1028: , - k_eq(18) *an_r(8)*pt*pt /const2
1029: d_eq(20,20) = d_eq(20,20) + an_r(21)*pt*pt/const2
1030: d_eq(20,21) = d_eq(20,21) + an_r(20)*pt*pt/const2
1032: !========
1033: !====for eq 21
1035: d_eq(21,1) = an_r(21)*an_r(23)*(pt**2)*(-2/const_cube)
1036: , - k_eq(19)*an_r(7)*an_r(8)*pt*pt * (-2/const_cube)
1037: do jj = 2,26
1038: d_eq(21,jj) = d_eq(21,1)
1039: enddo
1040: d_eq(21,7) = d_eq(21,7)
1041: , - k_eq(19) *an_r(8)*pt*pt /const2
1042: d_eq(21,8) = d_eq(21,8)
1043: , - k_eq(19) *an_r(7)*pt*pt /const2
1044: d_eq(21,21) = d_eq(21,21) + an_r(23)*pt*pt/const2
1045: d_eq(21,23) = d_eq(21,23) + an_r(21)*pt*pt/const2
1047: !========
1048: ! for 22
1049: d_eq(22,1) = an_r(5)*an_r(11)*(pt**2)*(-2/const_cube)
1050: , -k_eq(20)*an_r(21)*an_r(22)*pt*pt * (-2/const_cube)
1051: do jj = 2,26
1052: d_eq(22,jj) = d_eq(22,1)
1053: enddo
1054: d_eq(22,21) = d_eq(22,21)
1055: , - k_eq(20) *an_r(22)*pt*pt /const2
1056: d_eq(22,22) = d_eq(22,22)
1057: , - k_eq(20) *an_r(21)*pt*pt /const2
1058: d_eq(22,11) = d_eq(22,11) + an_r(5)*pt*pt/(const2)
1059: d_eq(22,5) = d_eq(22,5) + an_r(11)*pt*pt/(const2)
1063: !========
1064: ! for 23
1066: d_eq(23,1) = an_r(24)*(pt)*(-1/const2)
1067: , - k_eq(21)*an_r(21)*an_r(3)*pt*pt * (-2/const_cube)
1068: do jj = 2,26
1069: d_eq(23,jj) = d_eq(23,1)
1070: enddo
1071: d_eq(23,3) = d_eq(23,3)
1072: , - k_eq(21) *an_r(21)*pt*pt /const2
1073: d_eq(23,21) = d_eq(23,21)
1074: , - k_eq(21) *an_r(3)*pt*pt /const2
1075: d_eq(23,24) = d_eq(23,24) + pt/(an_t)
1077: !========
1078: ! for 24
1079: d_eq(24,1) = an_r(3)*an_r(25)*(pt**2)*(-2/const_cube)
1080: , - k_eq(22)*an_r(24)*an_r(8)*pt*pt * (-2/const_cube)
1081: do jj = 2,26
1082: d_eq(24,jj) = d_eq(24,1)
1083: enddo
1084: d_eq(24,8) = d_eq(24,8)
1085: , - k_eq(22) *an_r(24)*pt*pt /const2
1086: d_eq(24,24) = d_eq(24,24)
1087: , - k_eq(22) *an_r(8)*pt*pt /const2
1088: d_eq(24,3) = d_eq(24,3) + an_r(25)*pt*pt/const2
1089: d_eq(24,25) = d_eq(24,25) + an_r(3)*pt*pt/const2
1091: !========
1092: !for 25
1093:
1094: d_eq(25,1) = an_r(26)*(pt)*(-1/const2)
1095: , - k_eq(23)*an_r(21)*an_r(10)*pt*pt * (-2/const_cube)
1096: do jj = 2,26
1097: d_eq(25,jj) = d_eq(25,1)
1098: enddo
1099: d_eq(25,10) = d_eq(25,10)
1100: , - k_eq(23) *an_r(21)*pt*pt /const2
1101: d_eq(25,21) = d_eq(25,21)
1102: , - k_eq(23) *an_r(10)*pt*pt /const2
1103: d_eq(25,26) = d_eq(25,26) + pt/(an_t)
1105: !============
1106: ! for 26
1107: d_eq(26,20) = -1
1108: d_eq(26,22) = -1
1109: d_eq(26,23) = -1
1110: d_eq(26,21) = 1
1111: d_eq(26,24) = 1
1112: d_eq(26,25) = 1
1113: d_eq(26,26) = 1
1116: do j = 1,26
1117: do i = 1,26
1118: write(44,*)i,j,d_eq(i,j)
1119: enddo
1120: enddo
1122:
1123: return
1124: end
1126: subroutine ShashiPostCheck(ls,X,Y,W,c_Y,c_W,dummy)
1127: implicit none
1129: #include <petsc/finclude/petscsys.h>
1130: #include <petsc/finclude/petscvec.h>
1131: #include <petsc/finclude/petscsnes.h>
1132: #include <petsc/finclude/petscvec.h90>
1133: SNESLineSearch ls
1134: PetscErrorCode ierr
1135: Vec X,Y,W
1136: PetscObject dummy
1137: PetscBool c_Y,c_W
1138: PetscScalar,pointer :: xx(:)
1139: PetscInt i
1140: call VecGetArrayF90(W,xx,ierr)
1141: do i=1,26
1142: if (xx(i) < 0.0) then
1143: xx(i) = 0.0
1144: c_W = PETSC_TRUE
1145: endif
1146: if (xx(i) > 3.0) then
1147: xx(i) = 3.0
1148: endif
1149: enddo
1150: call VecRestoreArrayF90(W,xx,ierr)
1151: return
1152: end
1153: