xref: /petsc/src/ts/tutorials/ex1f.F90 (revision d8606c274c09e255c003062beb17b1be973467bc)
1*d8606c27SBarry Smith!
2*d8606c27SBarry Smith!   Solves the time dependent Bratu problem using pseudo-timestepping
3*d8606c27SBarry Smith!
4*d8606c27SBarry Smith!   This code demonstrates how one may solve a nonlinear problem
5*d8606c27SBarry Smith!   with pseudo-timestepping. In this simple example, the pseudo-timestep
6*d8606c27SBarry Smith!   is the same for all grid points, i.e., this is equivalent to using
7*d8606c27SBarry Smith!   the backward Euler method with a variable timestep.
8*d8606c27SBarry Smith!
9*d8606c27SBarry Smith!   Note: This example does not require pseudo-timestepping since it
10*d8606c27SBarry Smith!   is an easy nonlinear problem, but it is included to demonstrate how
11*d8606c27SBarry Smith!   the pseudo-timestepping may be done.
12*d8606c27SBarry Smith!
13*d8606c27SBarry Smith!   See snes/tutorials/ex4.c[ex4f.F] and
14*d8606c27SBarry Smith!   snes/tutorials/ex5.c[ex5f.F] where the problem is described
15*d8606c27SBarry Smith!   and solved using the method of Newton alone.
16*d8606c27SBarry Smith!
17*d8606c27SBarry Smith!
18*d8606c27SBarry Smith!23456789012345678901234567890123456789012345678901234567890123456789012
19*d8606c27SBarry Smith      program main
20*d8606c27SBarry Smith#include <petsc/finclude/petscts.h>
21*d8606c27SBarry Smith      use petscts
22*d8606c27SBarry Smith      implicit none
23*d8606c27SBarry Smith
24*d8606c27SBarry Smith!
25*d8606c27SBarry Smith!  Create an application context to contain data needed by the
26*d8606c27SBarry Smith!  application-provided call-back routines, FormJacobian() and
27*d8606c27SBarry Smith!  FormFunction(). We use a double precision array with three
28*d8606c27SBarry Smith!  entries indexed by param, lmx, lmy.
29*d8606c27SBarry Smith!
30*d8606c27SBarry Smith      PetscReal user(3)
31*d8606c27SBarry Smith      PetscInt          param,lmx,lmy,i5
32*d8606c27SBarry Smith      parameter (param = 1,lmx = 2,lmy = 3)
33*d8606c27SBarry Smith!
34*d8606c27SBarry Smith!   User-defined routines
35*d8606c27SBarry Smith!
36*d8606c27SBarry Smith      external FormJacobian,FormFunction
37*d8606c27SBarry Smith!
38*d8606c27SBarry Smith!   Data for problem
39*d8606c27SBarry Smith!
40*d8606c27SBarry Smith      TS                ts
41*d8606c27SBarry Smith      Vec               x,r
42*d8606c27SBarry Smith      Mat               J
43*d8606c27SBarry Smith      PetscInt           its,N,i1000,itmp
44*d8606c27SBarry Smith      PetscBool  flg
45*d8606c27SBarry Smith      PetscErrorCode      ierr
46*d8606c27SBarry Smith      PetscReal  param_max,param_min,dt
47*d8606c27SBarry Smith      PetscReal  tmax
48*d8606c27SBarry Smith      PetscReal  ftime
49*d8606c27SBarry Smith      TSConvergedReason reason
50*d8606c27SBarry Smith
51*d8606c27SBarry Smith      i5 = 5
52*d8606c27SBarry Smith      param_max = 6.81
53*d8606c27SBarry Smith      param_min = 0
54*d8606c27SBarry Smith
55*d8606c27SBarry Smith      PetscCallA(PetscInitialize(ierr))
56*d8606c27SBarry Smith      user(lmx)        = 4
57*d8606c27SBarry Smith      user(lmy)        = 4
58*d8606c27SBarry Smith      user(param)      = 6.0
59*d8606c27SBarry Smith
60*d8606c27SBarry Smith!
61*d8606c27SBarry Smith!     Allow user to set the grid dimensions and nonlinearity parameter at run-time
62*d8606c27SBarry Smith!
63*d8606c27SBarry Smith      PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-mx',user(lmx),flg,ierr))
64*d8606c27SBarry Smith      itmp = 4
65*d8606c27SBarry Smith      PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-my',itmp,flg,ierr))
66*d8606c27SBarry Smith      user(lmy) = itmp
67*d8606c27SBarry Smith      PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-param',user(param),flg,ierr))
68*d8606c27SBarry Smith      if (user(param) .ge. param_max .or. user(param) .le. param_min) then
69*d8606c27SBarry Smith        print*,'Parameter is out of range'
70*d8606c27SBarry Smith      endif
71*d8606c27SBarry Smith      if (user(lmx) .gt. user(lmy)) then
72*d8606c27SBarry Smith        dt = .5/user(lmx)
73*d8606c27SBarry Smith      else
74*d8606c27SBarry Smith        dt = .5/user(lmy)
75*d8606c27SBarry Smith      endif
76*d8606c27SBarry Smith      PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-dt',dt,flg,ierr))
77*d8606c27SBarry Smith      N          = int(user(lmx)*user(lmy))
78*d8606c27SBarry Smith
79*d8606c27SBarry Smith!
80*d8606c27SBarry Smith!      Create vectors to hold the solution and function value
81*d8606c27SBarry Smith!
82*d8606c27SBarry Smith      PetscCallA(VecCreateSeq(PETSC_COMM_SELF,N,x,ierr))
83*d8606c27SBarry Smith      PetscCallA(VecDuplicate(x,r,ierr))
84*d8606c27SBarry Smith
85*d8606c27SBarry Smith!
86*d8606c27SBarry Smith!    Create matrix to hold Jacobian. Preallocate 5 nonzeros per row
87*d8606c27SBarry Smith!    in the sparse matrix. Note that this is not the optimal strategy see
88*d8606c27SBarry Smith!    the Performance chapter of the users manual for information on
89*d8606c27SBarry Smith!    preallocating memory in sparse matrices.
90*d8606c27SBarry Smith!
91*d8606c27SBarry Smith      i5 = 5
92*d8606c27SBarry Smith      PetscCallA(MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,i5,PETSC_NULL_INTEGER,J,ierr))
93*d8606c27SBarry Smith
94*d8606c27SBarry Smith!
95*d8606c27SBarry Smith!     Create timestepper context
96*d8606c27SBarry Smith!
97*d8606c27SBarry Smith
98*d8606c27SBarry Smith      PetscCallA(TSCreate(PETSC_COMM_WORLD,ts,ierr))
99*d8606c27SBarry Smith      PetscCallA(TSSetProblemType(ts,TS_NONLINEAR,ierr))
100*d8606c27SBarry Smith
101*d8606c27SBarry Smith!
102*d8606c27SBarry Smith!     Tell the timestepper context where to compute solutions
103*d8606c27SBarry Smith!
104*d8606c27SBarry Smith
105*d8606c27SBarry Smith      PetscCallA(TSSetSolution(ts,x,ierr))
106*d8606c27SBarry Smith
107*d8606c27SBarry Smith!
108*d8606c27SBarry Smith!    Provide the call-back for the nonlinear function we are
109*d8606c27SBarry Smith!    evaluating. Thus whenever the timestepping routines need the
110*d8606c27SBarry Smith!    function they will call this routine. Note the final argument
111*d8606c27SBarry Smith!    is the application context used by the call-back functions.
112*d8606c27SBarry Smith!
113*d8606c27SBarry Smith
114*d8606c27SBarry Smith      PetscCallA(TSSetRHSFunction(ts,PETSC_NULL_VEC,FormFunction,user,ierr))
115*d8606c27SBarry Smith
116*d8606c27SBarry Smith!
117*d8606c27SBarry Smith!     Set the Jacobian matrix and the function used to compute
118*d8606c27SBarry Smith!     Jacobians.
119*d8606c27SBarry Smith!
120*d8606c27SBarry Smith
121*d8606c27SBarry Smith      PetscCallA(TSSetRHSJacobian(ts,J,J,FormJacobian,user,ierr))
122*d8606c27SBarry Smith
123*d8606c27SBarry Smith!
124*d8606c27SBarry Smith!       For the initial guess for the problem
125*d8606c27SBarry Smith!
126*d8606c27SBarry Smith
127*d8606c27SBarry Smith      PetscCallA(FormInitialGuess(x,user,ierr))
128*d8606c27SBarry Smith
129*d8606c27SBarry Smith!
130*d8606c27SBarry Smith!       This indicates that we are using pseudo timestepping to
131*d8606c27SBarry Smith!     find a steady state solution to the nonlinear problem.
132*d8606c27SBarry Smith!
133*d8606c27SBarry Smith
134*d8606c27SBarry Smith      PetscCallA(TSSetType(ts,TSPSEUDO,ierr))
135*d8606c27SBarry Smith
136*d8606c27SBarry Smith!
137*d8606c27SBarry Smith!       Set the initial time to start at (this is arbitrary for
138*d8606c27SBarry Smith!     steady state problems and the initial timestep given above
139*d8606c27SBarry Smith!
140*d8606c27SBarry Smith
141*d8606c27SBarry Smith      PetscCallA(TSSetTimeStep(ts,dt,ierr))
142*d8606c27SBarry Smith
143*d8606c27SBarry Smith!
144*d8606c27SBarry Smith!      Set a large number of timesteps and final duration time
145*d8606c27SBarry Smith!     to insure convergence to steady state.
146*d8606c27SBarry Smith!
147*d8606c27SBarry Smith      i1000 = 1000
148*d8606c27SBarry Smith      tmax  = 1.e12
149*d8606c27SBarry Smith      PetscCallA(TSSetMaxSteps(ts,i1000,ierr))
150*d8606c27SBarry Smith      PetscCallA(TSSetMaxTime(ts,tmax,ierr))
151*d8606c27SBarry Smith      PetscCallA(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER,ierr))
152*d8606c27SBarry Smith
153*d8606c27SBarry Smith!
154*d8606c27SBarry Smith!      Set any additional options from the options database. This
155*d8606c27SBarry Smith!     includes all options for the nonlinear and linear solvers used
156*d8606c27SBarry Smith!     internally the timestepping routines.
157*d8606c27SBarry Smith!
158*d8606c27SBarry Smith
159*d8606c27SBarry Smith      PetscCallA(TSSetFromOptions(ts,ierr))
160*d8606c27SBarry Smith
161*d8606c27SBarry Smith      PetscCallA(TSSetUp(ts,ierr))
162*d8606c27SBarry Smith
163*d8606c27SBarry Smith!
164*d8606c27SBarry Smith!      Perform the solve. This is where the timestepping takes place.
165*d8606c27SBarry Smith!
166*d8606c27SBarry Smith      PetscCallA(TSSolve(ts,x,ierr))
167*d8606c27SBarry Smith      PetscCallA(TSGetSolveTime(ts,ftime,ierr))
168*d8606c27SBarry Smith      PetscCallA(TSGetStepNumber(ts,its,ierr))
169*d8606c27SBarry Smith      PetscCallA(TSGetConvergedReason(ts,reason,ierr))
170*d8606c27SBarry Smith
171*d8606c27SBarry Smith      write(6,100) its,ftime,reason
172*d8606c27SBarry Smith 100  format('Number of pseudo time-steps ',i5,' final time ',1pe9.2,' reason ',i3)
173*d8606c27SBarry Smith
174*d8606c27SBarry Smith!
175*d8606c27SBarry Smith!     Free the data structures constructed above
176*d8606c27SBarry Smith!
177*d8606c27SBarry Smith
178*d8606c27SBarry Smith      PetscCallA(VecDestroy(x,ierr))
179*d8606c27SBarry Smith      PetscCallA(VecDestroy(r,ierr))
180*d8606c27SBarry Smith      PetscCallA(MatDestroy(J,ierr))
181*d8606c27SBarry Smith      PetscCallA(TSDestroy(ts,ierr))
182*d8606c27SBarry Smith      PetscCallA(PetscFinalize(ierr))
183*d8606c27SBarry Smith      end
184*d8606c27SBarry Smith
185*d8606c27SBarry Smith!
186*d8606c27SBarry Smith!  --------------------  Form initial approximation -----------------
187*d8606c27SBarry Smith!
188*d8606c27SBarry Smith      subroutine FormInitialGuess(X,user,ierr)
189*d8606c27SBarry Smith      use petscts
190*d8606c27SBarry Smith      implicit none
191*d8606c27SBarry Smith
192*d8606c27SBarry Smith      Vec              X
193*d8606c27SBarry Smith      PetscReal user(3)
194*d8606c27SBarry Smith      PetscInt  i,j,row,mx,my
195*d8606c27SBarry Smith      PetscErrorCode ierr
196*d8606c27SBarry Smith      PetscOffset      xidx
197*d8606c27SBarry Smith      PetscReal one,lambda
198*d8606c27SBarry Smith      PetscReal temp1,temp,hx,hy
199*d8606c27SBarry Smith      PetscScalar      xx(2)
200*d8606c27SBarry Smith      PetscInt          param,lmx,lmy
201*d8606c27SBarry Smith      parameter (param = 1,lmx = 2,lmy = 3)
202*d8606c27SBarry Smith
203*d8606c27SBarry Smith      one = 1.0
204*d8606c27SBarry Smith
205*d8606c27SBarry Smith      mx     = int(user(lmx))
206*d8606c27SBarry Smith      my     = int(user(lmy))
207*d8606c27SBarry Smith      lambda = user(param)
208*d8606c27SBarry Smith
209*d8606c27SBarry Smith      hy    = one / (my-1)
210*d8606c27SBarry Smith      hx    = one / (mx-1)
211*d8606c27SBarry Smith
212*d8606c27SBarry Smith      PetscCall(VecGetArray(X,xx,xidx,ierr))
213*d8606c27SBarry Smith      temp1 = lambda/(lambda + one)
214*d8606c27SBarry Smith      do 10, j=1,my
215*d8606c27SBarry Smith        temp = min(j-1,my-j)*hy
216*d8606c27SBarry Smith        do 20 i=1,mx
217*d8606c27SBarry Smith          row = i + (j-1)*mx
218*d8606c27SBarry Smith          if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then
219*d8606c27SBarry Smith            xx(row+xidx) = 0.0
220*d8606c27SBarry Smith          else
221*d8606c27SBarry Smith            xx(row+xidx) = temp1*sqrt(min(min(i-1,mx-i)*hx,temp))
222*d8606c27SBarry Smith          endif
223*d8606c27SBarry Smith 20     continue
224*d8606c27SBarry Smith 10   continue
225*d8606c27SBarry Smith      PetscCall(VecRestoreArray(X,xx,xidx,ierr))
226*d8606c27SBarry Smith      return
227*d8606c27SBarry Smith      end
228*d8606c27SBarry Smith!
229*d8606c27SBarry Smith!  --------------------  Evaluate Function F(x) ---------------------
230*d8606c27SBarry Smith!
231*d8606c27SBarry Smith      subroutine FormFunction(ts,t,X,F,user,ierr)
232*d8606c27SBarry Smith      use petscts
233*d8606c27SBarry Smith      implicit none
234*d8606c27SBarry Smith
235*d8606c27SBarry Smith      TS       ts
236*d8606c27SBarry Smith      PetscReal  t
237*d8606c27SBarry Smith      Vec               X,F
238*d8606c27SBarry Smith      PetscReal  user(3)
239*d8606c27SBarry Smith      PetscErrorCode     ierr
240*d8606c27SBarry Smith      PetscInt         i,j,row,mx,my
241*d8606c27SBarry Smith      PetscOffset       xidx,fidx
242*d8606c27SBarry Smith      PetscReal  two,lambda
243*d8606c27SBarry Smith      PetscReal  hx,hy,hxdhy,hydhx
244*d8606c27SBarry Smith      PetscScalar  ut,ub,ul,ur,u
245*d8606c27SBarry Smith      PetscScalar  uxx,uyy,sc
246*d8606c27SBarry Smith      PetscScalar  xx(2),ff(2)
247*d8606c27SBarry Smith      PetscInt     param,lmx,lmy
248*d8606c27SBarry Smith      parameter (param = 1,lmx = 2,lmy = 3)
249*d8606c27SBarry Smith
250*d8606c27SBarry Smith      two = 2.0
251*d8606c27SBarry Smith
252*d8606c27SBarry Smith      mx     = int(user(lmx))
253*d8606c27SBarry Smith      my     = int(user(lmy))
254*d8606c27SBarry Smith      lambda = user(param)
255*d8606c27SBarry Smith
256*d8606c27SBarry Smith      hx    = 1.0 / real(mx-1)
257*d8606c27SBarry Smith      hy    = 1.0 / real(my-1)
258*d8606c27SBarry Smith      sc    = hx*hy
259*d8606c27SBarry Smith      hxdhy = hx/hy
260*d8606c27SBarry Smith      hydhx = hy/hx
261*d8606c27SBarry Smith
262*d8606c27SBarry Smith      PetscCall(VecGetArrayRead(X,xx,xidx,ierr))
263*d8606c27SBarry Smith      PetscCall(VecGetArray(F,ff,fidx,ierr))
264*d8606c27SBarry Smith      do 10 j=1,my
265*d8606c27SBarry Smith        do 20 i=1,mx
266*d8606c27SBarry Smith          row = i + (j-1)*mx
267*d8606c27SBarry Smith          if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then
268*d8606c27SBarry Smith            ff(row+fidx) = xx(row+xidx)
269*d8606c27SBarry Smith          else
270*d8606c27SBarry Smith            u            = xx(row + xidx)
271*d8606c27SBarry Smith            ub           = xx(row - mx + xidx)
272*d8606c27SBarry Smith            ul           = xx(row - 1 + xidx)
273*d8606c27SBarry Smith            ut           = xx(row + mx + xidx)
274*d8606c27SBarry Smith            ur           = xx(row + 1 + xidx)
275*d8606c27SBarry Smith            uxx          = (-ur + two*u - ul)*hydhx
276*d8606c27SBarry Smith            uyy          = (-ut + two*u - ub)*hxdhy
277*d8606c27SBarry Smith            ff(row+fidx) = -uxx - uyy + sc*lambda*exp(u)
278*d8606c27SBarry Smith            u =  -uxx - uyy + sc*lambda*exp(u)
279*d8606c27SBarry Smith         endif
280*d8606c27SBarry Smith 20   continue
281*d8606c27SBarry Smith 10   continue
282*d8606c27SBarry Smith
283*d8606c27SBarry Smith      PetscCall(VecRestoreArrayRead(X,xx,xidx,ierr))
284*d8606c27SBarry Smith      PetscCall(VecRestoreArray(F,ff,fidx,ierr))
285*d8606c27SBarry Smith      return
286*d8606c27SBarry Smith      end
287*d8606c27SBarry Smith!
288*d8606c27SBarry Smith!  --------------------  Evaluate Jacobian of F(x) --------------------
289*d8606c27SBarry Smith!
290*d8606c27SBarry Smith      subroutine FormJacobian(ts,ctime,X,JJ,B,user,ierr)
291*d8606c27SBarry Smith      use petscts
292*d8606c27SBarry Smith      implicit none
293*d8606c27SBarry Smith
294*d8606c27SBarry Smith      TS               ts
295*d8606c27SBarry Smith      Vec              X
296*d8606c27SBarry Smith      Mat              JJ,B
297*d8606c27SBarry Smith      PetscReal user(3),ctime
298*d8606c27SBarry Smith      PetscErrorCode   ierr
299*d8606c27SBarry Smith      Mat              jac
300*d8606c27SBarry Smith      PetscOffset xidx
301*d8606c27SBarry Smith      PetscInt    i,j,row(1),mx,my
302*d8606c27SBarry Smith      PetscInt    col(5),i1,i5
303*d8606c27SBarry Smith      PetscScalar two,one,lambda
304*d8606c27SBarry Smith      PetscScalar v(5),sc,xx(2)
305*d8606c27SBarry Smith      PetscReal hx,hy,hxdhy,hydhx
306*d8606c27SBarry Smith
307*d8606c27SBarry Smith      PetscInt  param,lmx,lmy
308*d8606c27SBarry Smith      parameter (param = 1,lmx = 2,lmy = 3)
309*d8606c27SBarry Smith
310*d8606c27SBarry Smith      i1 = 1
311*d8606c27SBarry Smith      i5 = 5
312*d8606c27SBarry Smith      jac = B
313*d8606c27SBarry Smith      two = 2.0
314*d8606c27SBarry Smith      one = 1.0
315*d8606c27SBarry Smith
316*d8606c27SBarry Smith      mx     = int(user(lmx))
317*d8606c27SBarry Smith      my     = int(user(lmy))
318*d8606c27SBarry Smith      lambda = user(param)
319*d8606c27SBarry Smith
320*d8606c27SBarry Smith      hx    = 1.0 / real(mx-1)
321*d8606c27SBarry Smith      hy    = 1.0 / real(my-1)
322*d8606c27SBarry Smith      sc    = hx*hy
323*d8606c27SBarry Smith      hxdhy = hx/hy
324*d8606c27SBarry Smith      hydhx = hy/hx
325*d8606c27SBarry Smith
326*d8606c27SBarry Smith      PetscCall(VecGetArrayRead(X,xx,xidx,ierr))
327*d8606c27SBarry Smith      do 10 j=1,my
328*d8606c27SBarry Smith        do 20 i=1,mx
329*d8606c27SBarry Smith!
330*d8606c27SBarry Smith!      When inserting into PETSc matrices, indices start at 0
331*d8606c27SBarry Smith!
332*d8606c27SBarry Smith          row(1) = i - 1 + (j-1)*mx
333*d8606c27SBarry Smith          if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then
334*d8606c27SBarry Smith            PetscCall(MatSetValues(jac,i1,row,i1,row,one,INSERT_VALUES,ierr))
335*d8606c27SBarry Smith          else
336*d8606c27SBarry Smith            v(1)   = hxdhy
337*d8606c27SBarry Smith            col(1) = row(1) - mx
338*d8606c27SBarry Smith            v(2)   = hydhx
339*d8606c27SBarry Smith            col(2) = row(1) - 1
340*d8606c27SBarry Smith            v(3)   = -two*(hydhx+hxdhy)+sc*lambda*exp(xx(row(1)+1+xidx))
341*d8606c27SBarry Smith            col(3) = row(1)
342*d8606c27SBarry Smith            v(4)   = hydhx
343*d8606c27SBarry Smith            col(4) = row(1) + 1
344*d8606c27SBarry Smith            v(5)   = hxdhy
345*d8606c27SBarry Smith            col(5) = row(1) + mx
346*d8606c27SBarry Smith            PetscCall(MatSetValues(jac,i1,row,i5,col,v,INSERT_VALUES,ierr))
347*d8606c27SBarry Smith          endif
348*d8606c27SBarry Smith 20     continue
349*d8606c27SBarry Smith 10   continue
350*d8606c27SBarry Smith      PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr))
351*d8606c27SBarry Smith      PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr))
352*d8606c27SBarry Smith      PetscCall(VecRestoreArray(X,xx,xidx,ierr))
353*d8606c27SBarry Smith      return
354*d8606c27SBarry Smith      end
355*d8606c27SBarry Smith
356*d8606c27SBarry Smith!/*TEST
357*d8606c27SBarry Smith!
358*d8606c27SBarry Smith!    test:
359*d8606c27SBarry Smith!      TODO: broken
360*d8606c27SBarry Smith!      args: -ksp_gmres_cgs_refinement_type refine_always -snes_type newtonls -ts_monitor_pseudo -ts_max_snes_failures 3 -ts_pseudo_frtol 1.e-5  -snes_stol 1e-5
361*d8606c27SBarry Smith!
362*d8606c27SBarry Smith!TEST*/
363*d8606c27SBarry Smith
364