Actual source code: ex1f.F

  1: !
  2: !       Formatted test for TS routines.
  3: !
  4: !          Solves U_t = U_xx
  5: !     F(t,u) = (u_i+1 - 2u_i + u_i-1)/h^2
  6: !       using several different schemes.
  7: !
  8: !23456789012345678901234567890123456789012345678901234567890123456789012

 10: ! This example need rewritten using newly developed TS interface!

 12:       program main
 13:       implicit none
 14: #include <finclude/petscsys.h>
 15: #include <finclude/petscvec.h>
 16: #include <finclude/petscmat.h>
 17: #include <finclude/petscdmda.h>
 18: #include <finclude/petscpc.h>
 19: #include <finclude/petscksp.h>
 20: #include <finclude/petscsnes.h>
 21: #include <finclude/petscts.h>
 22: #include <finclude/petscdraw.h>
 23: #include <finclude/petscviewer.h>

 25: #include "ex1f.h"

 27:       integer   linear_no_matrix,linear_no_time,linear
 28:       integer   nonlinear_no_jacobian,nonlinear
 29:       parameter (linear_no_matrix = 0,linear_no_time = 1,linear = 2)
 30:       parameter (nonlinear_no_jacobian = 3,nonlinear = 4)

 32:       PetscErrorCode  ierr
 33:       PetscInt time_steps,steps
 34:       PetscMPIInt size
 35:       integer problem
 36:       Vec              local,global
 37:       double precision dt,ftime,zero,tmax
 38:       TS               ts
 39:       Mat              A
 40:       MatStructure     A_structure
 41:       TSProblemType    tsproblem
 42:       PetscDraw        draw
 43:       PetscViewer      viewer
 44:       character*(60)   type,tsinfo
 45:       character*(5)    etype
 46:       PetscInt         i1,i60
 47:       PetscBool  flg
 48:       PetscSizeT five

 50:       external Monitor,RHSFunctionHeat,RHSMatrixFree,Initial
 51:       external RHSMatrixHeat,RHSJacobianHeat
 52:       external TSComputeRHSFunctionLinear,TSComputeRHSJacobianConstant

 54:       i1         = 1
 55:       i60        = 60
 56:       zero       = 0.0
 57:       time_steps = 100
 58:       problem    = linear_no_matrix
 59:       A          = 0
 60:       tsproblem  = TS_LINEAR
 61: 
 62:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 63:       call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)

 65:       M = 60
 66:       call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-M',M,flg,ierr)
 67:       call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-time',time_steps,   &
 68:      &     flg,ierr)

 70:       call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-nox',flg,ierr)
 71:       if (flg) then
 72:         nox = 1
 73:       else
 74:         nox = 0
 75:       endif
 76:       nrm_2   = 0.0
 77:       nrm_max = 0.0

 79: !   Set up the ghost point communication pattern

 81:       call DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,M,i1,i1,            &
 82:      &     PETSC_NULL_INTEGER,da,ierr)
 83:       call DMCreateGlobalVector(da,global,ierr)
 84:       call VecGetLocalSize(global,m,ierr)
 85:       call DMCreateLocalVector(da,local,ierr)

 87: !   Set up display to show wave graph

 89:       call PetscViewerDrawOpen(PETSC_COMM_WORLD,PETSC_NULL_CHARACTER,        &
 90:      &     PETSC_NULL_CHARACTER,80,380,400,160,viewer1,ierr)
 91:       call PetscViewerDrawGetDraw(viewer1,0,draw,ierr)
 92:       call PetscDrawSetDoubleBuffer(draw,ierr)
 93:       call PetscViewerDrawOpen(PETSC_COMM_WORLD,PETSC_NULL_CHARACTER,        &
 94:      &     PETSC_NULL_CHARACTER,80,0,400,160,viewer2,ierr)
 95:       call PetscViewerDrawGetDraw(viewer2,0,draw,ierr)
 96:       call PetscDrawSetDoubleBuffer(draw,ierr)

 98: !   make work array for evaluating right hand side function

100:       call VecDuplicate(local,localwork,ierr)

102: !   make work array for storing exact solution

104:       call VecDuplicate(global,csolution,ierr)

106:       h = 1.0/(M-1.0)

