Actual source code: ex1f.F
1: !
2: ! "$Id: ex1f.F,v 1.26 2001/01/19 23:22:15 balay Exp $"
3: !
4: ! Solves the time dependent Bratu problem using pseudo-timestepping
5: !
6: ! Concepts: TS^pseudo-timestepping
7: ! Concepts: pseudo-timestepping
8: ! Concepts: nonlinear problems
9: ! Processors: 1
10: !
11: ! This code demonstrates how one may solve a nonlinear problem
12: ! with pseudo-timestepping. In this simple example, the pseudo-timestep
13: ! is the same for all grid points, i.e., this is equivalent to using
14: ! the backward Euler method with a variable timestep.
15: !
16: ! Note: This example does not require pseudo-timestepping since it
17: ! is an easy nonlinear problem, but it is included to demonstrate how
18: ! the pseudo-timestepping may be done.
19: !
20: ! See snes/examples/tutorials/ex4.c[ex4f.F] and
21: ! snes/examples/tutorials/ex5.c[ex5f.F] where the problem is described
22: ! and solved using the method of Newton alone.
23: !
24: ! Include "petscts.h" to use the PETSc timestepping routines,
25: ! "petsc.h" for basic PETSc operation,
26: ! "petscmat.h" for matrix operations,
27: ! "petscpc.h" for preconditions, and
28: ! "petscvec.h" for vector operations.
29: !
30: !23456789012345678901234567890123456789012345678901234567890123456789012
31: program main
32: implicit none
33: #include finclude/petsc.h
34: #include finclude/petscvec.h
35: #include finclude/petscmat.h
36: #include finclude/petscpc.h
37: #include finclude/petscts.h
38: !
39: ! Create an application context to contain data needed by the
40: ! application-provided call-back routines, FormJacobian() and
41: ! FormFunction(). We use a double precision array with three
42: ! entries indexed by param, lmx, lmy.
43: !
44: double precision user(3)
45: integer param,lmx,lmy
46: parameter (param = 1,lmx = 2,lmy = 3)
47: !
48: ! User-defined routines
49: !
50: external FormJacobian,FormFunction
51: !
52: ! Data for problem
53: !
54: TS ts
55: Vec x,r
56: Mat J
57: integer its
58: integer ierr,N,flg
59: double precision param_max,param_min,dt,tmax,zero
60: double precision ftime
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 PetscOptionsGetDouble(PETSC_NULL_CHARACTER,'-mx',user(lmx), &
74: & flg,ierr)
75: call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-my',user(lmy), &
76: & flg,ierr)
77: call PetscOptionsGetDouble(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 PetscOptionsGetDouble(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: !
104: call MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,5,PETSC_NULL_INTEGER, &
105: & J,ierr)
107: !
108: ! Create timestepper context
109: !
111: call TSCreate(PETSC_COMM_WORLD,TS_NONLINEAR,ts,ierr)
113: !
114: ! Tell the timestepper context where to compute solutions
115: !
117: call TSSetSolution(ts,x,ierr)
119: !
120: ! Provide the call-back for the nonlinear function we are
121: ! evaluating. Thus whenever the timestepping routines need the
122: ! function they will call this routine. Note the final argument
123: ! is the application context used by the call-back functions.
124: !
126: call TSSetRHSFunction(ts,FormFunction,user,ierr)
128: !
129: ! Set the Jacobian matrix and the function used to compute
130: ! Jacobians.
131: !
133: call TSSetRHSJacobian(ts,J,J,FormJacobian,user,ierr)
135: !
136: ! For the initial guess for the problem
137: !
139: call FormInitialGuess(x,user,ierr)
141: !
142: ! This indicates that we are using pseudo timestepping to
143: ! find a steady state solution to the nonlinear problem.
144: !
146: call TSSetType(ts,TS_PSEUDO,ierr)
148: !
149: ! Set the initial time to start at (this is arbitrary for
150: ! steady state problems and the initial timestep given above
151: !
153: zero = 0.0
154: call TSSetInitialTimeStep(ts,zero,dt,ierr)
156: !
157: ! Set a large number of timesteps and final duration time
158: ! to insure convergence to steady state.
159: !
160: tmax = 1.e12
161: call TSSetDuration(ts,1000,tmax,ierr)
163: !
164: ! Set any additional options from the options database. This
165: ! includes all options for the nonlinear and linear solvers used
166: ! internally the the timestepping routines.
167: !
169: call TSSetFromOptions(ts,ierr)
171: call TSSetUp(ts,ierr)
173: !
174: ! Perform the solve. This is where the timestepping takes place.
175: !
176:
177: call TSStep(ts,its,ftime,ierr)
178:
179: write(6,100) its,ftime
180: 100 format('Number of pseudo time-steps ',i5,' final time ',1pe8.2)
182: !
183: ! Free the data structures constructed above
184: !
186: call VecDestroy(x,ierr)
187: call VecDestroy(r,ierr)
188: call MatDestroy(J,ierr)
189: call TSDestroy(ts,ierr)
190: call PetscFinalize(ierr)
191: end
193: !
194: ! -------------------- Form initial approximation -----------------
195: !
196: subroutine FormInitialGuess(X,user,ierr)
197: implicit none
198: #include finclude/petsc.h
199: #include finclude/petscvec.h
200: #include finclude/petscmat.h
201: #include finclude/petscpc.h
202: #include finclude/petscts.h
203: Vec X
204: double precision user(3)
205: integer i,j,row,mx,my,ierr
206: PetscOffset xidx
207: double precision two,one,lambda
208: double precision temp1,temp,hx,hy,hxdhy,hydhx
209: double precision sc
210: Scalar xx(1)
211: integer param,lmx,lmy
212: parameter (param = 1,lmx = 2,lmy = 3)
214: two = 2.0
215: one = 1.0
217: mx = user(lmx)
218: my = user(lmy)
219: lambda = user(param)
221: hy = one / (my-1)
222: hx = one / (mx-1)
223: sc = hx*hy
224: hxdhy = hx/hy
225: hydhx = hy/hx
227: call VecGetArray(X,xx,xidx,ierr)
228: temp1 = lambda/(lambda + one)
229: do 10, j=1,my
230: temp = dble(min(j-1,my-j))*hy
231: do 20 i=1,mx
232: row = i + (j-1)*mx
233: if (i .eq. 1 .or. j .eq. 1 .or. &
234: & i .eq. mx .or. j .eq. my) then
235: xx(row+xidx) = 0.0
236: else
237: xx(row+xidx) = &
238: & temp1*sqrt(min(dble(min(i-1,mx-i))*hx,temp))
239: endif
240: 20 continue
241: 10 continue
242: call VecRestoreArray(X,xx,xidx,ierr)
243: return
244: end
245: !
246: ! -------------------- Evaluate Function F(x) ---------------------
247: !
248: subroutine FormFunction(ts,t,X,F,user,ierr)
249: implicit none
250: #include finclude/petsc.h
251: #include finclude/petscvec.h
252: #include finclude/petscmat.h
253: #include finclude/petscpc.h
254: #include finclude/petscts.h
255: TS ts
256: double precision t
257: Vec X,F
258: double precision user(3)
259: integer ierr,i,j,row,mx,my
260: PetscOffset xidx,fidx
261: double precision two,one,lambda
262: double precision hx,hy,hxdhy,hydhx
263: Scalar ut,ub,ul,ur,u,uxx,uyy,sc
264: Scalar xx(1),ff(1)
265: integer param,lmx,lmy
266: parameter (param = 1,lmx = 2,lmy = 3)
268: two = 2.0
269: one = 1.0
271: mx = user(lmx)
272: my = user(lmy)
273: lambda = user(param)
275: hx = 1.0 / dble(mx-1)
276: hy = 1.0 / dble(my-1)
277: sc = hx*hy
278: hxdhy = hx/hy
279: hydhx = hy/hx
281: call VecGetArray(X,xx,xidx,ierr)
282: call VecGetArray(F,ff,fidx,ierr)
283: do 10 j=1,my
284: do 20 i=1,mx
285: row = i + (j-1)*mx
286: if (i .eq. 1 .or. j .eq. 1 .or. &
287: & i .eq. mx .or. j .eq. my) then
288: ff(row+fidx) = xx(row+xidx)
289: else
290: u = xx(row + xidx)
291: ub = xx(row - mx + xidx)
292: ul = xx(row - 1 + xidx)
293: ut = xx(row + mx + xidx)
294: ur = xx(row + 1 + xidx)
295: uxx = (-ur + two*u - ul)*hydhx
296: uyy = (-ut + two*u - ub)*hxdhy
297: ff(row+fidx) = -uxx - uyy + sc*lambda*exp(u)
298: u = -uxx - uyy + sc*lambda*exp(u)
299: endif
300: 20 continue
301: 10 continue
303: call VecRestoreArray(X,xx,xidx,ierr)
304: call VecRestoreArray(F,ff,fidx,ierr)
305: return
306: end
307: !
308: ! -------------------- Evaluate Jacobian of F(x) --------------------
309: !
310: subroutine FormJacobian(ts,ctime,X,JJ,B,flag,user,ierr)
311: implicit none
312: #include finclude/petsc.h
313: #include finclude/petscvec.h
314: #include finclude/petscmat.h
315: #include finclude/petscpc.h
316: #include finclude/petscts.h
317: TS ts
318: Vec X
319: Mat JJ,B
320: MatStructure flag
321: double precision user(3),ctime
322: Mat jac
323: integer i,j,row,mx,my,col(5),ierr
324: PetscOffset xidx
325: Scalar two,one,lambda,v(5),sc,xx(1)
326: double precision hx,hy,hxdhy,hydhx
328: integer param,lmx,lmy
329: parameter (param = 1,lmx = 2,lmy = 3)
331: jac = B
332: two = 2.0
333: one = 1.0
335: mx = user(lmx)
336: my = user(lmy)
337: lambda = user(param)
339: hx = 1.0 / dble(mx-1)
340: hy = 1.0 / dble(my-1)
341: sc = hx*hy
342: hxdhy = hx/hy
343: hydhx = hy/hx
345: call VecGetArray(X,xx,xidx,ierr)
346: do 10 j=1,my
347: do 20 i=1,mx
348: !
349: ! When inserting into PETSc matrices, indices start at 0
350: !
351: ! call PetscTrValid(ierr)
352: row = i - 1 + (j-1)*mx
353: if (i .eq. 1 .or. j .eq. 1 .or. &
354: & i .eq. mx .or. j .eq. my) then
355: call MatSetValues(jac,1,row,1,row,one,INSERT_VALUES,ierr)
356: else
357: v(1) = hxdhy
358: col(1) = row - mx
359: v(2) = hydhx
360: col(2) = row - 1
361: v(3) = -two*(hydhx+hxdhy)+sc*lambda*exp(xx(row+1+xidx))
362: col(3) = row
363: v(4) = hydhx
364: col(4) = row + 1
365: v(5) = hxdhy
366: col(5) = row + mx
367: call MatSetValues(jac,1,row,5,col,v,INSERT_VALUES,ierr)
368: endif
369: 20 continue
370: 10 continue
371: call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
372: call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
373: call VecRestoreArray(X,xx,xidx,ierr)
374: flag = SAME_NONZERO_PATTERN
375: ! call PetscTrValid(ierr)
376: return
377: end