Actual source code: ex1f.F
1: !
2: ! "$Id: ex1f.F,v 1.38 2001/01/19 23:22:13 balay Exp $";
3: !
4: ! Formatted test for TS routines.
5: !
6: ! Solves U_t = U_xx
7: ! F(t,u) = (u_i+1 - 2u_i + u_i-1)/h^2
8: ! using several different schemes.
9: !
10: !23456789012345678901234567890123456789012345678901234567890123456789012
12: program main
13: implicit none
14: #include finclude/petsc.h
15: #include finclude/petscvec.h
16: #include finclude/petscmat.h
17: #include finclude/petscda.h
18: #include finclude/petscsys.h
19: #include finclude/petscpc.h
20: #include finclude/petscsles.h
21: #include finclude/petscsnes.h
22: #include finclude/petscts.h
23: #include finclude/petscdraw.h
24: #include finclude/petscviewer.h
26: integer linear_no_matrix,linear_no_time,linear
27: integer nonlinear_no_jacobian,nonlinear
28: parameter (linear_no_matrix = 0,linear_no_time = 1,linear = 2)
29: parameter (nonlinear_no_jacobian = 3,nonlinear = 4)
31: integer ierr,time_steps,steps,flg,size,problem
32: Vec local,global
33: double precision dt,ftime,zero,tmax
34: TS ts
35: Mat A
36: MatStructure A_structure
37: TSProblemType tsproblem
38: PetscDraw draw
39: PetscViewer viewer
40: character*(60) type,tsinfo
41: character*(5) etype
43: !
44: ! Application context data, stored in common block
45: !
46: Vec localwork,csolution
47: DA da
48: PetscViewer viewer1,viewer2
49: integer M
50: double precision h,nrm_2,nrm_max,nox
51: common /tsctx/ h,nrm_2,nrm_max,nox,localwork,csolution,da,viewer1 &
52: & ,viewer2,M
54: external Monitor,RHSFunctionHeat,RHSMatrixFree,Initial
55: external RHSMatrixHeat,RHSJacobianHeat
57: zero = 0.0
58: time_steps = 100
59: problem = linear_no_matrix
60: A = 0
61: tsproblem = TS_LINEAR
62:
63: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
64: call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
66: M = 60
67: call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-M',M,flg,ierr)
68: call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-time',time_steps, &
69: & flg,ierr)
71: call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-nox',flg,ierr)
72: if (flg .eq. 1) then
73: nox = 1
74: else
75: nox = 0
76: endif
77: nrm_2 = 0.0
78: nrm_max = 0.0
80: ! Set up the ghost point communication pattern
82: call DACreate1d(PETSC_COMM_WORLD,DA_NONPERIODIC,M,1,1, &
83: & PETSC_NULL_INTEGER,da,ierr)
84: call DACreateGlobalVector(da,global,ierr)
85: call VecGetLocalSize(global,m,ierr)
86: call DACreateLocalVector(da,local,ierr)
88: ! Set up display to show wave graph
90: call PetscViewerDrawOpen(PETSC_COMM_WORLD,PETSC_NULL_CHARACTER, &
91: & PETSC_NULL_CHARACTER,80,380,400,160,viewer1,ierr)
92: call PetscViewerDrawGetDraw(viewer1,0,draw,ierr)
93: call PetscDrawSetDoubleBuffer(draw,ierr)
94: call PetscViewerDrawOpen(PETSC_COMM_WORLD,PETSC_NULL_CHARACTER, &
95: & PETSC_NULL_CHARACTER,80,0,400,160,viewer2,ierr)
96: call PetscViewerDrawGetDraw(viewer2,0,draw,ierr)
97: call PetscDrawSetDoubleBuffer(draw,ierr)
99: ! make work array for evaluating right hand side function
101: call VecDuplicate(local,localwork,ierr)
103: ! make work array for storing exact solution
105: call VecDuplicate(global,csolution,ierr)
107: h = 1.0/(M-1.0)
109: ! set initial conditions
110:
111: call Initial(global,PETSC_NULL_OBJECT,ierr)
112:
113: !
114: ! This example is written to allow one to easily test parts
115: ! of TS, we do not expect users to generally need to use more
116: ! then a single TSProblemType
117: !
118: call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-linear_no_matrix',&
119: & flg,ierr)
120: if (flg .eq. 1) then
121: tsproblem = TS_LINEAR
122: problem = linear_no_matrix
123: endif
124: call PetscOptionsHasName(PETSC_NULL_CHARACTER, &
125: & '-linear_constant_matrix',flg,ierr)
126: if (flg .eq. 1) then
127: tsproblem = TS_LINEAR
128: problem = linear_no_time
129: endif
130: call PetscOptionsHasName(PETSC_NULL_CHARACTER, &
131: & '-linear_variable_matrix',flg,ierr)
132: if (flg .eq. 1) then
133: tsproblem = TS_LINEAR
134: problem = linear
135: endif
136: call PetscOptionsHasName(PETSC_NULL_CHARACTER, &
137: & '-nonlinear_no_jacobian',flg,ierr)
138: if (flg .eq. 1) then
139: tsproblem = TS_NONLINEAR
140: problem = nonlinear_no_jacobian
141: endif
142: call PetscOptionsHasName(PETSC_NULL_CHARACTER, &
143: & '-nonlinear_jacobian',flg,ierr)
144: if (flg .eq. 1) then
145: tsproblem = TS_NONLINEAR
146: problem = nonlinear
147: endif
148:
149: ! make timestep context
151: call TSCreate(PETSC_COMM_WORLD,tsproblem,ts,ierr)
152: call TSSetMonitor(ts,Monitor,PETSC_NULL_OBJECT, &
153: & PETSC_NULL_FUNCTION, ierr)
155: dt = h*h/2.01
157: if (problem .eq. linear_no_matrix) then
158: !
159: ! The user provides the RHS as a Shell matrix.
160: !
161: call MatCreateShell(PETSC_COMM_WORLD,m,M,M,M, &
162: & PETSC_NULL_OBJECT,A,ierr)
163: call MatShellSetOperation(A,MATOP_MULT,RHSMatrixFree,ierr)
164: call TSSetRHSMatrix(ts,A,A,PETSC_NULL_FUNCTION, &
165: & PETSC_NULL_OBJECT,ierr)
166: else if (problem .eq. linear_no_time) then
167: !
168: ! The user provides the RHS as a matrix
169: !
170: call MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,M,M, &
171: & A,ierr)
172: call MatSetFromOptions(A,ierr)
173: call RHSMatrixHeat(ts,zero,A,A,A_structure,PETSC_NULL_OBJECT, &
174: & ierr)
175: call TSSetRHSMatrix(ts,A,A,PETSC_NULL_FUNCTION, &
176: & PETSC_NULL_OBJECT,ierr)
177: else if (problem .eq. linear) then
178: !
179: ! The user provides the RHS as a time dependent matrix
180: !
181: call MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,M,M, &
182: & A,ierr)
183: call MatSetFromOptions(A,ierr)
184: call RHSMatrixHeat(ts,zero,A,A,A_structure,PETSC_NULL_OBJECT, &
185: & ierr)
186: call TSSetRHSMatrix(ts,A,A,RHSMatrixHeat,PETSC_NULL_OBJECT, &
187: & ierr)
188: else if (problem .eq. nonlinear_no_jacobian) then
189: !
190: ! The user provides the RHS and a Shell Jacobian
191: !
192: call TSSetRHSFunction(ts,RHSFunctionHeat,PETSC_NULL_OBJECT, &
193: & ierr)
194: call MatCreateShell(PETSC_COMM_WORLD,m,M,M,M, &
195: & PETSC_NULL_OBJECT,A,ierr)
196: call MatShellSetOperation(A,MATOP_MULT,RHSMatrixFree,ierr)
197: call TSSetRHSJacobian(ts,A,A,PETSC_NULL_FUNCTION, &
198: & PETSC_NULL_OBJECT,ierr)
199: else if (problem .eq. nonlinear) then
200: !
201: ! The user provides the RHS and Jacobian
202: !
203: call TSSetRHSFunction(ts,RHSFunctionHeat,PETSC_NULL_OBJECT, &
204: & ierr)
205: call MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,M,M,A &
206: & ,ierr)
207: call MatSetFromOptions(A,ierr)
208: call RHSMatrixHeat(ts,zero,A,A,A_structure,PETSC_NULL_OBJECT, &
209: & ierr)
210: call TSSetRHSJacobian(ts,A,A,RHSJacobianHeat, &
211: &PETSC_NULL_OBJECT,ierr)
212: endif
214: call TSSetFromOptions(ts,ierr)
216: call TSSetInitialTimeStep(ts,zero,dt,ierr)
217: tmax = 100.0
218: call TSSetDuration(ts,time_steps,tmax,ierr)
219: call TSSetSolution(ts,global,ierr)
221: call TSSetUp(ts,ierr)
222: call TSStep(ts,steps,ftime,ierr)
223: call PetscViewerStringOpen(PETSC_COMM_WORLD,tsinfo,60,viewer,ierr)
224: call TSView(ts,viewer,ierr)
225: call TSGetType(ts,type,ierr)
227: call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-test',flg,ierr)
228: if (flg .eq. 1) then
229: !
230: ! copy type to string of length 5 to ensure equality test
231: ! is done correctly
232: !
233: call PetscStrncpy(etype,type,5,ierr)
234: if (etype .eq. TS_EULER) then
235: if (abs(nrm_2/steps - 0.00257244) .gt. 1.e-4) then
236: print*,'Error in Euler method: 2-norm ',nrm_2/steps, &
237: & ' expecting: ',0.00257244
238: endif
239: else
240: if (abs(nrm_2/steps - 0.00506174) .gt. 1.e-4) then
241: print*,'Error in ',tsinfo,': 2-norm ',nrm_2/steps, &
242: & ' expecting: ',0.00506174
243: endif
244: endif
245: else
246: print*,size,' Procs Avg. error 2 norm ',nrm_2/steps, &
247: & nrm_max/steps,tsinfo
248: endif
250: call PetscViewerDestroy(viewer,ierr)
251: call TSDestroy(ts,ierr)
252: call PetscViewerDestroy(viewer1,ierr)
253: call PetscViewerDestroy(viewer2,ierr)
254: call VecDestroy(localwork,ierr)
255: call VecDestroy(csolution,ierr)
256: call VecDestroy(local,ierr)
257: call VecDestroy(global,ierr)
258: call DADestroy(da,ierr)
259: if (A .ne. 0) then
260: call MatDestroy(A,ierr)
261: endif
263: call PetscFinalize(ierr)
264: end
266: ! -------------------------------------------------------------------
267:
268: subroutine Initial(global,ctx,ierr)
269: implicit none
270: #include finclude/petsc.h
271: #include finclude/petscvec.h
272: #include finclude/petscmat.h
273: #include finclude/petscda.h
274: #include finclude/petscsys.h
275: #include finclude/petscpc.h
276: #include finclude/petscsles.h
277: #include finclude/petscsnes.h
278: #include finclude/petscts.h
279: #include finclude/petscviewer.h
280: Vec global
281: integer ctx
283: Scalar localptr(1)
284: integer i,mybase,myend,ierr
285: PetscOffset idx
287: !
288: ! Application context data, stored in common block
289: !
290: Vec localwork,csolution
291: DA da
292: PetscViewer viewer1,viewer2
293: integer M
294: double precision h,nrm_2,nrm_max,nox
295: common /tsctx/ h,nrm_2,nrm_max,nox,localwork,csolution,da,viewer1 &
296: & ,viewer2,M
298: ! determine starting point of each processor
300: call VecGetOwnershipRange(global,mybase,myend,ierr)
302: ! Initialize the array
304: call VecGetArray(global,localptr,idx,ierr)
305: do 10, i=mybase,myend-1
306: localptr(i-mybase+1+idx) = sin(PETSC_PI*i*6.*h) + &
307: & 3.*sin(PETSC_PI*i*2.*h)
308: 10 continue
309: call VecRestoreArray(global,localptr,idx,ierr)
310: return
311: end
313: ! ------------------------------------------------------------------------------
314: !
315: ! Exact solution
316: !
317: subroutine Solution(t,sol,ctx)
318: implicit none
319: #include finclude/petsc.h
320: #include finclude/petscvec.h
321: #include finclude/petscmat.h
322: #include finclude/petscda.h
323: #include finclude/petscsys.h
324: #include finclude/petscpc.h
325: #include finclude/petscsles.h
326: #include finclude/petscsnes.h
327: #include finclude/petscts.h
328: #include finclude/petscviewer.h
329: double precision t
330: Vec sol
331: integer ctx
333: Scalar localptr(1),ex1,ex2,sc1,sc2
334: integer i,mybase,myend,ierr
335: PetscOffset idx
337: !
338: ! Application context data, stored in common block
339: !
340: Vec localwork,csolution
341: DA da
342: PetscViewer viewer1,viewer2
343: integer M
344: double precision h,nrm_2,nrm_max,nox
345: common /tsctx/ h,nrm_2,nrm_max,nox,localwork,csolution,da,viewer1 &
346: & ,viewer2,M
348: ! determine starting point of each processor
350: call VecGetOwnershipRange(csolution,mybase,myend,ierr)
352: ex1 = exp(-36.*PETSC_PI*PETSC_PI*t)
353: ex2 = exp(-4.*PETSC_PI*PETSC_PI*t)
354: sc1 = PETSC_PI*6.*h
355: sc2 = PETSC_PI*2.*h
356: call VecGetArray(csolution,localptr,idx,ierr)
357: do 10, i=mybase,myend-1
358: localptr(i-mybase+1+idx) = sin(i*sc1)*ex1 + 3.*sin(i*sc2)*ex2
359: 10 continue
360: call VecRestoreArray(csolution,localptr,idx,ierr)
361: return
362: end
365: ! -----------------------------------------------------------------------------------
367: subroutine Monitor(ts,step,time,global,ctx,ierr)
368: implicit none
369: #include finclude/petsc.h
370: #include finclude/petscvec.h
371: #include finclude/petscmat.h
372: #include finclude/petscda.h
373: #include finclude/petscsys.h
374: #include finclude/petscpc.h
375: #include finclude/petscsles.h
376: #include finclude/petscsnes.h
377: #include finclude/petscts.h
378: #include finclude/petscdraw.h
379: #include finclude/petscviewer.h
380: TS ts
381: integer step,ctx,ierr
382: double precision time,lnrm_2,lnrm_max
383: Vec global
384: Scalar mone
386: !
387: ! Application context data, stored in common block
388: !
389: Vec localwork,csolution
390: DA da
391: PetscViewer viewer1,viewer2
392: integer M
393: double precision h,nrm_2,nrm_max,nox
394: common /tsctx/ h,nrm_2,nrm_max,nox,localwork,csolution,da,viewer1 &
395: & ,viewer2,M
397: call VecView(global,viewer1,ierr)
399: call Solution(time,csolution,ctx)
400: mone = -1.0
401: call VecAXPY(mone,global,csolution,ierr)
402: call VecNorm(csolution,NORM_2,lnrm_2,ierr)
403: lnrm_2 = sqrt(h)*lnrm_2
404: call VecNorm(csolution,NORM_MAX,lnrm_max,ierr)
406: if (nox .eq. 0) then
407: print*,'timestep ',step,' time ',time,' norm of error ', &
408: & lnrm_2,lnrm_max
409: endif
411: nrm_2 = nrm_2 + lnrm_2
412: nrm_max = nrm_max + lnrm_max
413: call VecView(csolution,viewer2,ierr)
415: return
416: end
418: ! -----------------------------------------------------------------------
420: subroutine RHSMatrixFree(mat,x,y,ierr)
421: implicit none
422: #include finclude/petsc.h
423: #include finclude/petscvec.h
424: #include finclude/petscmat.h
425: #include finclude/petscda.h
426: #include finclude/petscsys.h
427: #include finclude/petscpc.h
428: #include finclude/petscsles.h
429: #include finclude/petscsnes.h
430: #include finclude/petscts.h
431: #include finclude/petscviewer.h
432: Mat mat
433: Vec x,y
434: integer ierr
435: double precision zero
437: zero = 0.0
439: call RHSFunctionHeat(0,zero,x,y,0,ierr)
440: return
441: end
443: ! -------------------------------------------------------------------------
445: subroutine RHSFunctionHeat(ts,t,globalin,globalout,ctx,ierr)
446: implicit none
447: #include finclude/petsc.h
448: #include finclude/petscvec.h
449: #include finclude/petscmat.h
450: #include finclude/petscda.h
451: #include finclude/petscsys.h
452: #include finclude/petscpc.h
453: #include finclude/petscsles.h
454: #include finclude/petscsnes.h
455: #include finclude/petscts.h
456: #include finclude/petscviewer.h
457: TS ts
458: double precision t
459: Vec globalin,globalout
460: integer ctx
461: Vec local
462: integer ierr,i,localsize
463: PetscOffset ldx,cdx
464: Scalar copyptr(1),localptr(1),sc
465: !
466: ! Application context data, stored in common block
467: !
468: Vec localwork,csolution
469: DA da
470: PetscViewer viewer1,viewer2
471: integer M
472: double precision h,nrm_2,nrm_max,nox
473: common /tsctx/ h,nrm_2,nrm_max,nox,localwork,csolution,da,viewer1 &
474: & ,viewer2,M
476: ! Extract local array
478: call DACreateLocalVector(da,local,ierr)
479: call DAGlobalToLocalBegin(da,globalin,INSERT_VALUES,local,ierr)
480: call DAGlobalToLocalEnd(da,globalin,INSERT_VALUES,local,ierr)
481: call VecGetArray(local,localptr,ldx,ierr)
483: ! Extract work vector
485: call VecGetArray(localwork,copyptr,cdx,ierr)
487: ! Update Locally - Make array of new values
488: ! Note: For the first and last entry I copy the value
489: ! if this is an interior node it is irrelevant
491: sc = 1.0/(h*h)
492: call VecGetLocalSize(local,localsize,ierr)
493: copyptr(1+cdx) = localptr(1+ldx)
494: do 10, i=1,localsize-2
495: copyptr(i+1+cdx) = sc * (localptr(i+2+ldx) + localptr(i+ldx) - &
496: & 2.0*localptr(i+1+ldx))
497: 10 continue
498: copyptr(localsize-1+1+cdx) = localptr(localsize-1+1+ldx)
499: call VecRestoreArray(localwork,copyptr,cdx,ierr)
500: call VecRestoreArray(local,localptr,ldx,ierr)
501: call VecDestroy(local,ierr)
503: ! Local to Global
505: call DALocalToGlobal(da,localwork,INSERT_VALUES,globalout,ierr)
506: return
507: end
509: ! ---------------------------------------------------------------------
511: subroutine RHSMatrixHeat(ts,t,AA,BB,str,ctx,ierr)
512: implicit none
513: #include finclude/petsc.h
514: #include finclude/petscvec.h
515: #include finclude/petscmat.h
516: #include finclude/petscda.h
517: #include finclude/petscsys.h
518: #include finclude/petscpc.h
519: #include finclude/petscsles.h
520: #include finclude/petscsnes.h
521: #include finclude/petscts.h
522: #include finclude/petscviewer.h
523: Mat AA,BB
524: double precision t
525: TS ts
526: MatStructure str
527: integer ctx,ierr
529: Mat A
530: integer i,mstart,mend,rank,size,idx(3)
531: Scalar v(3),stwo,sone
533: !
534: ! Application context data, stored in common block
535: !
536: Vec localwork,csolution
537: DA da
538: PetscViewer viewer1,viewer2
539: integer M
540: double precision h,nrm_2,nrm_max,nox
541: common /tsctx/ h,nrm_2,nrm_max,nox,localwork,csolution,da,viewer1 &
542: & ,viewer2,M
544: A = AA
545: stwo = -2./(h*h)
546: sone = -.5*stwo
547: str = SAME_NONZERO_PATTERN
549: call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
550: call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
552: call MatGetOwnershipRange(A,mstart,mend,ierr)
553: if (mstart .eq. 0) then
554: v(1) = 1.0
555: call MatSetValues(A,1,mstart,1,mstart,v,INSERT_VALUES,ierr)
556: mstart = mstart + 1
557: endif
558: if (mend .eq. M) then
559: mend = mend - 1
560: v(1) = 1.0
561: call MatSetValues(A,1,mend,1,mend,v,INSERT_VALUES,ierr)
562: endif
564: !
565: ! Construct matrice one row at a time
566: !
567: v(1) = sone
568: v(2) = stwo
569: v(3) = sone
570: do 10, i=mstart,mend-1
571: idx(1) = i-1
572: idx(2) = i
573: idx(3) = i+1
574: call MatSetValues(A,1,i,3,idx,v,INSERT_VALUES,ierr)
575: 10 continue
577: call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
578: call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)
579: return
580: end
582: ! --------------------------------------------------------------------------------------
584: subroutine RHSJacobianHeat(ts,t,x,AA,BB,str,ctx,ierr)
585: implicit none
586: #include finclude/petsc.h
587: #include finclude/petscvec.h
588: #include finclude/petscmat.h
589: #include finclude/petscda.h
590: #include finclude/petscsys.h
591: #include finclude/petscpc.h
592: #include finclude/petscsles.h
593: #include finclude/petscsnes.h
594: #include finclude/petscts.h
595: #include finclude/petscviewer.h
596: TS ts
597: double precision t
598: Vec x
599: Mat AA,BB
600: MatStructure str
601: integer ctx,ierr
603: call RHSMatrixHeat(ts,t,AA,BB,str,ctx,ierr)
604: return
605: end