108: !   set initial conditions
109: 
110:       call Initial(global,PETSC_NULL_OBJECT,ierr)
111: 
112: !
113: !     This example is written to allow one to easily test parts
114: !    of TS, we do not expect users to generally need to use more
115: !    then a single TSProblemType
116: !
117:       call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-linear_no_matrix',&
118:      &                    flg,ierr)
119:       if (flg) then
120:         tsproblem = TS_LINEAR
121:         problem   = linear_no_matrix
122:       endif
123:       call PetscOptionsHasName(PETSC_NULL_CHARACTER,                    &
124:      &     '-linear_constant_matrix',flg,ierr)
125:       if (flg) then
126:         tsproblem = TS_LINEAR
127:         problem   = linear_no_time
128:       endif
129:       call PetscOptionsHasName(PETSC_NULL_CHARACTER,                    &
130:      &     '-nonlinear_no_jacobian',flg,ierr)
131:       if (flg) then
132:         tsproblem = TS_NONLINEAR
133:         problem   = nonlinear_no_jacobian
134:       endif
135:       call PetscOptionsHasName(PETSC_NULL_CHARACTER,                    &
136:      &     '-nonlinear_jacobian',flg,ierr)
137:       if (flg) then
138:         tsproblem = TS_NONLINEAR
139:         problem   = nonlinear
140:       endif
141: 
142: !   make timestep context

144:       call TSCreate(PETSC_COMM_WORLD,ts,ierr)
145:       call TSSetProblemType(ts,tsproblem,ierr)
146:       call TSMonitorSet(ts,Monitor,PETSC_NULL_OBJECT,                   &
147:      &                  PETSC_NULL_FUNCTION, ierr)

149:       dt = h*h/2.01

151:       if (problem .eq. linear_no_matrix) then
152: !
153: !         The user provides the RHS as a Shell matrix.
154: !
155:          call MatCreateShell(PETSC_COMM_WORLD,m,M,M,M,                  &
156:      &        PETSC_NULL_OBJECT,A,ierr)
157:          call MatShellSetOperation(A,MATOP_MULT,RHSMatrixFree,ierr)
158:          call TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,     &
159:      &        PETSC_NULL_OBJECT,ierr)
160:          call TSSetRHSFunction(ts,PETSC_NULL_OBJECT,                    &
161:      &        TSComputeRHSFunctionLinear,PETSC_NULL_OBJECT,ierr)
162:       else if (problem .eq. linear_no_time) then
163: !
164: !         The user provides the RHS as a matrix
165: !
166:          call MatCreate(PETSC_COMM_WORLD,A,ierr)
167:          call MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,M,ierr)
168:          call MatSetFromOptions(A,ierr)
169:          call RHSMatrixHeat(ts,zero,A,A,A_structure,PETSC_NULL_OBJECT,  &
170:      &        ierr)
171:          call TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,     &
172:      &        PETSC_NULL_OBJECT,ierr)
173:          call TSSetRHSFunction(ts,PETSC_NULL_OBJECT,                    &
174:      &        TSComputeRHSFunctionLinear,PETSC_NULL_OBJECT,ierr)
175:       else if (problem .eq. nonlinear_no_jacobian) then
176: !
177: !         The user provides the RHS and a Shell Jacobian
178: !
179:          call TSSetRHSFunction(ts,PETSC_NULL_OBJECT,RHSFunctionHeat,    &
180:      &        PETSC_NULL_OBJECT,ierr)
181:          call MatCreateShell(PETSC_COMM_WORLD,m,M,M,M,                  &
182:      &        PETSC_NULL_OBJECT,A,ierr)
183:          call MatShellSetOperation(A,MATOP_MULT,RHSMatrixFree,ierr)
184:          call TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,     &
185:      &        PETSC_NULL_OBJECT,ierr)
186:       else if (problem .eq. nonlinear) then
187: !
188: !         The user provides the RHS and Jacobian
189: !
190:          call TSSetRHSFunction(ts,PETSC_NULL_OBJECT,RHSFunctionHeat,    &
191:      &        PETSC_NULL_OBJECT,ierr)
192:          call MatCreate(PETSC_COMM_WORLD,A,ierr)
193:          call MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,M,ierr)
194:          call MatSetFromOptions(A,ierr)
195:          call RHSMatrixHeat(ts,zero,A,A,A_structure,PETSC_NULL_OBJECT,   &
196:      &        ierr)
197:          call TSSetRHSJacobian(ts,A,A,RHSJacobianHeat,                   &
198:      &        PETSC_NULL_OBJECT,ierr)
199:       endif

201:       call TSSetFromOptions(ts,ierr)

