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:       program main
 11:       implicit none
 12:  #include finclude/petsc.h
 13:  #include finclude/petscvec.h
 14:  #include finclude/petscmat.h
 15:  #include finclude/petscda.h
 16:  #include finclude/petscsys.h
 17:  #include finclude/petscpc.h
 18:  #include finclude/petscksp.h
 19:  #include finclude/petscsnes.h
 20:  #include finclude/petscts.h
 21:  #include finclude/petscdraw.h
 22:  #include finclude/petscviewer.h

 24: #include "ex1f.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:       PetscErrorCode  ierr
 32:       PetscInt time_steps,steps
 33:       PetscMPIInt size
 34:       integer problem
 35:       Vec              local,global
 36:       double precision dt,ftime,zero,tmax
 37:       TS               ts
 38:       Mat              A
 39:       MatStructure     A_structure
 40:       TSProblemType    tsproblem
 41:       PetscDraw        draw
 42:       PetscViewer      viewer
 43:       character*(60)   type,tsinfo
 44:       character*(5)    etype
 45:       PetscInt         i1,i60
 46:       PetscTruth flg

 48:       external Monitor,RHSFunctionHeat,RHSMatrixFree,Initial
 49:       external RHSMatrixHeat,RHSJacobianHeat

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

 62:       M = 60
 63:       call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-M',M,flg,ierr)
 64:       call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-time',time_steps,   &
 65:      &     flg,ierr)

 67:       call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-nox',flg,ierr)
 68:       if (flg .eq. 1) then
 69:         nox = 1
 70:       else
 71:         nox = 0
 72:       endif
 73:       nrm_2   = 0.0
 74:       nrm_max = 0.0

 76: !   Set up the ghost point communication pattern

 78:       call DACreate1d(PETSC_COMM_WORLD,DA_NONPERIODIC,M,i1,i1,            &
 79:      &     PETSC_NULL_INTEGER,da,ierr)
 80:       call DACreateGlobalVector(da,global,ierr)
 81:       call VecGetLocalSize(global,m,ierr)
 82:       call DACreateLocalVector(da,local,ierr)

 84: !   Set up display to show wave graph

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

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

 97:       call VecDuplicate(local,localwork,ierr)

 99: !   make work array for storing exact solution

101:       call VecDuplicate(global,csolution,ierr)

103:       h = 1.0/(M-1.0)

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

147:       call TSCreate(PETSC_COMM_WORLD,ts,ierr)
148:       call TSSetProblemType(ts,tsproblem,ierr)
149:       call TSSetMonitor(ts,Monitor,PETSC_NULL_OBJECT,                   &
150:      &                  PETSC_NULL_FUNCTION, ierr)

152:       dt = h*h/2.01

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

211:       call TSSetFromOptions(ts,ierr)

213:       call TSSetInitialTimeStep(ts,zero,dt,ierr)
214:       tmax = 100.0
215:       call TSSetDuration(ts,time_steps,tmax,ierr)
216:       call TSSetSolution(ts,global,ierr)

218:       call TSSetUp(ts,ierr)
219:       call TSStep(ts,steps,ftime,ierr)
220:       call PetscViewerStringOpen(PETSC_COMM_WORLD,tsinfo,i60,viewer,         &
221:      &                           ierr)
222:       call TSView(ts,viewer,ierr)
223:       call TSGetType(ts,type,ierr)

225:       call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-test',flg,ierr)
226:       if (flg .eq. 1) then
227: !
228: !         copy type to string of length 5 to ensure equality test
229: !         is done correctly
230: !
231:         call PetscStrncpy(etype,type,5,ierr)
232:         if (etype .eq. TS_EULER) then
233:           if (abs(nrm_2/steps - 0.00257244) .gt. 1.e-4) then
234:             print*,'Error in Euler method: 2-norm ',nrm_2/steps,        &
235:      &            ' expecting: ',0.00257244
236:           endif
237:         else
238:           if (abs(nrm_2/steps - 0.00506174) .gt. 1.e-4) then
239:             print*,'Error in ',tsinfo,': 2-norm ',nrm_2/steps,          &
240:      &             ' expecting: ',0.00506174
241:           endif
242:         endif
243:       else
244:         print*,size,' Procs Avg. error 2 norm ',nrm_2/steps,            &
245:      &              nrm_max/steps,tsinfo
246:       endif

248:       call PetscViewerDestroy(viewer,ierr)
249:       call TSDestroy(ts,ierr)
250:       call PetscViewerDestroy(viewer1,ierr)
251:       call PetscViewerDestroy(viewer2,ierr)
252:       call VecDestroy(localwork,ierr)
253:       call VecDestroy(csolution,ierr)
254:       call VecDestroy(local,ierr)
255:       call VecDestroy(global,ierr)
256:       call DADestroy(da,ierr)
257:       if (A .ne. 0) then
258:         call MatDestroy(A,ierr)
259:       endif

261:       call PetscFinalize(ierr)
262:       end

264: !  -------------------------------------------------------------------
265: 
266:       subroutine Initial(global,ctx,ierr)
267:       implicit none
268:  #include finclude/petsc.h
269:  #include finclude/petscvec.h
270:  #include finclude/petscmat.h
271:  #include finclude/petscda.h
272:  #include finclude/petscsys.h
273:  #include finclude/petscpc.h
274:  #include finclude/petscksp.h
275:  #include finclude/petscsnes.h
276:  #include finclude/petscts.h
277:  #include finclude/petscviewer.h

279: #include "ex1f.h"

281:       Vec         global
282:       PetscObject    ctx

284:       PetscScalar localptr(1)
285:       PetscInt     i,mybase,myend
286:       PetscErrorCode ierr
287:       PetscOffset idx

289: !   determine starting point of each processor

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

293: !   Initialize the array

295:       call VecGetArray(global,localptr,idx,ierr)
296:       do 10, i=mybase,myend-1
297:         localptr(i-mybase+1+idx) = sin(PETSC_PI*i*6.*h) +               &
298:      &                             3.*sin(PETSC_PI*i*2.*h)
299:  10   continue
300:       call VecRestoreArray(global,localptr,idx,ierr)
301:       return
302:       end

304: ! ------------------------------------------------------------------------------
305: !
306: !       Exact solution
307: !
308:       subroutine Solution(t,sol,ctx)
309:       implicit none
310:  #include finclude/petsc.h
311:  #include finclude/petscvec.h
312:  #include finclude/petscmat.h
313:  #include finclude/petscda.h
314:  #include finclude/petscsys.h
315:  #include finclude/petscpc.h
316:  #include finclude/petscksp.h
317:  #include finclude/petscsnes.h
318:  #include finclude/petscts.h
319:  #include finclude/petscviewer.h

321: #include "ex1f.h"

323:       double precision  t
324:       Vec               sol
325:       integer           ctx

327:       PetscScalar       localptr(1),ex1,ex2,sc1,sc2
328:       PetscInt           i,mybase,myend
329:       PetscErrorCode ierr
330:       PetscOffset       idx

332: !   determine starting point of each processor

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

336:       ex1 = exp(-36.*PETSC_PI*PETSC_PI*t)
337:       ex2 = exp(-4.*PETSC_PI*PETSC_PI*t)
338:       sc1 = PETSC_PI*6.*h
339:       sc2 = PETSC_PI*2.*h
340:       call VecGetArray(csolution,localptr,idx,ierr)
341:       do 10, i=mybase,myend-1
342:         localptr(i-mybase+1+idx) = sin(i*sc1)*ex1 + 3.*sin(i*sc2)*ex2
343:  10   continue
344:       call VecRestoreArray(csolution,localptr,idx,ierr)
345:       return
346:       end


349: ! -----------------------------------------------------------------------------------

351:       subroutine Monitor(ts,step,time,global,ctx,ierr)
352:       implicit none
353:  #include finclude/petsc.h
354:  #include finclude/petscvec.h
355:  #include finclude/petscmat.h
356:  #include finclude/petscda.h
357:  #include finclude/petscsys.h
358:  #include finclude/petscpc.h
359:  #include finclude/petscksp.h
360:  #include finclude/petscsnes.h
361:  #include finclude/petscts.h
362:  #include finclude/petscdraw.h
363:  #include finclude/petscviewer.h

365: #include "ex1f.h"
366:       TS                ts
367:       PetscInt           step
368:       integer ctx
369:       PetscErrorCode ierr
370:       double precision  time,lnrm_2,lnrm_max
371:       Vec               global
372:       PetscScalar       mone

374:       call VecView(global,viewer1,ierr)

376:       call Solution(time,csolution,ctx)
377:       mone = -1.0
378:       call VecAXPY(csolution,mone,global,ierr)
379:       call VecNorm(csolution,NORM_2,lnrm_2,ierr)
380:       lnrm_2 = sqrt(h)*lnrm_2
381:       call VecNorm(csolution,NORM_MAX,lnrm_max,ierr)

383:       if (nox .eq. 0) then
384:         print*,'timestep ',step,' time ',time,' norm of error ',        &
385:      &         lnrm_2,lnrm_max
386:       endif

388:       nrm_2   = nrm_2 + lnrm_2
389:       nrm_max = nrm_max +  lnrm_max
390:       call VecView(csolution,viewer2,ierr)

392:       return
393:       end

395: !  -----------------------------------------------------------------------

397:       subroutine RHSMatrixFree(mat,x,y,ierr)
398:       implicit none
399:  #include finclude/petsc.h
400:  #include finclude/petscvec.h
401:  #include finclude/petscmat.h
402:  #include finclude/petscda.h
403:  #include finclude/petscsys.h
404:  #include finclude/petscpc.h
405:  #include finclude/petscksp.h
406:  #include finclude/petscsnes.h
407:  #include finclude/petscts.h
408:  #include finclude/petscviewer.h
409:       Mat              mat
410:       Vec              x,y
411:       PetscErrorCode          ierr
412:       double precision zero
413:       TS               ts0

415:       zero = 0.0

417:       ts0 = PETSC_NULL_OBJECT

419:       call RHSFunctionHeat(ts0,zero,x,y,PETSC_NULL_OBJECT,ierr)
420:       return
421:       end

423: ! -------------------------------------------------------------------------

425:       subroutine RHSFunctionHeat(ts,t,globalin,globalout,ctx,ierr)
426:       implicit none
427:  #include finclude/petsc.h
428:  #include finclude/petscvec.h
429:  #include finclude/petscmat.h
430:  #include finclude/petscda.h
431:  #include finclude/petscsys.h
432:  #include finclude/petscpc.h
433:  #include finclude/petscksp.h
434:  #include finclude/petscsnes.h
435:  #include finclude/petscts.h
436:  #include finclude/petscviewer.h

438: #include "ex1f.h"
439:       TS               ts
440:       double precision t
441:       Vec              globalin,globalout
442:       PetscObject      ctx
443:       Vec              local
444:       PetscErrorCode  ierr
445:       PetscInt i,localsize
446:       PetscOffset      ldx,cdx
447:       PetscScalar      copyptr(1),localptr(1),sc

449: !  Extract local array

451:       call DACreateLocalVector(da,local,ierr)
452:       call DAGlobalToLocalBegin(da,globalin,INSERT_VALUES,local,ierr)
453:       call DAGlobalToLocalEnd(da,globalin,INSERT_VALUES,local,ierr)
454:       call VecGetArray(local,localptr,ldx,ierr)

456: !  Extract work vector

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

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

464:       sc = 1.0/(h*h)
465:       call VecGetLocalSize(local,localsize,ierr)
466:       copyptr(1+cdx) = localptr(1+ldx)
467:       do 10, i=1,localsize-2
468:         copyptr(i+1+cdx) = sc * (localptr(i+2+ldx) + localptr(i+ldx) -  &
469:      &                     2.0*localptr(i+1+ldx))
470:  10   continue
471:       copyptr(localsize-1+1+cdx) = localptr(localsize-1+1+ldx)
472:       call VecRestoreArray(localwork,copyptr,cdx,ierr)
473:       call VecRestoreArray(local,localptr,ldx,ierr)
474:       call VecDestroy(local,ierr)

476: !   Local to Global

478:       call DALocalToGlobal(da,localwork,INSERT_VALUES,globalout,ierr)
479:       return
480:       end

482: !  ---------------------------------------------------------------------

484:       subroutine RHSMatrixHeat(ts,t,AA,BB,str,ctx,ierr)
485:       implicit none
486:  #include finclude/petsc.h
487:  #include finclude/petscvec.h
488:  #include finclude/petscmat.h
489:  #include finclude/petscda.h
490:  #include finclude/petscsys.h
491:  #include finclude/petscpc.h
492:  #include finclude/petscksp.h
493:  #include finclude/petscsnes.h
494:  #include finclude/petscts.h
495:  #include finclude/petscviewer.h

497: #include "ex1f.h"
498:       Mat              AA,BB
499:       double precision t
500:       TS               ts
501:       MatStructure     str
502:       PetscObject          ctx
503:       PetscErrorCode ierr

505:       Mat              A
506:       PetscInt       i,mstart,mend,idx(3)
507:       PetscMPIInt rank,size
508:       PetscScalar      v(3),stwo,sone
509:       PetscInt          i1,i3

511:       i1 = 1
512:       i3 = 3
513:       A    = AA
514:       stwo = -2./(h*h)
515:       sone = -.5*stwo
516:       str  = SAME_NONZERO_PATTERN

518:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
519:       call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)

521:       call MatGetOwnershipRange(A,mstart,mend,ierr)
522:       if (mstart .eq. 0) then
523:         v(1) = 1.0
524:         call MatSetValues(A,i1,mstart,i1,mstart,v,INSERT_VALUES,ierr)
525:         mstart = mstart + 1
526:       endif
527:       if (mend .eq. M) then
528:         mend = mend - 1
529:         v(1) = 1.0
530:         call MatSetValues(A,i1,mend,i1,mend,v,INSERT_VALUES,ierr)
531:       endif

533: !
534: !     Construct matrice one row at a time
535: !
536:       v(1) = sone
537:       v(2) = stwo
538:       v(3) = sone
539:       do 10, i=mstart,mend-1
540:         idx(1) = i-1
541:         idx(2) = i
542:         idx(3) = i+1
543:         call MatSetValues(A,i1,i,i3,idx,v,INSERT_VALUES,ierr)
544:  10   continue

546:       call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
547:       call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)
548:       return
549:       end

551: ! --------------------------------------------------------------------------------------

553:       subroutine RHSJacobianHeat(ts,t,x,AA,BB,str,ctx,ierr)
554:       implicit none
555:  #include finclude/petsc.h
556:  #include finclude/petscvec.h
557:  #include finclude/petscmat.h
558:  #include finclude/petscda.h
559:  #include finclude/petscsys.h
560:  #include finclude/petscpc.h
561:  #include finclude/petscksp.h
562:  #include finclude/petscsnes.h
563:  #include finclude/petscts.h
564:  #include finclude/petscviewer.h
565:       TS               ts
566:       double precision t
567:       Vec              x
568:       Mat              AA,BB
569:       MatStructure     str
570:       PetscObject      ctx
571:       PetscErrorCode ierr

573:       call RHSMatrixHeat(ts,t,AA,BB,str,ctx,ierr)
574:       return
575:       end