Actual source code: ex5f.F
1: ! "$Id: ex5f.F,v 1.80 2001/08/24 16:23:36 bsmith Exp $";
2: !
3: ! Description: This example solves a nonlinear system in parallel with SNES.
4: ! We solve the Bratu (SFI - solid fuel ignition) problem in a 2D rectangular
5: ! domain, using distributed arrays (DAs) to partition the parallel grid.
6: ! The command line options include:
7: ! -par <param>, where <param> indicates the nonlinearity of the problem
8: ! problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81)
9: !
10: ! Program usage: mpirun -np <procs> ex5f [-help] [all PETSc options]
11: !
12: !/*T
13: ! Concepts: SNES^parallel Bratu example
14: ! Concepts: DA^using distributed arrays;
15: ! Processors: n
16: !T*/
17: !
18: ! --------------------------------------------------------------------------
19: !
20: ! Solid Fuel Ignition (SFI) problem. This problem is modeled by
21: ! the partial differential equation
22: !
23: ! -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1,
24: !
25: ! with boundary conditions
26: !
27: ! u = 0 for x = 0, x = 1, y = 0, y = 1.
28: !
29: ! A finite difference approximation with the usual 5-point stencil
30: ! is used to discretize the boundary value problem to obtain a nonlinear
31: ! system of equations.
32: !
33: ! --------------------------------------------------------------------------
35: program main
36: implicit none
37: !
38: ! We place common blocks, variable declarations, and other include files
39: ! needed for this code in the single file ex5f.h. We then need to include
40: ! only this file throughout the various routines in this program. See
41: ! additional comments in the file ex5f.h.
42: !
43: #include "ex5f.h"
45: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
46: ! Variable declarations
47: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48: !
49: ! Variables:
50: ! snes - nonlinear solver
51: ! x, r - solution, residual vectors
52: ! J - Jacobian matrix
53: ! its - iterations for convergence
54: !
55: ! See additional variable declarations in the file ex5f.h
56: !
57: SNES snes
58: Vec x,r
59: Mat J,A
60: integer its,flg,ierr
61: double precision lambda_max,lambda_min
62: ISColoring coloring
63: PetscTruth adifor_jacobian,adiformf_jacobian
65: ! Note: Any user-defined Fortran routines (such as FormJacobianLocal)
66: ! MUST be declared as external.
68: external FormInitialGuess
69: external FormFunctionLocal,FormJacobianLocal
70: #if defined(PETSC_HAVE_ADIFOR)
71: external g_FormFunctionLocal,m_FormFunctionLocal
72: #endif
74: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75: ! Initialize program
76: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
78: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
79: call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
80: call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
82: ! Initialize problem parameters
84: lambda_max = 6.81
85: lambda_min = 0.0
86: lambda = 6.0
87: call PetscOptionsGetReal(PETSC_NULL_CHARACTER,'-par',lambda, &
88: & flg,ierr)
89: if (lambda .ge. lambda_max .or. lambda .le. lambda_min) then
90: if (rank .eq. 0) write(6,*) 'Lambda is out of range'
91: SETERRQ(1,' ',ierr)
92: endif
94: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95: ! Create nonlinear solver context
96: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98: call SNESCreate(PETSC_COMM_WORLD,snes,ierr)
100: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101: ! Create vector data structures; set function evaluation routine
102: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104: ! Create distributed array (DA) to manage parallel grid and vectors
106: ! This really needs only the star-type stencil, but we use the box
107: ! stencil temporarily.
108: call DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_BOX,-4, &
109: & -4,PETSC_DECIDE,PETSC_DECIDE,1,1,PETSC_NULL_INTEGER, &
110: & PETSC_NULL_INTEGER,da,ierr)
112: ! Extract global and local vectors from DA; then duplicate for remaining
113: ! vectors that are the same types
115: call DACreateGlobalVector(da,x,ierr)
116: call VecDuplicate(x,r,ierr)
118: ! Get local grid boundaries (for 2-dimensional DA)
120: call DAGetInfo(da,PETSC_NULL_INTEGER,mx,my,PETSC_NULL_INTEGER, &
121: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
122: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
123: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
124: & PETSC_NULL_INTEGER,ierr)
125: call DAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym, &
126: & PETSC_NULL_INTEGER,ierr)
127: call DAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,gxm,gym, &
128: & PETSC_NULL_INTEGER,ierr)
130: ! Here we shift the starting indices up by one so that we can easily
131: ! use the Fortran convention of 1-based indices (rather 0-based indices).
133: xs = xs+1
134: ys = ys+1
135: gxs = gxs+1
136: gys = gys+1
138: ye = ys+ym-1
139: xe = xs+xm-1
140: gye = gys+gym-1
141: gxe = gxs+gxm-1
143: ! Set function evaluation routine and vector
145: call DASetLocalFunction(da,FormFunctionLocal,ierr)
146: call DASetLocalJacobian(da,FormJacobianLocal,ierr)
147: #if defined(PETSC_HAVE_ADIFOR)
148: call DASetLocalAdiforFunction(da, &
149: & g_FormFunctionLocal,ierr)
150: #endif
151: call SNESSetFunction(snes,r,SNESDAFormFunction,da,ierr)
153: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
154: ! Create matrix data structure; set Jacobian evaluation routine
155: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157: ! Set Jacobian matrix data structure and default Jacobian evaluation
158: ! routine. User can override with:
159: ! -snes_mf : matrix-free Newton-Krylov method with no preconditioning
160: ! (unless user explicitly sets preconditioner)
161: ! -snes_mf_operator : form preconditioning matrix as set by the user,
162: ! but use matrix-free approx for Jacobian-vector
163: ! products within Newton-Krylov method
164: !
166: call DAGetMatrix(da,MATMPIAIJ,J,ierr)
168: #if defined(PETSC_HAVE_ADIFOR)
169: call PetscOptionsGetLogical(PETSC_NULL_CHARACTER &
170: & ,'-adiformf_jacobian', &
171: & adiformf_jacobian,PETSC_NULL_INTEGER,ierr)
172: if (adiformf_jacobian .eq. 1) then
173: call DASetLocalAdiforMFFunction(da, &
174: & m_FormFunctionLocal,ierr)
175: call MatRegisterDAAD(ierr)
176: call MatCreateDAAD(da,A,ierr)
177: call MatDAADSetSNES(A,snes,ierr)
178: else
179: A = J
180: endif
181: #else
182: A = J
183: #endif
185: call SNESSetJacobian(snes,A,J,SNESDAComputeJacobian, &
186: & da,ierr)
188: #if defined(PETSC_HAVE_ADIFOR)
189: call PetscOptionsGetLogical(PETSC_NULL_CHARACTER &
190: & ,'-adifor_jacobian', &
191: & adifor_jacobian,PETSC_NULL_INTEGER,ierr)
192: if (adifor_jacobian .eq. 1) then
193: call SNESSetJacobian(snes,A,J,SNESDAComputeJacobianWithAdifor, &
194: & da,ierr)
195: endif
196: #endif
199: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
200: ! Customize nonlinear solver; set runtime options
201: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
203: ! Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
205: call SNESSetFromOptions(snes,ierr)
207: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
208: ! Evaluate initial guess; then solve nonlinear system.
209: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
211: ! Note: The user should initialize the vector, x, with the initial guess
212: ! for the nonlinear solver prior to calling SNESSolve(). In particular,
213: ! to employ an initial guess of zero, the user should explicitly set
214: ! this vector to zero by calling VecSet().
216: call FormInitialGuess(x,ierr)
217: call SNESSolve(snes,x,its,ierr)
218: if (rank .eq. 0) then
219: write(6,100) its
220: endif
221: 100 format('Number of Newton iterations = ',i5)
224: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
225: ! Free work space. All PETSc objects should be destroyed when they
226: ! are no longer needed.
227: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
229: if (A .ne. J) call MatDestroy(A,ierr)
230: call MatDestroy(J,ierr)
231: call VecDestroy(x,ierr)
232: call VecDestroy(r,ierr)
233: call SNESDestroy(snes,ierr)
234: call DADestroy(da,ierr)
235: call PetscFinalize(ierr)
236: end
238: ! ---------------------------------------------------------------------
239: !
240: ! FormInitialGuess - Forms initial approximation.
241: !
242: ! Input Parameters:
243: ! X - vector
244: !
245: ! Output Parameter:
246: ! X - vector
247: !
248: ! Notes:
249: ! This routine serves as a wrapper for the lower-level routine
250: ! "ApplicationInitialGuess", where the actual computations are
251: ! done using the standard Fortran style of treating the local
252: ! vector data as a multidimensional array over the local mesh.
253: ! This routine merely handles ghost point scatters and accesses
254: ! the local vector data via VecGetArray() and VecRestoreArray().
255: !
256: subroutine FormInitialGuess(X,ierr)
257: implicit none
259: #include "ex5f.h"
261: ! Input/output variables:
262: Vec X
263: integer ierr
265: ! Declarations for use with local arrays:
266: PetscScalar lx_v(0:1)
267: PetscOffset lx_i
268: Vec localX
270: 0
272: ! Get a pointer to vector data.
273: ! - For default PETSc vectors, VecGetArray() returns a pointer to
274: ! the data array. Otherwise, the routine is implementation dependent.
275: ! - You MUST call VecRestoreArray() when you no longer need access to
276: ! the array.
277: ! - Note that the Fortran interface to VecGetArray() differs from the
278: ! C version. See the users manual for details.
280: call DAGetLocalVector(da,localX,ierr)
281: call VecGetArray(localX,lx_v,lx_i,ierr)
283: ! Compute initial guess over the locally owned part of the grid
285: call InitialGuessLocal(lx_v(lx_i),ierr)
287: ! Restore vector
289: call VecRestoreArray(localX,lx_v,lx_i,ierr)
291: ! Insert values into global vector
293: call DALocalToGlobal(da,localX,INSERT_VALUES,X,ierr)
294: call DARestoreLocalVector(da,localX,ierr)
295: return
296: end
298: ! ---------------------------------------------------------------------
299: !
300: ! InitialGuessLocal - Computes initial approximation, called by
301: ! the higher level routine FormInitialGuess().
302: !
303: ! Input Parameter:
304: ! x - local vector data
305: !
306: ! Output Parameters:
307: ! x - local vector data
308: ! ierr - error code
309: !
310: ! Notes:
311: ! This routine uses standard Fortran-style computations over a 2-dim array.
312: !
313: subroutine InitialGuessLocal(x,ierr)
314: implicit none
316: #include "ex5f.h"
318: ! Input/output variables:
319: PetscScalar x(gxs:gxe,gys:gye)
320: integer ierr
322: ! Local variables:
323: integer i,j,hxdhy,hydhx
324: PetscScalar temp1,temp,hx,hy,sc,one
326: ! Set parameters
328: ierr = 0
329: one = 1.0
330: hx = one/(dble(mx-1))
331: hy = one/(dble(my-1))
332: sc = hx*hy*lambda
333: hxdhy = hx/hy
334: hydhx = hy/hx
335: temp1 = lambda/(lambda + one)
337: do 20 j=ys,ye
338: temp = dble(min(j-1,my-j))*hy
339: do 10 i=xs,xe
340: if (i .eq. 1 .or. j .eq. 1 &
341: & .or. i .eq. mx .or. j .eq. my) then
342: x(i,j) = 0.0
343: else
344: x(i,j) = temp1 * &
345: & sqrt(min(dble(min(i-1,mx-i)*hx),dble(temp)))
346: endif
347: 10 continue
348: 20 continue
350: return
351: end
353: ! ---------------------------------------------------------------------
354: !
355: ! FormFunctionLocal - Computes nonlinear function, called by
356: ! the higher level routine FormFunction().
357: !
358: ! Input Parameter:
359: ! x - local vector data
360: !
361: ! Output Parameters:
362: ! f - local vector data, f(x)
363: ! ierr - error code
364: !
365: ! Notes:
366: ! This routine uses standard Fortran-style computations over a 2-dim array.
367: !
368: ! Process adifor: FormFunctionLocal
369: !
370: subroutine FormFunctionLocal(info,x,f,dummy,ierr)
372: implicit none
374: #include "ex5f.h"
376: ! Input/output variables:
377: DALocalInfo info(DA_LOCAL_INFO_SIZE)
378: PetscScalar x(gxs:gxe,gys:gye)
379: PetscScalar f(xs:xe,ys:ye)
380: integer ierr
381: PetscObject dummy
383: ! Local variables:
384: PetscScalar two,one,hx,hy,hxdhy,hydhx,sc
385: PetscScalar u,uxx,uyy
386: integer i,j
389: xs = info(DA_LOCAL_INFO_XS)+1
390: xe = xs+info(DA_LOCAL_INFO_XM)-1
391: ys = info(DA_LOCAL_INFO_YS)+1
392: ye = ys+info(DA_LOCAL_INFO_YM)-1
393: mx = info(DA_LOCAL_INFO_MX)
394: my = info(DA_LOCAL_INFO_MY)
396: one = 1.0
397: two = 2.0
398: hx = one/dble(mx-1)
399: hy = one/dble(my-1)
400: sc = hx*hy*lambda
401: hxdhy = hx/hy
402: hydhx = hy/hx
404: ! Compute function over the locally owned part of the grid
406: do 20 j=ys,ye
407: do 10 i=xs,xe
408: if (i .eq. 1 .or. j .eq. 1 &
409: & .or. i .eq. mx .or. j .eq. my) then
410: f(i,j) = x(i,j)
411: else
412: u = x(i,j)
413: uxx = hydhx * (two*u &
414: & - x(i-1,j) - x(i+1,j))
415: uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1))
416: f(i,j) = uxx + uyy - sc*exp(u)
417: endif
418: 10 continue
419: 20 continue
421: call PetscLogFlops(11*ym*xm,ierr)
423: return
424: end
426: ! ---------------------------------------------------------------------
427: !
428: ! FormJacobianLocal - Computes Jacobian matrix, called by
429: ! the higher level routine FormJacobian().
430: !
431: ! Input Parameters:
432: ! x - local vector data
433: !
434: ! Output Parameters:
435: ! jac - Jacobian matrix
436: ! jac_prec - optionally different preconditioning matrix (not used here)
437: ! ierr - error code
438: !
439: ! Notes:
440: ! This routine uses standard Fortran-style computations over a 2-dim array.
441: !
442: ! Notes:
443: ! Due to grid point reordering with DAs, we must always work
444: ! with the local grid points, and then transform them to the new
445: ! global numbering with the "ltog" mapping (via DAGetGlobalIndices()).
446: ! We cannot work directly with the global numbers for the original
447: ! uniprocessor grid!
448: !
449: ! Two methods are available for imposing this transformation
450: ! when setting matrix entries:
451: ! (A) MatSetValuesLocal(), using the local ordering (including
452: ! ghost points!)
453: ! - Use DAGetGlobalIndices() to extract the local-to-global map
454: ! - Associate this map with the matrix by calling
455: ! MatSetLocalToGlobalMapping() once
456: ! - Set matrix entries using the local ordering
457: ! by calling MatSetValuesLocal()
458: ! (B) MatSetValues(), using the global ordering
459: ! - Use DAGetGlobalIndices() to extract the local-to-global map
460: ! - Then apply this map explicitly yourself
461: ! - Set matrix entries using the global ordering by calling
462: ! MatSetValues()
463: ! Option (A) seems cleaner/easier in many cases, and is the procedure
464: ! used in this example.
465: !
466: subroutine FormJacobianLocal(info,x,jac,ctx,ierr)
467: implicit none
469: #include "ex5f.h"
471: ! Input/output variables:
472: PetscScalar x(gxs:gxe,gys:gye)
473: Mat jac
474: integer ierr,ctx
475: DALocalInfo info(DA_LOCAL_INFO_SIZE)
476:
478: ! Local variables:
479: integer row,col(5),i,j
480: PetscScalar two,one,hx,hy,hxdhy,hydhx,sc,v(5)
482: ! Set parameters
484: one = 1.0
485: two = 2.0
486: hx = one/dble(mx-1)
487: hy = one/dble(my-1)
488: sc = hx*hy
489: hxdhy = hx/hy
490: hydhx = hy/hx
492: ! Compute entries for the locally owned part of the Jacobian.
493: ! - Currently, all PETSc parallel matrix formats are partitioned by
494: ! contiguous chunks of rows across the processors.
495: ! - Each processor needs to insert only elements that it owns
496: ! locally (but any non-local elements will be sent to the
497: ! appropriate processor during matrix assembly).
498: ! - Here, we set all entries for a particular row at once.
499: ! - We can set matrix entries either using either
500: ! MatSetValuesLocal() or MatSetValues(), as discussed above.
501: ! - Note that MatSetValues() uses 0-based row and column numbers
502: ! in Fortran as well as in C.
504: do 20 j=ys,ye
505: row = (j - gys)*gxm + xs - gxs - 1
506: do 10 i=xs,xe
507: row = row + 1
508: ! boundary points
509: if (i .eq. 1 .or. j .eq. 1 &
510: & .or. i .eq. mx .or. j .eq. my) then
511: call MatSetValuesLocal(jac,1,row,1,row,one, &
512: & INSERT_VALUES,ierr)
513: ! interior grid points
514: else
515: v(1) = -hxdhy
516: v(2) = -hydhx
517: v(3) = two*(hydhx + hxdhy) &
518: & - sc*lambda*exp(x(i,j))
519: v(4) = -hydhx
520: v(5) = -hxdhy
521: col(1) = row - gxm
522: col(2) = row - 1
523: col(3) = row
524: col(4) = row + 1
525: col(5) = row + gxm
526: call MatSetValuesLocal(jac,1,row,5,col,v, &
527: & INSERT_VALUES,ierr)
528: endif
529: 10 continue
530: 20 continue
531: call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
532: call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
533: return
534: end