203:       call TSSetInitialTimeStep(ts,zero,dt,ierr)
204:       tmax = 100.0
205:       call TSSetDuration(ts,time_steps,tmax,ierr)
206:       call TSSetSolution(ts,global,ierr)

208:       call TSSetUp(ts,ierr)
209:       call TSStep(ts,steps,ftime,ierr)
210:       call PetscViewerStringOpen(PETSC_COMM_WORLD,tsinfo,i60,viewer,         &
211:      &                           ierr)
212:       call TSView(ts,viewer,ierr)
213:       call TSGetType(ts,type,ierr)

215:       call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-test',flg,ierr)
216:       if (flg) then
217: !
218: !         copy type to string of length 5 to ensure equality test
219: !         is done correctly
220: !
221:         five = 5
222:         call PetscStrncpy(etype,type,five,ierr)
223:         if (etype .eq. TSEULER) then
224:           if (abs(nrm_2/steps - 0.00257244) .gt. 1.e-4) then
225:             print*,'Error in Euler method: 2-norm ',nrm_2/steps,        &
226:      &            ' expecting: ',0.00257244
227:           endif
228:         else
229:           if (abs(nrm_2/steps - 0.00506174) .gt. 1.e-4) then
230:             print*,'Error in ',tsinfo,': 2-norm ',nrm_2/steps,          &
231:      &             ' expecting: ',0.00506174
232:           endif
233:         endif
234:       else
235:         print*,size,' Procs Avg. error 2 norm ',nrm_2/steps,            &
236:      &              nrm_max/steps,tsinfo
237:       endif

239:       call PetscViewerDestroy(viewer,ierr)
240:       call TSDestroy(ts,ierr)
241:       call PetscViewerDestroy(viewer1,ierr)
242:       call PetscViewerDestroy(viewer2,ierr)
243:       call VecDestroy(localwork,ierr)
244:       call VecDestroy(csolution,ierr)
245:       call VecDestroy(local,ierr)
246:       call VecDestroy(global,ierr)
247:       call DMDestroy(da,ierr)
248:       if (A .ne. 0) then
249:         call MatDestroy(A,ierr)
250:       endif

252:       call PetscFinalize(ierr)
253:       end

255: !  -------------------------------------------------------------------
256: 
257:       subroutine Initial(global,ctx,ierr)
258:       implicit none
259: #include <finclude/petscsys.h>
260: #include <finclude/petscvec.h>
261: #include <finclude/petscmat.h>
262: #include <finclude/petscdmda.h>
263: #include <finclude/petscpc.h>
264: #include <finclude/petscksp.h>
265: #include <finclude/petscsnes.h>
266: #include <finclude/petscts.h>
267: #include <finclude/petscviewer.h>

269: #include "ex1f.h"

271:       Vec         global
272:       PetscObject    ctx

274:       PetscScalar localptr(1)
275:       PetscInt     i,mybase,myend
276:       PetscErrorCode ierr
277:       PetscOffset idx

279: !   determine starting point of each processor

281:       call VecGetOwnershipRange(global,mybase,myend,ierr)

283: !   Initialize the array

285:       call VecGetArray(global,localptr,idx,ierr)
286:       do 10, i=mybase,myend-1
287:         localptr(i-mybase+1+idx) = sin(PETSC_PI*i*6.*h) +               &
288:      &                             3.*sin(PETSC_PI*i*2.*h)
289:  10   continue
290:       call VecRestoreArray(global,localptr,idx,ierr)
291:       return
292:       end

294: ! ------------------------------------------------------------------------------
295: !
296: !       Exact solution
297: !
298:       subroutine Solution(t,sol,ctx)
299:       implicit none
300: #include <finclude/petscsys.h>
301: #include <finclude/petscvec.h>
302: #include <finclude/petscmat.h>
303: #include <finclude/petscdmda.h>
304: #include <finclude/petscpc.h>
305: #include <finclude/petscksp.h>
306: #include <finclude/petscsnes.h>
307: #include <finclude/petscts.h>
308: #include <finclude/petscviewer.h>

310: #include "ex1f.h"

312:       double precision  t
313:       Vec               sol
314:       PetscObject       ctx

316:       PetscScalar localptr(1),ex1
317:       PetscScalar ex2,sc1,sc2
318:       PetscInt    i,mybase,myend
319:       PetscErrorCode ierr
320:       PetscOffset       idx

322: !   determine starting point of each processor

324:       call VecGetOwnershipRange(csolution,mybase,myend,ierr)

326:       ex1 = exp(-36.*PETSC_PI*PETSC_PI*t)
327:       ex2 = exp(-4.*PETSC_PI*PETSC_PI*t)
328:       sc1 = PETSC_PI*6.*h
329:       sc2 = PETSC_PI*2.*h
330:       call VecGetArray(csolution,localptr,idx,ierr)
331:       do 10, i=mybase,myend-1
332:         localptr(i-mybase+1+idx) = sin(i*sc1)*ex1 + 3.*sin(i*sc2)*ex2
333:  10   continue
334:       call VecRestoreArray(csolution,localptr,idx,ierr)
335:       return
336:       end


339: ! -----------------------------------------------------------------------------------

341:       subroutine Monitor(ts,step,time,global,ctx,ierr)
342:       implicit none
343: #include <finclude/petscsys.h>
344: #include <finclude/petscvec.h>
345: #include <finclude/petscmat.h>
346: #include <finclude/petscdmda.h>
347: #include <finclude/petscpc.h>
348: #include <finclude/petscksp.h>
349: #include <finclude/petscsnes.h>
350: #include <finclude/petscts.h>
351: #include <finclude/petscdraw.h>
352: #include <finclude/petscviewer.h>

354: #include "ex1f.h"
355:       TS                ts
356:       PetscInt           step
357:       PetscObject     ctx
358:       PetscErrorCode ierr
359:       double precision  time,lnrm_2,lnrm_max
360:       Vec               global
361:       PetscScalar       mone

363:       call VecView(global,viewer1,ierr)

365:       call Solution(time,csolution,ctx)
366:       mone = -1.0
367:       call VecAXPY(csolution,mone,global,ierr)
368:       call VecNorm(csolution,NORM_2,lnrm_2,ierr)
369:       lnrm_2 = sqrt(h)*lnrm_2
370:       call VecNorm(csolution,NORM_MAX,lnrm_max,ierr)

372:       if (nox .eq. 0) then
373:         print*,'timestep ',step,' time ',time,' norm of error ',        &
374:      &         lnrm_2,lnrm_max
375:       endif

377:       nrm_2   = nrm_2 + lnrm_2
378:       nrm_max = nrm_max +  lnrm_max
379:       call VecView(csolution,viewer2,ierr)

381:       return
382:       end

384: !  -----------------------------------------------------------------------

386:       subroutine RHSMatrixFree(mat,x,y,ierr)
387:       implicit none
388: #include <finclude/petscsys.h>
389: #include <finclude/petscvec.h>
390: #include <finclude/petscmat.h>
391: #include <finclude/petscdmda.h>
392: #include <finclude/petscpc.h>
393: #include <finclude/petscksp.h>
394: #include <finclude/petscsnes.h>
395: #include <finclude/petscts.h>
396: #include <finclude/petscviewer.h>
397:       Mat              mat
398:       Vec              x,y
399:       PetscErrorCode          ierr
400:       double precision zero
401:       TS               ts0

403:       zero = 0.0

405:       ts0 = PETSC_NULL_OBJECT

407:       call RHSFunctionHeat(ts0,zero,x,y,PETSC_NULL_OBJECT,ierr)
408:       return
409:       end

411: ! -------------------------------------------------------------------------

413:       subroutine RHSFunctionHeat(ts,t,globalin,globalout,ctx,ierr)
414:       implicit none
415: #include <finclude/petscsys.h>
416: #include <finclude/petscvec.h>
417: #include <finclude/petscmat.h>
418: #include <finclude/petscdmda.h>
419: #include <finclude/petscpc.h>
420: #include <finclude/petscksp.h>
421: #include <finclude/petscsnes.h>
422: #include <finclude/petscts.h>
423: #include <finclude/petscviewer.h>

425: #include "ex1f.h"
426:       TS               ts
427:       double precision t
428:       Vec globalin,globalout
429:       PetscObject ctx
430:       Vec local
431:       PetscErrorCode  ierr
432:       PetscInt i,localsize
433:       PetscOffset ldx,cdx
434:       PetscScalar copyptr(1),localptr(1)
435:       PetscScalar sc

437: !  Extract local array

439:       call DMCreateLocalVector(da,local,ierr)
440:       call DMGlobalToLocalBegin(da,globalin,INSERT_VALUES,local,ierr)
441:       call DMGlobalToLocalEnd(da,globalin,INSERT_VALUES,local,ierr)
442:       call VecGetArray(local,localptr,ldx,ierr)

