Actual source code: shashi.F
petsc-3.7.5 2017-01-01
2:
3: program main
4: implicit none
6: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
7: ! Include files
8: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: !
10: ! The following include statements are generally used in SNES Fortran
11: ! programs:
12: ! petscsys.h - base PETSc routines
13: ! petscvec.h - vectors
14: ! petscmat.h - matrices
15: ! petscksp.h - Krylov subspace methods
16: ! petscpc.h - preconditioners
17: ! petscsnes.h - SNES interface
18: ! Other include statements may be needed if using additional PETSc
19: ! routines in a Fortran program, e.g.,
20: ! petscviewer.h - viewers
21: ! petscis.h - index sets
22: !
23: #include <petsc/finclude/petscsys.h>
24: #include <petsc/finclude/petscvec.h>
25: #include <petsc/finclude/petscmat.h>
26: #include <petsc/finclude/petscksp.h>
27: #include <petsc/finclude/petscpc.h>
28: #include <petsc/finclude/petscsnes.h>
29: #include <petsc/finclude/petscvec.h90>
30: !
31: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
32: ! Variable declarations
33: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
34: !
35: ! Variables:
36: ! snes - nonlinear solver
37: ! ksp - linear solver
38: ! pc - preconditioner context
39: ! ksp - Krylov subspace method context
40: ! x, r - solution, residual vectors
41: ! J - Jacobian matrix
42: ! its - iterations for convergence
43: !
44: SNES snes
45: PC pc
46: KSP ksp
47: Vec x,r,lb,ub
48: Mat J
49: SNESLineSearch linesearch
50: PetscErrorCode ierr
51: PetscInt its,i2,i20
52: PetscMPIInt size,rank
53: PetscScalar pfive
54: PetscReal tol
55: PetscBool setls
56: PetscScalar,pointer :: xx(:)
57: PetscScalar zero,big
58: SNESLineSearch ls
59:
60: ! Note: Any user-defined Fortran routines (such as FormJacobian)
61: ! MUST be declared as external.
63: external FormFunction, FormJacobian
64: external ShashiPostCheck
66: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67: ! Macro definitions
68: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
69: !
70: ! Macros to make clearer the process of setting values in vectors and
71: ! getting values from vectors. These vectors are used in the routines
72: ! FormFunction() and FormJacobian().
73: ! - The element lx_a(ib) is element ib in the vector x
74: !
75: #define lx_a(ib) lx_v(lx_i + (ib))
76: #define lf_a(ib) lf_v(lf_i + (ib))
77: !
78: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79: ! Beginning of program
80: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
82: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
83: call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
84: call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
85: if (size .ne. 1) then
86: if (rank .eq. 0) then
87: write(6,*) 'This is a uniprocessor example only!'
88: endif
89: SETERRQ(PETSC_COMM_SELF,1,' ',ierr)
90: endif
92: big = 2.88
93: big = PETSC_INFINITY
94: zero = 0.0
95: i2 = 26
96: ! - - - - - - - - - -- - - - - - - - - - - - - - - - - - - - - - - - - -
97: ! Create nonlinear solver context
98: ! - - - - - - - - - -- - - - - - - - - - - - - - - - - - - - - - - - - -
100: call SNESCreate(PETSC_COMM_WORLD,snes,ierr)
102: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103: ! Create matrix and vector data structures; set corresponding routines
104: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106: ! Create vectors for solution and nonlinear function
108: call VecCreateSeq(PETSC_COMM_SELF,i2,x,ierr)
109: call VecDuplicate(x,r,ierr)
112: ! Create Jacobian matrix data structure
114: call MatCreateDense(PETSC_COMM_SELF,26,26,26,26, &
115: & PETSC_NULL_SCALAR,J,ierr)
117: ! Set function evaluation routine and vector
119: call SNESSetFunction(snes,r,FormFunction,PETSC_NULL_OBJECT,ierr)
121: ! Set Jacobian matrix data structure and Jacobian evaluation routine
123: call SNESSetJacobian(snes,J,J,FormJacobian,PETSC_NULL_OBJECT, &
124: & ierr)
126: call VecDuplicate(x,lb,ierr)
127: call VecDuplicate(x,ub,ierr)
128: call VecSet(lb,zero,ierr)
129: ! call VecGetArrayF90(lb,xx,ierr)
130: ! call ShashiLowerBound(xx)
131: ! call VecRestoreArrayF90(lb,xx,ierr)
132: call VecSet(ub,big,ierr)
133: ! call SNESVISetVariableBounds(snes,lb,ub,ierr)
135: call SNESGetLineSearch(snes,ls,ierr)
136: call SNESLineSearchSetPostCheck(ls,ShashiPostCheck, &
137: , PETSC_NULL_OBJECT,ierr)
138: call SNESSetType(snes,SNESVINEWTONRSLS,ierr)
140: call SNESSetFromOptions(snes,ierr)
142: ! set initial guess
144: call VecGetArrayF90(x,xx,ierr)
145: call ShashiInitialGuess(xx)
146: call VecRestoreArrayF90(x,xx,ierr)
148: call SNESSolve(snes,PETSC_NULL_OBJECT,x,ierr)
150: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151: ! Free work space. All PETSc objects should be destroyed when they
152: ! are no longer needed.
153: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
154: call VecDestroy(lb,ierr)
155: call VecDestroy(ub,ierr)
156: call VecDestroy(x,ierr)
157: call VecDestroy(r,ierr)
158: call MatDestroy(J,ierr)
159: call SNESDestroy(snes,ierr)
160: call PetscFinalize(ierr)
161: end
162: !
163: ! ------------------------------------------------------------------------
164: !
165: ! FormFunction - Evaluates nonlinear function, F(x).
166: !
167: ! Input Parameters:
168: ! snes - the SNES context
169: ! x - input vector
170: ! dummy - optional user-defined context (not used here)
171: !
172: ! Output Parameter:
173: ! f - function vector
174: !
175: subroutine FormFunction(snes,x,f,dummy,ierr)
176: implicit none
178: #include <petsc/finclude/petscsys.h>
179: #include <petsc/finclude/petscvec.h>
180: #include <petsc/finclude/petscsnes.h>
182: SNES snes
183: Vec x,f
184: PetscErrorCode ierr
185: integer dummy(*)
187: ! Declarations for use with local arrays
189: PetscScalar lx_v(2),lf_v(2)
190: PetscOffset lx_i,lf_i
192: ! Get pointers to vector data.
193: ! - For default PETSc vectors, VecGetArray() returns a pointer to
194: ! the data array. Otherwise, the routine is implementation dependent.
195: ! - You MUST call VecRestoreArray() when you no longer need access to
196: ! the array.
197: ! - Note that the Fortran interface to VecGetArray() differs from the
198: ! C version. See the Fortran chapter of the users manual for details.
200: call VecGetArrayRead(x,lx_v,lx_i,ierr)
201: call VecGetArray(f,lf_v,lf_i,ierr)
202: call ShashiFormFunction(lx_a(1),lf_a(1))
203: call VecRestoreArrayRead(x,lx_v,lx_i,ierr)
204: call VecRestoreArray(f,lf_v,lf_i,ierr)
206: return
207: end
209: ! ---------------------------------------------------------------------
210: !
211: ! FormJacobian - Evaluates Jacobian matrix.
212: !
213: ! Input Parameters:
214: ! snes - the SNES context
215: ! x - input vector
216: ! dummy - optional user-defined context (not used here)
217: !
218: ! Output Parameters:
219: ! A - Jacobian matrix
220: ! B - optionally different preconditioning matrix
221: ! flag - flag indicating matrix structure
222: !
223: subroutine FormJacobian(snes,X,jac,B,dummy,ierr)
224: implicit none
226: #include <petsc/finclude/petscsys.h>
227: #include <petsc/finclude/petscvec.h>
228: #include <petsc/finclude/petscmat.h>
229: #include <petsc/finclude/petscpc.h>
230: #include <petsc/finclude/petscsnes.h>
231: #include <petsc/finclude/petscmat.h90>
232:
233: SNES snes
234: Vec X
235: Mat jac,B
236: PetscScalar A(4)
237: PetscErrorCode ierr
238: PetscInt idx(2),i2
239: integer dummy(*)
240:
241: ! Declarations for use with local arrays
243: PetscScalar lx_v(1),lf_v(1)
244: PetscOffset lx_i,lf_i
246: ! Get pointer to vector data
248: call VecGetArrayRead(x,lx_v,lx_i,ierr)
249: call MatDenseGetArray(B,lf_v,lf_i,ierr)
250: call ShashiFormJacobian(lx_a(1),lf_a(1))
251: call MatDenseRestoreArray(B,lf_v,lf_i,ierr)
252: call VecRestoreArrayRead(x,lx_v,lx_i,ierr)
253:
254: ! Assemble matrix
256: call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
257: call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
258:
259: return
260: end
262: subroutine ShashiLowerBound(an_r)
263: ! implicit PetscScalar (a-h,o-z)
264: implicit none
265: PetscScalar an_r(26)
266: PetscInt i
268: do i=2,26
269: an_r(i) = 1000.0/6.023D+23
270: enddo
271: return
272: end
274: subroutine ShashiInitialGuess(an_r)
275: ! implicit PetscScalar (a-h,o-z)
276: implicit none
277: PetscScalar an_c_additive
278: PetscScalar an_h_additive
279: PetscScalar an_o_additive
280: PetscScalar atom_c_init
281: PetscScalar atom_h_init
282: PetscScalar atom_n_init
283: PetscScalar atom_o_init
284: PetscScalar h_init
285: PetscScalar p_init
286: PetscInt nfuel
287: PetscScalar temp,pt
288: PetscScalar an_r(26),k_eq(23),f_eq(26)
289: PetscScalar d_eq(26,26),H_molar(26)
290: PetscInt an_h(1),an_c(1)
291: PetscScalar part_p(26)
292: PetscInt i_cc,i_hh,i_h2o
293: PetscInt i_pwr2,i_co2_h2o
294:
295: pt = 0.1
296: atom_c_init =6.7408177364816552D-022
297: atom_h_init = 2.0
298: atom_o_init = 1.0
299: atom_n_init = 3.76
300: nfuel = 1
301: an_c(1) = 1
302: an_h(1) = 4
303: an_c_additive = 2
304: an_h_additive = 6
305: an_o_additive = 1
306: h_init = 128799.7267952987
307: p_init = 0.1
308: temp = 1500
311: an_r( 1) = 1.66000D-24
312: an_r( 2) = 1.66030D-22
313: an_r( 3) = 5.00000D-01
314: an_r( 4) = 1.66030D-22
315: an_r( 5) = 1.66030D-22
316: an_r( 6) = 1.88000D+00
317: an_r( 7) = 1.66030D-22
318: an_r( 8) = 1.66030D-22
319: an_r( 9) = 1.66030D-22
320: an_r(10) = 1.66030D-22
321: an_r(11) = 1.66030D-22
322: an_r(12) = 1.66030D-22
323: an_r(13) = 1.66030D-22
324: an_r(14) = 1.00000D+00
325: an_r(15) = 1.66030D-22
326: an_r(16) = 1.66030D-22
327: an_r(17) = 1.66000D-24
328: an_r(18) = 1.66030D-24
329: an_r(19) = 1.66030D-24
330: an_r(20) = 1.66030D-24
331: an_r(21) = 1.66030D-24
332: an_r(22) = 1.66030D-24
333: an_r(23) = 1.66030D-24
334: an_r(24) = 1.66030D-24
335: an_r(25) = 1.66030D-24
336: an_r(26) = 1.66030D-24
338: an_r = 0
339: an_r( 3) = .5
340: an_r( 6) = 1.88000
341: an_r(14) = 1.
343:
344: #if defined(solution)
345: an_r( 1 ) = 3.802208D-33
346: an_r( 2 ) = 1.298287D-29
347: an_r( 3 ) = 2.533067D-04
348: an_r( 4 ) = 6.865078D-22
349: an_r( 5 ) = 9.993125D-01
350: an_r( 6 ) = 1.879964D+00
351: an_r( 7 ) = 4.449489D-13
352: an_r( 8 ) = 3.428687D-07
353: an_r( 9 ) = 7.105138D-05
354: an_r(10 ) = 1.094368D-04
355: an_r(11 ) = 2.362305D-06
356: an_r(12 ) = 1.107145D-09
357: an_r(13 ) = 1.276162D-24
358: an_r(14 ) = 6.315538D-04
359: an_r(15 ) = 2.356540D-09
360: an_r(16 ) = 2.048248D-09
361: an_r(17 ) = 1.966187D-22
362: an_r(18 ) = 7.856497D-29
363: an_r(19 ) = 1.987840D-36
364: an_r(20 ) = 8.182441D-22
365: an_r(21 ) = 2.684880D-16
366: an_r(22 ) = 2.680473D-16
367: an_r(23 ) = 6.594967D-18
368: an_r(24 ) = 2.509714D-21
369: an_r(25 ) = 3.096459D-21
370: an_r(26 ) = 6.149551D-18
371: #endif
372:
373: return
374: end
378: subroutine ShashiFormFunction(an_r,f_eq)
379: ! implicit PetscScalar (a-h,o-z)
380: implicit none
381: PetscScalar an_c_additive
382: PetscScalar an_h_additive
383: PetscScalar an_o_additive
384: PetscScalar atom_c_init
385: PetscScalar atom_h_init
386: PetscScalar atom_n_init
387: PetscScalar atom_o_init
388: PetscScalar h_init
389: PetscScalar p_init
390: PetscInt nfuel
391: PetscScalar temp,pt
392: PetscScalar an_r(26),k_eq(23),f_eq(26)
393: PetscScalar d_eq(26,26),H_molar(26)
394: PetscInt an_h(1),an_c(1)
395: PetscScalar part_p(26),idiff
396: PetscInt i_cc,i_hh,i_h2o
397: PetscInt i_pwr2,i_co2_h2o
398: PetscScalar an_t,sum_h,pt_cubed,pt_five
399: PetscScalar pt_four,pt_val1,pt_val2
400: PetscScalar a_io2
401: PetscInt i,ip
402: pt = 0.1
403: atom_c_init =6.7408177364816552D-022
404: atom_h_init = 2.0
405: atom_o_init = 1.0
406: atom_n_init = 3.76
407: nfuel = 1
408: an_c(1) = 1
409: an_h(1) = 4
410: an_c_additive = 2
411: an_h_additive = 6
412: an_o_additive = 1
413: h_init = 128799.7267952987
414: p_init = 0.1
415: temp = 1500
416:
417: k_eq( 1) = 1.75149D-05
418: k_eq( 2) = 4.01405D-06
419: k_eq( 3) = 6.04663D-14
420: k_eq( 4) = 2.73612D-01
421: k_eq( 5) = 3.25592D-03
422: k_eq( 6) = 5.33568D+05
423: k_eq( 7) = 2.07479D+05
424: k_eq( 8) = 1.11841D-02
425: k_eq( 9) = 1.72684D-03
426: k_eq(10) = 1.98588D-07
427: k_eq(11) = 7.23600D+27
428: k_eq(12) = 5.73926D+49
429: k_eq(13) = 1.00000D+00
430: k_eq(14) = 1.64493D+16
431: k_eq(15) = 2.73837D-29
432: k_eq(16) = 3.27419D+50
433: k_eq(17) = 1.72447D-23
434: k_eq(18) = 4.24657D-06
435: k_eq(19) = 1.16065D-14
436: k_eq(20) = 3.28020D+25
437: k_eq(21) = 1.06291D+00
438: k_eq(22) = 9.11507D+02
439: k_eq(23) = 6.02837D+03
441: H_molar( 1) = 3.26044D+03
442: H_molar( 2) = -8.00407D+04
443: H_molar( 3) = 4.05873D+04
444: H_molar( 4) = -3.31849D+05
445: H_molar( 5) = -1.93654D+05
446: H_molar( 6) = 3.84035D+04
447: H_molar( 7) = 4.97589D+05
448: H_molar( 8) = 2.74483D+05
449: H_molar( 9) = 1.30022D+05
450: H_molar(10) = 7.58429D+04
451: H_molar(11) = 2.42948D+05
452: H_molar(12) = 1.44588D+05
453: H_molar(13) = -7.16891D+04
454: H_molar(14) = 3.63075D+04
455: H_molar(15) = 9.23880D+04
456: H_molar(16) = 6.50477D+04
457: H_molar(17) = 3.04310D+05
458: H_molar(18) = 7.41707D+05
459: H_molar(19) = 6.32767D+05
460: H_molar(20) = 8.90624D+05
461: H_molar(21) = 2.49805D+04
462: H_molar(22) = 6.43473D+05
463: H_molar(23) = 1.02861D+06
464: H_molar(24) = -6.07503D+03
465: H_molar(25) = 1.27020D+05
466: H_molar(26) = -1.07011D+05
467: !=============
468: an_t = 0
469: sum_h = 0
471: do i = 1,26
472: an_t = an_t + an_r(i)
473: enddo
475: f_eq(1) = atom_h_init
476: , - (an_h(1)*an_r(1) + an_h_additive*an_r(2)
477: , + 2*an_r(5) + an_r(10) + an_r(11) + 2*an_r(14)
478: , + an_r(16)+ 2*an_r(17) + an_r(19)
479: , +an_r(20) + 3*an_r(22)+an_r(26))
482: f_eq(2) = atom_o_init
483: , - (an_o_additive*an_r(2) + 2*an_r(3)
484: , + 2*an_r(4) + an_r(5)
485: , + an_r(8) + an_r(9) + an_r(10) + an_r(12) + an_r(13)
486: , + 2*an_r(15) + 2*an_r(16)+ an_r(20) + an_r(22)
487: , + an_r(23) + 2*an_r(24) + 1*an_r(25)+an_r(26))
490: f_eq(3) = an_r(2)-1.0d-150
492: f_eq(4) = atom_c_init
493: , - (an_c(1)*an_r(1) + an_c_additive * an_r(2)
494: , + an_r(4) + an_r(13)+ 2*an_r(17) + an_r(18)
495: , + an_r(19) + an_r(20))
499:
500: do ip = 1,26
501: part_p(ip) = (an_r(ip)/an_t) * pt
502: enddo
504: i_cc = an_c(1)
505: i_hh = an_h(1)
506: a_io2 = i_cc + i_hh/4.0
507: i_h2o = i_hh/2
508: idiff = (i_cc + i_h2o) - (a_io2 + 1)
510: f_eq(5) = k_eq(11)*an_r(1)*an_r(3)**a_io2
511: , - (an_r(4)**i_cc)*(an_r(5)**i_h2o)*((pt/an_t)**idiff)
512: ! write(6,*)f_eq(5),an_r(1), an_r(3), an_r(4),an_r(5),' II'
513: ! stop
514: f_eq(6) = atom_n_init
515: , - (2*an_r(6) + an_r(7) + an_r(9) + 2*an_r(12)
516: , + an_r(15)
517: , + an_r(23))
522: f_eq( 7) = part_p(11)
523: , - (k_eq( 1) * sqrt(part_p(14)+1d-23))
524: f_eq( 8) = part_p( 8)
525: , - (k_eq( 2) * sqrt(part_p( 3)+1d-23))
527: f_eq( 9) = part_p( 7)
528: , - (k_eq( 3) * sqrt(part_p( 6)+1d-23))
530: f_eq(10) = part_p(10)
531: , - (k_eq( 4) * sqrt(part_p( 3)+1d-23))
532: , * sqrt(part_p(14))
534: f_eq(11) = part_p( 9)
535: , - (k_eq( 5) * sqrt(part_p( 3)+1d-23))
536: , * sqrt(part_p( 6)+1d-23)
537: f_eq(12) = part_p( 5)
538: , - (k_eq( 6) * sqrt(part_p( 3)+1d-23))
539: , * (part_p(14))
542: f_eq(13) = part_p( 4)
543: , - (k_eq( 7) * sqrt(part_p(3)+1.0d-23))
544: , * (part_p(13))
546: f_eq(14) = part_p(15)
547: , - (k_eq( 8) * sqrt(part_p(3)+1.0d-50))
548: , * (part_p( 9))
549: f_eq(15) = part_p(16)
550: , - (k_eq( 9) * part_p( 3))
551: , * sqrt(part_p(14)+1d-23)
553: f_eq(16) = part_p(12)
554: , - (k_eq(10) * sqrt(part_p( 3)+1d-23))
555: , * (part_p( 6))
557: f_eq(17) = part_p(14)*part_p(18)**2
558: , - (k_eq(15) * part_p(17))
560: f_eq(18) = (part_p(13)**2)
561: , - (k_eq(16) * part_p(3)*part_p(18)**2)
562: print*,f_eq(18),part_p(3),part_p(18),part_p(13),k_eq(16)
564: f_eq(19) = part_p(19)*part_p(3) - k_eq(17)*part_p(13)*part_p(10)
566: f_eq(20) = part_p(21)*part_p(20) - k_eq(18)*part_p(19)*part_p(8)
568: f_eq(21) = part_p(21)*part_p(23) - k_eq(19)*part_p(7)*part_p(8)
571: f_eq(22) = part_p(5)*part_p(11) - k_eq(20)*part_p(21)*part_p(22)
574: f_eq(23) = part_p(24) - k_eq(21)*part_p(21)*part_p(3)
576: f_eq(24) = part_p(3)*part_p(25) - k_eq(22)*part_p(24)*part_p(8)
578: f_eq(25) = part_p(26) - k_eq(23)*part_p(21)*part_p(10)
580: f_eq(26) = -(an_r(20) + an_r(22) + an_r(23))
581: , +(an_r(21) + an_r(24) + an_r(25) + an_r(26))
583: do i = 1,26
584: write(44,*)i,f_eq(i)
585: enddo
586:
587: return
588: end
590: subroutine ShashiFormJacobian(an_r,d_eq)
591: ! implicit PetscScalar (a-h,o-z)
592: implicit none
593: PetscScalar an_c_additive
594: PetscScalar an_h_additive
595: PetscScalar an_o_additive
596: PetscScalar atom_c_init
597: PetscScalar atom_h_init
598: PetscScalar atom_n_init
599: PetscScalar atom_o_init
600: PetscScalar h_init
601: PetscScalar p_init
602: PetscInt nfuel
603: PetscScalar temp,pt
604: PetscScalar an_t,ai_o2,sum_h
605: PetscScalar an_tot1_d,an_tot1
606: PetscScalar an_tot2_d,an_tot2
607: PetscScalar const5,const3,const2
608: PetscScalar const_cube
609: PetscScalar const_five
610: PetscScalar const_four
611: PetscScalar const_six
612: PetscInt jj,jb,ii3,id,ib,ip,i
613: PetscScalar pt2,pt1
614: PetscScalar an_r(26),k_eq(23),f_eq(26)
615: PetscScalar d_eq(26,26),H_molar(26)
616: PetscInt an_h(1),an_c(1)
617: PetscScalar ai_pwr1,part_p(26),idiff
618: PetscInt i_cc,i_hh
619: PetscInt i_h2o,i_pwr2,i_co2_h2o
620: PetscScalar pt_cube,pt_five
621: PetscScalar pt_four
622: PetscScalar pt_val1,pt_val2,a_io2
623: PetscInt j
625: pt = 0.1
626: atom_c_init =6.7408177364816552D-022
627: atom_h_init = 2.0
628: atom_o_init = 1.0
629: atom_n_init = 3.76
630: nfuel = 1
631: an_c(1) = 1
632: an_h(1) = 4
633: an_c_additive = 2
634: an_h_additive = 6
635: an_o_additive = 1
636: h_init = 128799.7267952987
637: p_init = 0.1
638: temp = 1500
640: k_eq( 1) = 1.75149D-05
641: k_eq( 2) = 4.01405D-06
642: k_eq( 3) = 6.04663D-14
643: k_eq( 4) = 2.73612D-01
644: k_eq( 5) = 3.25592D-03
645: k_eq( 6) = 5.33568D+05
646: k_eq( 7) = 2.07479D+05
647: k_eq( 8) = 1.11841D-02
648: k_eq( 9) = 1.72684D-03
649: k_eq(10) = 1.98588D-07
650: k_eq(11) = 7.23600D+27
651: k_eq(12) = 5.73926D+49
652: k_eq(13) = 1.00000D+00
653: k_eq(14) = 1.64493D+16
654: k_eq(15) = 2.73837D-29
655: k_eq(16) = 3.27419D+50
656: k_eq(17) = 1.72447D-23
657: k_eq(18) = 4.24657D-06
658: k_eq(19) = 1.16065D-14
659: k_eq(20) = 3.28020D+25
660: k_eq(21) = 1.06291D+00
661: k_eq(22) = 9.11507D+02
662: k_eq(23) = 6.02837D+03
664: H_molar( 1) = 3.26044D+03
665: H_molar( 2) = -8.00407D+04
666: H_molar( 3) = 4.05873D+04
667: H_molar( 4) = -3.31849D+05
668: H_molar( 5) = -1.93654D+05
669: H_molar( 6) = 3.84035D+04
670: H_molar( 7) = 4.97589D+05
671: H_molar( 8) = 2.74483D+05
672: H_molar( 9) = 1.30022D+05
673: H_molar(10) = 7.58429D+04
674: H_molar(11) = 2.42948D+05
675: H_molar(12) = 1.44588D+05
676: H_molar(13) = -7.16891D+04
677: H_molar(14) = 3.63075D+04
678: H_molar(15) = 9.23880D+04
679: H_molar(16) = 6.50477D+04
680: H_molar(17) = 3.04310D+05
681: H_molar(18) = 7.41707D+05
682: H_molar(19) = 6.32767D+05
683: H_molar(20) = 8.90624D+05
684: H_molar(21) = 2.49805D+04
685: H_molar(22) = 6.43473D+05
686: H_molar(23) = 1.02861D+06
687: H_molar(24) = -6.07503D+03
688: H_molar(25) = 1.27020D+05
689: H_molar(26) = -1.07011D+05
690:
691: !
692: !=======
693: do jb = 1,26
694: do ib = 1,26
695: d_eq(ib,jb) = 0.0d0
696: end do
697: end do
699: an_t = 0.0
700: do id = 1,26
701: an_t = an_t + an_r(id)
702: enddo
704: !====
705: !====
706: d_eq(1,1) = -an_h(1)
707: d_eq(1,2) = -an_h_additive
708: d_eq(1,5) = -2
709: d_eq(1,10) = -1
710: d_eq(1,11) = -1
711: d_eq(1,14) = -2
712: d_eq(1,16) = -1
713: d_eq(1,17) = -2
714: d_eq(1,19) = -1
715: d_eq(1,20) = -1
716: d_eq(1,22) = -3
717: d_eq(1,26) = -1
719: d_eq(2,2) = -1*an_o_additive
720: d_eq(2,3) = -2
721: d_eq(2,4) = -2
722: d_eq(2,5) = -1
723: d_eq(2,8) = -1
724: d_eq(2,9) = -1
725: d_eq(2,10) = -1
726: d_eq(2,12) = -1
727: d_eq(2,13) = -1
728: d_eq(2,15) = -2
729: d_eq(2,16) = -2
730: d_eq(2,20) = -1
731: d_eq(2,22) = -1
732: d_eq(2,23) = -1
733: d_eq(2,24) = -2
734: d_eq(2,25) = -1
735: d_eq(2,26) = -1
739: d_eq(6,6) = -2
740: d_eq(6,7) = -1
741: d_eq(6,9) = -1
742: d_eq(6,12) = -2
743: d_eq(6,15) = -1
744: d_eq(6,23) = -1
748: d_eq(4,1) = -an_c(1)
749: d_eq(4,2) = -an_c_additive
750: d_eq(4,4) = -1
751: d_eq(4,13) = -1
752: d_eq(4,17) = -2
753: d_eq(4,18) = -1
754: d_eq(4,19) = -1
755: d_eq(4,20) = -1
758: !----------
759: const2 = an_t*an_t
760: const3 = (an_t)*sqrt(an_t)
761: const5 = an_t*const3
764: const_cube = an_t*an_t*an_t
765: const_four = const2*const2
766: const_five = const2*const_cube
767: const_six = const_cube*const_cube
768: pt_cube = pt*pt*pt
769: pt_four = pt_cube*pt
770: pt_five = pt_cube*pt*pt
772: i_cc = an_c(1)
773: i_hh = an_h(1)
774: ai_o2 = i_cc + float(i_hh)/4.0
775: i_co2_h2o = i_cc + i_hh/2
776: i_h2o = i_hh/2
777: ai_pwr1 = 1 + i_cc + float(i_hh)/4.0
778: i_pwr2 = i_cc + i_hh/2
779: idiff = (i_cc + i_h2o) - (ai_o2 + 1)
781: pt1 = pt**(ai_pwr1)
782: an_tot1 = an_t**(ai_pwr1)
783: pt_val1 = (pt/an_t)**(ai_pwr1)
784: an_tot1_d = an_tot1*an_t
786: pt2 = pt**(i_pwr2)
787: an_tot2 = an_t**(i_pwr2)
788: pt_val2 = (pt/an_t)**(i_pwr2)
789: an_tot2_d = an_tot2*an_t
792: d_eq(5,1) =
793: , -(an_r(4)**i_cc)*(an_r(5)**i_h2o)
794: , *((pt/an_t)**idiff) *(-idiff/an_t)
797: do jj = 2,26
798: d_eq(5,jj) = d_eq(5,1)
799: enddo
801: d_eq(5,1) = d_eq(5,1) + k_eq(11)* (an_r(3) **ai_o2)
803: d_eq(5,3) = d_eq(5,3) + k_eq(11)* (ai_o2*an_r(3)**(ai_o2-1))
804: , * an_r(1)
806: d_eq(5,4) = d_eq(5,4) - (i_cc*an_r(4)**(i_cc-1))*
807: , (an_r(5)**i_h2o)* ((pt/an_t)**idiff)
808: d_eq(5,5) = d_eq(5,5)
809: , - (i_h2o*(an_r(5)**(i_h2o-1)))
810: , * (an_r(4)**i_cc)* ((pt/an_t)**idiff)
814: d_eq(3,1) = -(an_r(4)**2)*(an_r(5)**3)*(pt/an_t)*(-1.0/an_t)
815: do jj = 2,26
816: d_eq(3,jj) = d_eq(3,1)
817: enddo
820: d_eq(3,2) = d_eq(3,2) + k_eq(12)* (an_r(3)**3)
824: d_eq(3,3) = d_eq(3,3) + k_eq(12)* (3*an_r(3)**2)*an_r(2)
826: d_eq(3,4) = d_eq(3,4) - 2*an_r(4)*(an_r(5)**3)*(pt/an_t)
829: d_eq(3,5) = d_eq(3,5) - 3*(an_r(5)**2)*(an_r(4)**2)*(pt/an_t)
830: ! , *(pt_five/const_five)
832: do ii3 = 1,26
833: d_eq(3,ii3) = 0.0d0
834: enddo
835: d_eq(3,2) = 1.0d0
839: d_eq(7,1) = pt*an_r(11)*(-1.0)/const2
840: , -k_eq(1)*sqrt(pt)*sqrt(an_r(14)+1d-50)*(-0.5/const3)
842: do jj = 2,26
843: d_eq(7,jj) = d_eq(7,1)
844: enddo
846: d_eq(7,11) = d_eq(7,11) + pt/an_t
847: d_eq(7,14) = d_eq(7,14)
848: , - k_eq(1)*sqrt(pt)*(0.5/(sqrt((an_r(14)+1d-50)*an_t)))
851: d_eq(8,1) = pt*an_r(8)*(-1.0)/const2
852: , -k_eq(2)*sqrt(pt)*sqrt(an_r(3)+1.0d-50)*(-0.5/const3)
854: do jj = 2,26
855: d_eq(8,jj) = d_eq(8,1)
856: enddo
858: d_eq(8,3) = d_eq(8,3)
859: , -k_eq(2)*sqrt(pt)*(0.5/(sqrt((an_r(3)+1.0d-50)*an_t)))
860: d_eq(8,8) = d_eq(8,8) + pt/an_t
863: d_eq(9,1) = pt*an_r(7)*(-1.0)/const2
864: , -k_eq(3)*sqrt(pt)*sqrt(an_r(6))*(-0.5/const3)
866: do jj = 2,26
867: d_eq(9,jj) = d_eq(9,1)
868: enddo
870: d_eq(9,7) = d_eq(9,7) + pt/an_t
871: d_eq(9,6) = d_eq(9,6)
872: , -k_eq(3)*sqrt(pt)*(0.5/(sqrt(an_r(6)*an_t)))
875: d_eq(10,1) = pt*an_r(10)*(-1.0)/const2
876: , -k_eq(4)*(pt)*sqrt((an_r(3)+1.0d-50)
877: , *an_r(14))*(-1.0/const2)
878: do jj = 2,26
879: d_eq(10,jj) = d_eq(10,1)
880: enddo
882: d_eq(10,3) = d_eq(10,3)
883: , -k_eq(4)*(pt)*sqrt(an_r(14))
884: , *(0.5/(sqrt(an_r(3)+1.0d-50)*an_t))
885: d_eq(10,10) = d_eq(10,10) + pt/an_t
886: d_eq(10,14) = d_eq(10,14)
887: , -k_eq(4)*(pt)*sqrt(an_r(3)+1.0d-50)
888: , *(0.5/(sqrt(an_r(14)+1.0d-50)*an_t))
890: d_eq(11,1) = pt*an_r(9)*(-1.0)/const2
891: , -k_eq(5)*(pt)*sqrt((an_r(3)+1.0d-50)*an_r(6))
892: , *(-1.0/const2)
894: do jj = 2,26
895: d_eq(11,jj) = d_eq(11,1)
896: enddo
898: d_eq(11,3) = d_eq(11,3)
899: , -k_eq(5)*(pt)*sqrt(an_r(6))*(0.5/
900: , (sqrt(an_r(3)+1.0d-50)*an_t))
901: d_eq(11,6) = d_eq(11,6)
902: , -k_eq(5)*(pt)*sqrt(an_r(3)+1.0d-50)
903: , *(0.5/(sqrt(an_r(6))*an_t))
904: d_eq(11,9) = d_eq(11,9) + pt/an_t
908: d_eq(12,1) = pt*an_r(5)*(-1.0)/const2
909: , -k_eq(6)*(pt**1.5)*sqrt(an_r(3)+1.0d-50)
910: , *(an_r(14))*(-1.5/const5)
913: do jj = 2,26
914: d_eq(12,jj) = d_eq(12,1)
915: enddo
917: d_eq(12,3) = d_eq(12,3)
918: , -k_eq(6)*(pt**1.5)*((an_r(14)+1.0d-50)/const3)
919: , *(0.5/sqrt(an_r(3)+1.0d-50))
921: d_eq(12,5) = d_eq(12,5) + pt/an_t
922: d_eq(12,14) = d_eq(12,14)
923: , -k_eq(6)*(pt**1.5)*(sqrt(an_r(3)+1.0d-50)/const3)
926: d_eq(13,1) = pt*an_r(4)*(-1.0)/const2
927: , -k_eq(7)*(pt**1.5)*sqrt(an_r(3)+1.0d-50)
928: , *(an_r(13))*(-1.5/const5)
930: do jj = 2,26
931: d_eq(13,jj) = d_eq(13,1)
932: enddo
934: d_eq(13,3) = d_eq(13,3)
935: , -k_eq(7)*(pt**1.5)*(an_r(13)/const3)
936: , *(0.5/sqrt(an_r(3)+1.0d-50))
938: d_eq(13,4) = d_eq(13,4) + pt/an_t
939: d_eq(13,13) = d_eq(13,13)
940: , -k_eq(7)*(pt**1.5)*(sqrt(an_r(3)+1.0d-50)/const3)
942:
945: d_eq(14,1) = pt*an_r(15)*(-1.0)/const2
946: , -k_eq(8)*(pt**1.5)*sqrt(an_r(3)+1.0d-50)
947: , *(an_r(9))*(-1.5/const5)
949: do jj = 2,26
950: d_eq(14,jj) = d_eq(14,1)
951: enddo
953: d_eq(14,3) = d_eq(14,3)
954: , -k_eq(8)*(pt**1.5)*(an_r(9)/const3)
955: , *(0.5/sqrt(an_r(3)+1.0d-50))
956: d_eq(14,9) = d_eq(14,9)
957: , -k_eq(8)*(pt**1.5)*(sqrt(an_r(3)+1.0d-50)/const3)
958: d_eq(14,15) = d_eq(14,15)+ pt/an_t
962: d_eq(15,1) = pt*an_r(16)*(-1.0)/const2
963: , -k_eq(9)*(pt**1.5)*sqrt(an_r(14)+1.0d-50)
964: , *(an_r(3))*(-1.5/const5)
966: do jj = 2,26
967: d_eq(15,jj) = d_eq(15,1)
968: enddo
970: d_eq(15,3) = d_eq(15,3)
971: , -k_eq(9)*(pt**1.5)*(sqrt(an_r(14)+1.0d-50)/const3)
972: d_eq(15,14) = d_eq(15,14)
973: , -k_eq(9)*(pt**1.5)*(an_r(3)/const3)
974: , *(0.5/sqrt(an_r(14)+1.0d-50))
975: d_eq(15,16) = d_eq(15,16) + pt/an_t
978: d_eq(16,1) = pt*an_r(12)*(-1.0)/const2
979: , -k_eq(10)*(pt**1.5)*sqrt(an_r(3)+1.0d-50)
980: , *(an_r(6))*(-1.5/const5)
982: do jj = 2,26
983: d_eq(16,jj) = d_eq(16,1)
984: enddo
986: d_eq(16,3) = d_eq(16,3)
987: , -k_eq(10)*(pt**1.5)*(an_r(6)/const3)
988: , *(0.5/sqrt(an_r(3)+1.0d-50))
990: d_eq(16,6) = d_eq(16,6)
991: , -k_eq(10)*(pt**1.5)*(sqrt(an_r(3)+1.0d-50)/const3)
992: d_eq(16,12) = d_eq(16,12) + pt/an_t
998: const_cube = an_t*an_t*an_t
999: const_four = const2*const2
1002: d_eq(17,1) = an_r(14)*an_r(18)*an_r(18)*(pt**3)*(-3/const_four)
1003: , - k_eq(15) * an_r(17)*pt * (-1/const2)
1004: do jj = 2,26
1005: d_eq(17,jj) = d_eq(17,1)
1006: enddo
1007: d_eq(17,14) = d_eq(17,14) + an_r(18)*an_r(18)*(pt**3)/const_cube
1008: d_eq(17,17) = d_eq(17,17) - k_eq(15)*pt/an_t
1009: d_eq(17,18) = d_eq(17,18) + 2*an_r(18)*an_r(14)
1010: , *(pt**3)/const_cube
1013: d_eq(18,1) = an_r(13)*an_r(13)*(pt**2)*(-2/const_cube)
1014: , - k_eq(16) * an_r(3)*an_r(18)*an_r(18)
1015: , * (pt*pt*pt) * (-3/const_four)
1016: do jj = 2,26
1017: d_eq(18,jj) = d_eq(18,1)
1018: enddo
1019: d_eq(18,3) = d_eq(18,3)
1020: , - k_eq(16) *an_r(18)* an_r(18)*pt*pt*pt /const_cube
1021: d_eq(18,13) = d_eq(18,13)
1022: , + 2* an_r(13)*pt*pt /const2
1023: d_eq(18,18) = d_eq(18,18) -k_eq(16)*an_r(3)
1024: , * 2*an_r(18)*pt*pt*pt/const_cube
1028: !====for eq 19
1030: d_eq(19,1) = an_r(3)*an_r(19)*(pt**2)*(-2/const_cube)
1031: , - k_eq(17)*an_r(13)*an_r(10)*pt*pt * (-2/const_cube)
1032: do jj = 2,26
1033: d_eq(19,jj) = d_eq(19,1)
1034: enddo
1035: d_eq(19,13) = d_eq(19,13)
1036: , - k_eq(17) *an_r(10)*pt*pt /const2
1037: d_eq(19,10) = d_eq(19,10)
1038: , - k_eq(17) *an_r(13)*pt*pt /const2
1039: d_eq(19,3) = d_eq(19,3) + an_r(19)*pt*pt/const2
1040: d_eq(19,19) = d_eq(19,19) + an_r(3)*pt*pt/const2
1041: !====for eq 20
1043: d_eq(20,1) = an_r(21)*an_r(20)*(pt**2)*(-2/const_cube)
1044: , - k_eq(18) * an_r(19)*an_r(8)*pt*pt * (-2/const_cube)
1045: do jj = 2,26
1046: d_eq(20,jj) = d_eq(20,1)
1047: enddo
1048: d_eq(20,8) = d_eq(20,8)
1049: , - k_eq(18) *an_r(19)*pt*pt /const2
1050: d_eq(20,19) = d_eq(20,19)
1051: , - k_eq(18) *an_r(8)*pt*pt /const2
1052: d_eq(20,20) = d_eq(20,20) + an_r(21)*pt*pt/const2
1053: d_eq(20,21) = d_eq(20,21) + an_r(20)*pt*pt/const2
1055: !========
1056: !====for eq 21
1058: d_eq(21,1) = an_r(21)*an_r(23)*(pt**2)*(-2/const_cube)
1059: , - k_eq(19)*an_r(7)*an_r(8)*pt*pt * (-2/const_cube)
1060: do jj = 2,26
1061: d_eq(21,jj) = d_eq(21,1)
1062: enddo
1063: d_eq(21,7) = d_eq(21,7)
1064: , - k_eq(19) *an_r(8)*pt*pt /const2
1065: d_eq(21,8) = d_eq(21,8)
1066: , - k_eq(19) *an_r(7)*pt*pt /const2
1067: d_eq(21,21) = d_eq(21,21) + an_r(23)*pt*pt/const2
1068: d_eq(21,23) = d_eq(21,23) + an_r(21)*pt*pt/const2
1070: !========
1071: ! for 22
1072: d_eq(22,1) = an_r(5)*an_r(11)*(pt**2)*(-2/const_cube)
1073: , -k_eq(20)*an_r(21)*an_r(22)*pt*pt * (-2/const_cube)
1074: do jj = 2,26
1075: d_eq(22,jj) = d_eq(22,1)
1076: enddo
1077: d_eq(22,21) = d_eq(22,21)
1078: , - k_eq(20) *an_r(22)*pt*pt /const2
1079: d_eq(22,22) = d_eq(22,22)
1080: , - k_eq(20) *an_r(21)*pt*pt /const2
1081: d_eq(22,11) = d_eq(22,11) + an_r(5)*pt*pt/(const2)
1082: d_eq(22,5) = d_eq(22,5) + an_r(11)*pt*pt/(const2)
1086: !========
1087: ! for 23
1089: d_eq(23,1) = an_r(24)*(pt)*(-1/const2)
1090: , - k_eq(21)*an_r(21)*an_r(3)*pt*pt * (-2/const_cube)
1091: do jj = 2,26
1092: d_eq(23,jj) = d_eq(23,1)
1093: enddo
1094: d_eq(23,3) = d_eq(23,3)
1095: , - k_eq(21) *an_r(21)*pt*pt /const2
1096: d_eq(23,21) = d_eq(23,21)
1097: , - k_eq(21) *an_r(3)*pt*pt /const2
1098: d_eq(23,24) = d_eq(23,24) + pt/(an_t)
1100: !========
1101: ! for 24
1102: d_eq(24,1) = an_r(3)*an_r(25)*(pt**2)*(-2/const_cube)
1103: , - k_eq(22)*an_r(24)*an_r(8)*pt*pt * (-2/const_cube)
1104: do jj = 2,26
1105: d_eq(24,jj) = d_eq(24,1)
1106: enddo
1107: d_eq(24,8) = d_eq(24,8)
1108: , - k_eq(22) *an_r(24)*pt*pt /const2
1109: d_eq(24,24) = d_eq(24,24)
1110: , - k_eq(22) *an_r(8)*pt*pt /const2
1111: d_eq(24,3) = d_eq(24,3) + an_r(25)*pt*pt/const2
1112: d_eq(24,25) = d_eq(24,25) + an_r(3)*pt*pt/const2
1114: !========
1115: !for 25
1116:
1117: d_eq(25,1) = an_r(26)*(pt)*(-1/const2)
1118: , - k_eq(23)*an_r(21)*an_r(10)*pt*pt * (-2/const_cube)
1119: do jj = 2,26
1120: d_eq(25,jj) = d_eq(25,1)
1121: enddo
1122: d_eq(25,10) = d_eq(25,10)
1123: , - k_eq(23) *an_r(21)*pt*pt /const2
1124: d_eq(25,21) = d_eq(25,21)
1125: , - k_eq(23) *an_r(10)*pt*pt /const2
1126: d_eq(25,26) = d_eq(25,26) + pt/(an_t)
1128: !============
1129: ! for 26
1130: d_eq(26,20) = -1
1131: d_eq(26,22) = -1
1132: d_eq(26,23) = -1
1133: d_eq(26,21) = 1
1134: d_eq(26,24) = 1
1135: d_eq(26,25) = 1
1136: d_eq(26,26) = 1
1139: do j = 1,26
1140: do i = 1,26
1141: write(44,*)i,j,d_eq(i,j)
1142: enddo
1143: enddo
1145:
1146: return
1147: end
1149: subroutine ShashiPostCheck(ls,X,Y,W,c_Y,c_W,dummy)
1150: implicit none
1152: #include <petsc/finclude/petscsys.h>
1153: #include <petsc/finclude/petscvec.h>
1154: #include <petsc/finclude/petscsnes.h>
1155: #include <petsc/finclude/petscvec.h90>
1156: SNESLineSearch ls
1157: PetscErrorCode ierr
1158: Vec X,Y,W
1159: PetscObject dummy
1160: PetscBool c_Y,c_W
1161: PetscScalar,pointer :: xx(:)
1162: PetscInt i
1163: call VecGetArrayF90(W,xx,ierr)
1164: do i=1,26
1165: if (xx(i) < 0.0) then
1166: xx(i) = 0.0
1167: c_W = PETSC_TRUE
1168: endif
1169: if (xx(i) > 3.0) then
1170: xx(i) = 3.0
1171: endif
1172: enddo
1173: call VecRestoreArrayF90(W,xx,ierr)
1174: return
1175: end
1176: