Actual source code: ex1f.F

  1: !
  2: !   Solves the time dependent Bratu problem using pseudo-timestepping
  3: !
  4: !   Concepts: TS^pseudo-timestepping
  5: !   Concepts: pseudo-timestepping
  6: !   Concepts: nonlinear problems
  7: !   Processors: 1
  8: !
  9: !   This code demonstrates how one may solve a nonlinear problem
 10: !   with pseudo-timestepping. In this simple example, the pseudo-timestep
 11: !   is the same for all grid points, i.e., this is equivalent to using
 12: !   the backward Euler method with a variable timestep.
 13: !
 14: !   Note: This example does not require pseudo-timestepping since it
 15: !   is an easy nonlinear problem, but it is included to demonstrate how
 16: !   the pseudo-timestepping may be done.
 17: !
 18: !   See snes/examples/tutorials/ex4.c[ex4f.F] and
 19: !   snes/examples/tutorials/ex5.c[ex5f.F] where the problem is described
 20: !   and solved using the method of Newton alone.
 21: !
 22: !   Include "petscts.h"    to use the PETSc timestepping routines,
 23: !           "petsc.h" for basic PETSc operation,
 24: !           "petscmat.h"   for matrix operations,
 25: !           "petscpc.h"    for preconditions, and
 26: !           "petscvec.h"   for vector operations.
 27: !
 28: !23456789012345678901234567890123456789012345678901234567890123456789012
 29:       program main
 30:       implicit none
 31:  #include finclude/petsc.h
 32:  #include finclude/petscvec.h
 33:  #include finclude/petscmat.h
 34:  #include finclude/petscpc.h
 35:  #include finclude/petscts.h
 36: !
 37: !  Create an application context to contain data needed by the
 38: !  application-provided call-back routines, FormJacobian() and
 39: !  FormFunction(). We use a double precision array with three
 40: !  entries indexed by param, lmx, lmy.
 41: !
 42:       double precision user(3)
 43:       PetscInt          param,lmx,lmy,i5
 44:       parameter (param = 1,lmx = 2,lmy = 3)
 45: !
 46: !   User-defined routines
 47: !
 48:       external FormJacobian,FormFunction
 49: !
 50: !   Data for problem
 51: !
 52:       TS                ts
 53:       Vec               x,r
 54:       Mat               J
 55:       PetscInt           its,N,i1000
 56:       PetscTruth flg
 57:       PetscErrorCode      ierr
 58:       double precision  param_max,param_min,dt,tmax,zero
 59:       double precision  ftime

 61:       i5 = 5
 62:       param_max = 6.81
 63:       param_min = 0

 65:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 66:       user(lmx)        = 4
 67:       user(lmy)        = 4
 68:       user(param)      = 6.0
 69: 
 70: !
 71: !     Allow user to set the grid dimensions and nonlinearity parameter at run-time
 72: !
 73:       call PetscOptionsGetReal(PETSC_NULL_CHARACTER,'-mx',user(lmx),    &
 74:      &     flg,ierr)
 75:       call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-my',user(lmy),     &
 76:      &     flg,ierr)
 77:       call PetscOptionsGetReal(PETSC_NULL_CHARACTER,'-param',           &
 78:      &     user(param),flg,ierr)
 79:       if (user(param) .ge. param_max .or.                               &
 80:      &                                user(param) .le. param_min) then
 81:         print*,'Parameter is out of range'
 82:       endif
 83:       if (user(lmx) .gt. user(lmy)) then
 84:         dt = .5/user(lmx)
 85:       else
 86:         dt = .5/user(lmy)
 87:       endif
 88:       call PetscOptionsGetReal(PETSC_NULL_CHARACTER,'-dt',dt,flg,ierr)
 89:       N          = user(lmx)*user(lmy)
 90: 
 91: !
 92: !      Create vectors to hold the solution and function value
 93: !
 94:       call VecCreateSeq(PETSC_COMM_SELF,N,x,ierr)
 95:       call VecDuplicate(x,r,ierr)

 97: !
 98: !    Create matrix to hold Jacobian. Preallocate 5 nonzeros per row
 99: !    in the sparse matrix. Note that this is not the optimal strategy see
100: !    the Performance chapter of the users manual for information on
101: !    preallocating memory in sparse matrices.
102: !
103:       i5 = 5
104:       call MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,i5,PETSC_NULL_INTEGER,    &
105:      &     J,ierr)

107: !
108: !     Create timestepper context
109: !

111:       call TSCreate(PETSC_COMM_WORLD,ts,ierr)
112:       call TSSetProblemType(ts,TS_NONLINEAR,ierr)

114: !
115: !     Tell the timestepper context where to compute solutions
116: !

118:       call TSSetSolution(ts,x,ierr)

120: !
121: !    Provide the call-back for the nonlinear function we are
122: !    evaluating. Thus whenever the timestepping routines need the
123: !    function they will call this routine. Note the final argument
124: !    is the application context used by the call-back functions.
125: !

127:       call TSSetRHSFunction(ts,FormFunction,user,ierr)

129: !
130: !     Set the Jacobian matrix and the function used to compute
131: !     Jacobians.
132: !

134:       call TSSetRHSJacobian(ts,J,J,FormJacobian,user,ierr)

136: !
137: !       For the initial guess for the problem
138: !

140:       call FormInitialGuess(x,user,ierr)

142: !
143: !       This indicates that we are using pseudo timestepping to
144: !     find a steady state solution to the nonlinear problem.
145: !

147:       call TSSetType(ts,TS_PSEUDO,ierr)

149: !
150: !       Set the initial time to start at (this is arbitrary for
151: !     steady state problems and the initial timestep given above
152: !

154:       zero = 0.0
155:       call TSSetInitialTimeStep(ts,zero,dt,ierr)

157: !
158: !      Set a large number of timesteps and final duration time
159: !     to insure convergence to steady state.
160: !
161:       i1000 = 1000
162:       tmax  = 1.e12
163:       call TSSetDuration(ts,i1000,tmax,ierr)

165: !
166: !      Set any additional options from the options database. This
167: !     includes all options for the nonlinear and linear solvers used
168: !     internally the the timestepping routines.
169: !

171:       call TSSetFromOptions(ts,ierr)

173:       call TSSetUp(ts,ierr)

175: !
176: !      Perform the solve. This is where the timestepping takes place.
177: !
178: 
179:       call TSStep(ts,its,ftime,ierr)
180: 
181:       write(6,100) its,ftime
182:  100  format('Number of pseudo time-steps ',i5,' final time ',1pe8.2)

184: !
185: !     Free the data structures constructed above
186: !

188:       call VecDestroy(x,ierr)
189:       call VecDestroy(r,ierr)
190:       call MatDestroy(J,ierr)
191:       call TSDestroy(ts,ierr)
192:       call PetscFinalize(ierr)
193:       end

195: !
196: !  --------------------  Form initial approximation -----------------
197: !
198:       subroutine FormInitialGuess(X,user,ierr)
199:       implicit none
200:  #include finclude/petsc.h
201:  #include finclude/petscvec.h
202:  #include finclude/petscmat.h
203:  #include finclude/petscpc.h
204:  #include finclude/petscts.h
205:       Vec              X
206:       double precision user(3)
207:       PetscInt  i,j,row,mx,my
208:       PetscErrorCode ierr
209:       PetscOffset      xidx
210:       double precision two,one,lambda
211:       double precision temp1,temp,hx,hy,hxdhy,hydhx
212:       double precision sc
213:       PetscScalar      xx(1)
214:       PetscInt          param,lmx,lmy
215:       parameter (param = 1,lmx = 2,lmy = 3)

217:       two = 2.0
218:       one = 1.0

220:       mx     = user(lmx)
221:       my     = user(lmy)
222:       lambda = user(param)

224:       hy    = one / (my-1)
225:       hx    = one / (mx-1)
226:       sc    = hx*hy
227:       hxdhy = hx/hy
228:       hydhx = hy/hx

230:       call VecGetArray(X,xx,xidx,ierr)
231:       temp1 = lambda/(lambda + one)
232:       do 10, j=1,my
233:         temp = dble(min(j-1,my-j))*hy
234:         do 20 i=1,mx
235:           row = i + (j-1)*mx
236:           if (i .eq. 1 .or. j .eq. 1 .or.                               &
237:      &        i .eq. mx .or. j .eq. my) then
238:             xx(row+xidx) = 0.0
239:           else
240:             xx(row+xidx) =                                              &
241:      &        temp1*sqrt(min(dble(min(i-1,mx-i))*hx,temp))
242:           endif
243:  20     continue
244:  10   continue
245:       call VecRestoreArray(X,xx,xidx,ierr)
246:       return
247:       end
248: !
249: !  --------------------  Evaluate Function F(x) ---------------------
250: !
251:       subroutine FormFunction(ts,t,X,F,user,ierr)
252:       implicit none
253:  #include finclude/petsc.h
254:  #include finclude/petscvec.h
255:  #include finclude/petscmat.h
256:  #include finclude/petscpc.h
257:  #include finclude/petscts.h
258:       TS       ts
259:       double precision  t
260:       Vec               X,F
261:       double precision  user(3)
262:       PetscErrorCode     ierr
263:       PetscInt         i,j,row,mx,my
264:       PetscOffset       xidx,fidx
265:       double precision  two,one,lambda
266:       double precision  hx,hy,hxdhy,hydhx
267:       PetscScalar  ut,ub,ul,ur,u
268:       PetscScalar  uxx,uyy,sc
269:       PetscScalar  xx(1),ff(1)
270:       PetscInt     param,lmx,lmy
271:       parameter (param = 1,lmx = 2,lmy = 3)

273:       two = 2.0
274:       one = 1.0

276:       mx     = user(lmx)
277:       my     = user(lmy)
278:       lambda = user(param)

280:       hx    = 1.0 / dble(mx-1)
281:       hy    = 1.0 / dble(my-1)
282:       sc    = hx*hy
283:       hxdhy = hx/hy
284:       hydhx = hy/hx

286:       call VecGetArray(X,xx,xidx,ierr)
287:       call VecGetArray(F,ff,fidx,ierr)
288:       do 10 j=1,my
289:         do 20 i=1,mx
290:           row = i + (j-1)*mx
291:           if (i .eq. 1 .or. j .eq. 1 .or.                               &
292:      &        i .eq. mx .or. j .eq. my) then
293:             ff(row+fidx) = xx(row+xidx)
294:           else
295:             u            = xx(row + xidx)
296:             ub           = xx(row - mx + xidx)
297:             ul           = xx(row - 1 + xidx)
298:             ut           = xx(row + mx + xidx)
299:             ur           = xx(row + 1 + xidx)
300:             uxx          = (-ur + two*u - ul)*hydhx
301:             uyy          = (-ut + two*u - ub)*hxdhy
302:             ff(row+fidx) = -uxx - uyy + sc*lambda*exp(u)
303:             u =  -uxx - uyy + sc*lambda*exp(u)
304:          endif
305:  20   continue
306:  10   continue

308:       call VecRestoreArray(X,xx,xidx,ierr)
309:       call VecRestoreArray(F,ff,fidx,ierr)
310:       return
311:       end
312: !
313: !  --------------------  Evaluate Jacobian of F(x) --------------------
314: !
315:       subroutine FormJacobian(ts,ctime,X,JJ,B,flag,user,ierr)
316:       implicit none
317:  #include finclude/petsc.h
318:  #include finclude/petscvec.h
319:  #include finclude/petscmat.h
320:  #include finclude/petscpc.h
321:  #include finclude/petscts.h
322:       TS               ts
323:       Vec              X
324:       Mat              JJ,B
325:       MatStructure     flag
326:       double precision user(3),ctime
327:       PetscErrorCode   ierr
328:       Mat              jac
329:       PetscOffset xidx
330:       PetscInt    i,j,row(1),mx,my
331:       PetscInt    col(5),i1,i5
332:       PetscScalar two,one,lambda
333:       PetscScalar v(5),sc,xx(1)
334:       double precision hx,hy,hxdhy,hydhx

336:       PetscInt  param,lmx,lmy
337:       parameter (param = 1,lmx = 2,lmy = 3)

339:       i1 = 1
340:       i5 = 5
341:       jac = B
342:       two = 2.0
343:       one = 1.0

345:       mx     = user(lmx)
346:       my     = user(lmy)
347:       lambda = user(param)

349:       hx    = 1.0 / dble(mx-1)
350:       hy    = 1.0 / dble(my-1)
351:       sc    = hx*hy
352:       hxdhy = hx/hy
353:       hydhx = hy/hx

355:       call VecGetArray(X,xx,xidx,ierr)
356:       do 10 j=1,my
357:         do 20 i=1,mx
358: !
359: !      When inserting into PETSc matrices, indices start at 0
360: !
361:           row(1) = i - 1 + (j-1)*mx
362:           if (i .eq. 1 .or. j .eq. 1 .or.                               &
363:      &        i .eq. mx .or. j .eq. my) then
364:             call MatSetValues(jac,i1,row,i1,row,one,INSERT_VALUES,ierr)
365:           else
366:             v(1)   = hxdhy
367:             col(1) = row(1) - mx
368:             v(2)   = hydhx
369:             col(2) = row(1) - 1
370:             v(3)   = -two*(hydhx+hxdhy)+sc*lambda*exp(xx(row(1)+1+xidx))
371:             col(3) = row(1)
372:             v(4)   = hydhx
373:             col(4) = row(1) + 1
374:             v(5)   = hxdhy
375:             col(5) = row(1) + mx
376:             call MatSetValues(jac,i1,row,i5,col,v,INSERT_VALUES,ierr)
377:           endif
378:  20     continue
379:  10   continue
380:       call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
381:       call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
382:       call VecRestoreArray(X,xx,xidx,ierr)
383:       flag = SAME_NONZERO_PATTERN
384:       return
385:       end