444: !  Extract work vector

446:       call VecGetArray(localwork,copyptr,cdx,ierr)

448: !   Update Locally - Make array of new values
449: !   Note: For the first and last entry I copy the value
450: !   if this is an interior node it is irrelevant

452:       sc = 1.0/(h*h)
453:       call VecGetLocalSize(local,localsize,ierr)
454:       copyptr(1+cdx) = localptr(1+ldx)
455:       do 10, i=1,localsize-2
456:         copyptr(i+1+cdx) = sc * (localptr(i+2+ldx) + localptr(i+ldx) -  &
457:      &                     2.0*localptr(i+1+ldx))
458:  10   continue
459:       copyptr(localsize-1+1+cdx) = localptr(localsize-1+1+ldx)
460:       call VecRestoreArray(localwork,copyptr,cdx,ierr)
461:       call VecRestoreArray(local,localptr,ldx,ierr)
462:       call VecDestroy(local,ierr)

464: !   Local to Global

466:       call DMLocalToGlobalBegin(da,localwork,INSERT_VALUES,globalout,    &
467:      &             ierr)
468:       call DMLocalToGlobalEnd(da,localwork,INSERT_VALUES,globalout,      &
469:      &            ierr)
470:       return
471:       end

473: !  ---------------------------------------------------------------------

475:       subroutine RHSMatrixHeat(ts,t,AA,BB,str,ctx,ierr)
476:       implicit none
477: #include <finclude/petscsys.h>
478: #include <finclude/petscvec.h>
479: #include <finclude/petscmat.h>
480: #include <finclude/petscdmda.h>
481: #include <finclude/petscpc.h>
482: #include <finclude/petscksp.h>
483: #include <finclude/petscsnes.h>
484: #include <finclude/petscts.h>
485: #include <finclude/petscviewer.h>

487: #include "ex1f.h"
488:       Mat              AA,BB
489:       double precision t
490:       TS               ts
491:       MatStructure     str
492:       PetscObject          ctx
493:       PetscErrorCode ierr

495:       Mat              A
496:       PetscInt    i,mstart(1),mend(1),idx(3)
497:       PetscMPIInt rank,size
498:       PetscScalar v(3),stwo,sone
499:       PetscInt    i1,i3

501:       i1 = 1
502:       i3 = 3
503:       A    = AA
504:       stwo = -2./(h*h)
505:       sone = -.5*stwo
506:       str  = SAME_NONZERO_PATTERN

508:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
509:       call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)

511:       call MatGetOwnershipRange(A,mstart,mend,ierr)
512:       if (mstart(1) .eq. 0) then
513:         v(1) = 1.0
514:         call MatSetValues(A,i1,mstart(1),i1,mstart(1),v,INSERT_VALUES,  &
515:      &       ierr)
516:         mstart(1) = mstart(1) + 1
517:       endif
518:       if (mend(1) .eq. M) then
519:         mend(1) = mend(1) - 1
520:         v(1) = 1.0
521:         call MatSetValues(A,i1,mend,i1,mend,v,INSERT_VALUES,ierr)
522:       endif

524: !
525: !     Construct matrice one row at a time
526: !
527:       v(1) = sone
528:       v(2) = stwo
529:       v(3) = sone
530:       do 10, i=mstart(1),mend(1)-1
531:         idx(1) = i-1
532:         idx(2) = i
533:         idx(3) = i+1
534:         call MatSetValues(A,i1,i,i3,idx,v,INSERT_VALUES,ierr)
535:  10   continue

537:       call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
538:       call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)
539:       return
540:       end

542: ! --------------------------------------------------------------------------------------

544:       subroutine RHSJacobianHeat(ts,t,x,AA,BB,str,ctx,ierr)
545:       implicit none
546: #include <finclude/petscsys.h>
547: #include <finclude/petscvec.h>
548: #include <finclude/petscmat.h>
549: #include <finclude/petscdmda.h>
550: #include <finclude/petscpc.h>
551: #include <finclude/petscksp.h>
552: #include <finclude/petscsnes.h>
553: #include <finclude/petscts.h>
554: #include <finclude/petscviewer.h>
555:       TS               ts
556:       double precision t
557:       Vec              x
558:       Mat              AA,BB
559:       MatStructure     str
560:       PetscObject      ctx
561:       PetscErrorCode ierr

563:       call RHSMatrixHeat(ts,t,AA,BB,str,ctx,ierr)
564:       return
565:       end