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_NONLINEAR_EQUATIONS, &
99: & snes,ierr)
101: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102: ! Create vector data structures; set function evaluation routine
103: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105: ! Create distributed array (DA) to manage parallel grid and vectors
107: ! This really needs only the star-type stencil, but we use the box
108: ! stencil temporarily.
109: call DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_BOX,-4, &
110: & -4,PETSC_DECIDE,PETSC_DECIDE,1,1,PETSC_NULL_INTEGER, &
111: & PETSC_NULL_INTEGER,da,ierr)
113: ! Extract global and local vectors from DA; then duplicate for remaining
114: ! vectors that are the same types
116: call DACreateGlobalVector(da,x,ierr)
117: call VecDuplicate(x,r,ierr)
119: ! Get local grid boundaries (for 2-dimensional DA)
121: call DAGetInfo(da,PETSC_NULL_INTEGER,mx,my,PETSC_NULL_INTEGER, &
122: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
123: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
124: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
125: & PETSC_NULL_INTEGER,ierr)
126: call DAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym, &
127: & PETSC_NULL_INTEGER,ierr)
128: call DAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,gxm,gym, &
129: & PETSC_NULL_INTEGER,ierr)
131: ! Here we shift the starting indices up by one so that we can easily
132: ! use the Fortran convention of 1-based indices (rather 0-based indices).
134: xs = xs+1
135: ys = ys+1
136: gxs = gxs+1
137: gys = gys+1
139: ye = ys+ym-1
140: xe = xs+xm-1
141: gye = gys+gym-1
142: gxe = gxs+gxm-1
144: ! Set function evaluation routine and vector
146: call DASetLocalFunction(da,FormFunctionLocal,ierr)
147: call DASetLocalJacobian(da,FormJacobianLocal,ierr)
148: #if defined(PETSC_HAVE_ADIFOR)
149: call DASetLocalAdiforFunction(da, &
150: & g_FormFunctionLocal,ierr)
151: #endif
152: call SNESSetFunction(snes,r,SNESDAFormFunction,da,ierr)
154: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
155: ! Create matrix data structure; set Jacobian evaluation routine
156: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158: ! Set Jacobian matrix data structure and default Jacobian evaluation
159: ! routine. User can override with:
160: ! -snes_mf : matrix-free Newton-Krylov method with no preconditioning
161: ! (unless user explicitly sets preconditioner)
162: ! -snes_mf_operator : form preconditioning matrix as set by the user,
163: ! but use matrix-free approx for Jacobian-vector
164: ! products within Newton-Krylov method
165: !
167: call DAGetMatrix(da,MATMPIAIJ,J,ierr)
169: #if defined(PETSC_HAVE_ADIFOR)
170: call PetscOptionsGetLogical(PETSC_NULL_CHARACTER &
171: & ,'-adiformf_jacobian', &
172: & adiformf_jacobian,PETSC_NULL_INTEGER,ierr)
173: if (adiformf_jacobian .eq. 1) then
174: call DASetLocalAdiforMFFunction(da, &
175: & m_FormFunctionLocal,ierr)
176: call MatRegisterDAAD(ierr)
177: call MatCreateDAAD(da,A,ierr)
178: call MatDAADSetSNES(A,snes,ierr)
179: else
180: A = J
181: endif
182: #else
183: A = J
184: #endif
186: call SNESSetJacobian(snes,A,J,SNESDAComputeJacobian, &
187: & da,ierr)
189: #if defined(PETSC_HAVE_ADIFOR)
190: call PetscOptionsGetLogical(PETSC_NULL_CHARACTER &
191: & ,'-adifor_jacobian', &
192: & adifor_jacobian,PETSC_NULL_INTEGER,ierr)
193: if (adifor_jacobian .eq. 1) then
194: call SNESSetJacobian(snes,A,J,SNESDAComputeJacobianWithAdifor, &
195: & da,ierr)
196: endif
197: #endif
200: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
201: ! Customize nonlinear solver; set runtime options
202: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204: ! Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
206: call SNESSetFromOptions(snes,ierr)
208: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
209: ! Evaluate initial guess; then solve nonlinear system.
210: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
212: ! Note: The user should initialize the vector, x, with the initial guess
213: ! for the nonlinear solver prior to calling SNESSolve(). In particular,
214: ! to employ an initial guess of zero, the user should explicitly set
215: ! this vector to zero by calling VecSet().
217: call FormInitialGuess(x,ierr)
218: call SNESSolve(snes,x,its,ierr)
219: if (rank .eq. 0) then
220: write(6,100) its
221: endif
222: 100 format('Number of Newton iterations = ',i5)
225: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226: ! Free work space. All PETSc objects should be destroyed when they
227: ! are no longer needed.
228: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
230: if (A .ne. J) call MatDestroy(A,ierr)
231: call MatDestroy(J,ierr)
232: call VecDestroy(x,ierr)
233: call VecDestroy(r,ierr)
234: call SNESDestroy(snes,ierr)
235: call DADestroy(da,ierr)
236: call PetscFinalize(ierr)
237: end
239: ! ---------------------------------------------------------------------
240: !
241: ! FormInitialGuess - Forms initial approximation.
242: !
243: ! Input Parameters:
244: ! X - vector
245: !
246: ! Output Parameter:
247: ! X - vector
248: !
249: ! Notes:
250: ! This routine serves as a wrapper for the lower-level routine
251: ! "ApplicationInitialGuess", where the actual computations are
252: ! done using the standard Fortran style of treating the local
253: ! vector data as a multidimensional array over the local mesh.
254: ! This routine merely handles ghost point scatters and accesses
255: ! the local vector data via VecGetArray() and VecRestoreArray().
256: !
257: subroutine FormInitialGuess(X,ierr)
258: implicit none
260: #include "ex5f.h"
262: ! Input/output variables:
263: Vec X
264: integer ierr
266: ! Declarations for use with local arrays:
267: PetscScalar lx_v(0:1)
268: PetscOffset lx_i
269: Vec localX
271: 0
273: ! Get a pointer to vector data.
274: ! - For default PETSc vectors, VecGetArray() returns a pointer to
275: ! the data array. Otherwise, the routine is implementation dependent.
276: ! - You MUST call VecRestoreArray() when you no longer need access to
277: ! the array.
278: ! - Note that the Fortran interface to VecGetArray() differs from the
279: ! C version. See the users manual for details.
281: call DAGetLocalVector(da,localX,ierr)
282: call VecGetArray(localX,lx_v,lx_i,ierr)
284: ! Compute initial guess over the locally owned part of the grid
286: call InitialGuessLocal(lx_v(lx_i),ierr)
288: ! Restore vector
290: call VecRestoreArray(localX,lx_v,lx_i,ierr)
292: ! Insert values into global vector
294: call DALocalToGlobal(da,localX,INSERT_VALUES,X,ierr)
295: call DARestoreLocalVector(da,localX,ierr)
296: return
297: end
299: ! ---------------------------------------------------------------------
300: !
301: ! InitialGuessLocal - Computes initial approximation, called by
302: ! the higher level routine FormInitialGuess().
303: !
304: ! Input Parameter:
305: ! x - local vector data
306: !
307: ! Output Parameters:
308: ! x - local vector data
309: ! ierr - error code
310: !
311: ! Notes:
312: ! This routine uses standard Fortran-style computations over a 2-dim array.
313: !
314: subroutine InitialGuessLocal(x,ierr)
315: implicit none
317: #include "ex5f.h"
319: ! Input/output variables:
320: PetscScalar x(gxs:gxe,gys:gye)
321: integer ierr
323: ! Local variables:
324: integer i,j,hxdhy,hydhx
325: PetscScalar temp1,temp,hx,hy,sc,one
327: ! Set parameters
329: ierr = 0
330: one = 1.0
331: hx = one/(dble(mx-1))
332: hy = one/(dble(my-1))
333: sc = hx*hy*lambda
334: hxdhy = hx/hy
335: hydhx = hy/hx
336: temp1 = lambda/(lambda + one)
338: do 20 j=ys,ye
339: temp = dble(min(j-1,my-j))*hy
340: do 10 i=xs,xe
341: if (i .eq. 1 .or. j .eq. 1 &
342: & .or. i .eq. mx .or. j .eq. my) then
343: x(i,j) = 0.0
344: else
345: x(i,j) = temp1 * &
346: & sqrt(min(dble(min(i-1,mx-i)*hx),dble(temp)))
347: endif
348: 10 continue
349: 20 continue
351: return
352: end
354: ! ---------------------------------------------------------------------
355: !
356: ! FormFunctionLocal - Computes nonlinear function, called by
357: ! the higher level routine FormFunction().
358: !
359: ! Input Parameter:
360: ! x - local vector data
361: !
362: ! Output Parameters:
363: ! f - local vector data, f(x)
364: ! ierr - error code
365: !
366: ! Notes:
367: ! This routine uses standard Fortran-style computations over a 2-dim array.
368: !
369: ! Process adifor: FormFunctionLocal
370: !
371: subroutine FormFunctionLocal(info,x,f,dummy,ierr)
373: implicit none
375: #include "ex5f.h"
377: ! Input/output variables:
378: DALocalInfo info(DA_LOCAL_INFO_SIZE)
379: PetscScalar x(gxs:gxe,gys:gye)
380: PetscScalar f(xs:xe,ys:ye)
381: integer ierr
382: PetscObject dummy
384: ! Local variables:
385: PetscScalar two,one,hx,hy,hxdhy,hydhx,sc
386: PetscScalar u,uxx,uyy
387: integer i,j
390: xs = info(DA_LOCAL_INFO_XS)+1
391: xe = xs+info(DA_LOCAL_INFO_XM)-1
392: ys = info(DA_LOCAL_INFO_YS)+1
393: ye = ys+info(DA_LOCAL_INFO_YM)-1
394: mx = info(DA_LOCAL_INFO_MX)
395: my = info(DA_LOCAL_INFO_MY)
397: one = 1.0
398: two = 2.0
399: hx = one/dble(mx-1)
400: hy = one/dble(my-1)
401: sc = hx*hy*lambda
402: hxdhy = hx/hy
403: hydhx = hy/hx
405: ! Compute function over the locally owned part of the grid
407: do 20 j=ys,ye
408: do 10 i=xs,xe
409: if (i .eq. 1 .or. j .eq. 1 &
410: & .or. i .eq. mx .or. j .eq. my) then
411: f(i,j) = x(i,j)
412: else
413: u = x(i,j)
414: uxx = hydhx * (two*u &
415: & - x(i-1,j) - x(i+1,j))
416: uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1))
417: f(i,j) = uxx + uyy - sc*exp(u)
418: endif
419: 10 continue
420: 20 continue
422: call PetscLogFlops(11*ym*xm,ierr)
424: return
425: end
427: ! ---------------------------------------------------------------------
428: !
429: ! FormJacobianLocal - Computes Jacobian matrix, called by
430: ! the higher level routine FormJacobian().
431: !
432: ! Input Parameters:
433: ! x - local vector data
434: !
435: ! Output Parameters:
436: ! jac - Jacobian matrix
437: ! jac_prec - optionally different preconditioning matrix (not used here)
438: ! ierr - error code
439: !
440: ! Notes:
441: ! This routine uses standard Fortran-style computations over a 2-dim array.
442: !
443: ! Notes:
444: ! Due to grid point reordering with DAs, we must always work
445: ! with the local grid points, and then transform them to the new
446: ! global numbering with the "ltog" mapping (via DAGetGlobalIndices()).
447: ! We cannot work directly with the global numbers for the original
448: ! uniprocessor grid!
449: !
450: ! Two methods are available for imposing this transformation
451: ! when setting matrix entries:
452: ! (A) MatSetValuesLocal(), using the local ordering (including
453: ! ghost points!)
454: ! - Use DAGetGlobalIndices() to extract the local-to-global map
455: ! - Associate this map with the matrix by calling
456: ! MatSetLocalToGlobalMapping() once
457: ! - Set matrix entries using the local ordering
458: ! by calling MatSetValuesLocal()
459: ! (B) MatSetValues(), using the global ordering
460: ! - Use DAGetGlobalIndices() to extract the local-to-global map
461: ! - Then apply this map explicitly yourself
462: ! - Set matrix entries using the global ordering by calling
463: ! MatSetValues()
464: ! Option (A) seems cleaner/easier in many cases, and is the procedure
465: ! used in this example.
466: !
467: subroutine FormJacobianLocal(info,x,jac,ctx,ierr)
468: implicit none
470: #include "ex5f.h"
472: ! Input/output variables:
473: PetscScalar x(gxs:gxe,gys:gye)
474: Mat jac
475: integer ierr,ctx
476: DALocalInfo info(DA_LOCAL_INFO_SIZE)
477:
479: ! Local variables:
480: integer row,col(5),i,j
481: PetscScalar two,one,hx,hy,hxdhy,hydhx,sc,v(5)
483: ! Set parameters
485: one = 1.0
486: two = 2.0
487: hx = one/dble(mx-1)
488: hy = one/dble(my-1)
489: sc = hx*hy
490: hxdhy = hx/hy
491: hydhx = hy/hx
493: ! Compute entries for the locally owned part of the Jacobian.
494: ! - Currently, all PETSc parallel matrix formats are partitioned by
495: ! contiguous chunks of rows across the processors.
496: ! - Each processor needs to insert only elements that it owns
497: ! locally (but any non-local elements will be sent to the
498: ! appropriate processor during matrix assembly).
499: ! - Here, we set all entries for a particular row at once.
500: ! - We can set matrix entries either using either
501: ! MatSetValuesLocal() or MatSetValues(), as discussed above.
502: ! - Note that MatSetValues() uses 0-based row and column numbers
503: ! in Fortran as well as in C.
505: do 20 j=ys,ye
506: row = (j - gys)*gxm + xs - gxs - 1
507: do 10 i=xs,xe
508: row = row + 1
509: ! boundary points
510: if (i .eq. 1 .or. j .eq. 1 &
511: & .or. i .eq. mx .or. j .eq. my) then
512: call MatSetValuesLocal(jac,1,row,1,row,one, &
513: & INSERT_VALUES,ierr)
514: ! interior grid points
515: else
516: v(1) = -hxdhy
517: v(2) = -hydhx
518: v(3) = two*(hydhx + hxdhy) &
519: & - sc*lambda*exp(x(i,j))
520: v(4) = -hydhx
521: v(5) = -hxdhy
522: col(1) = row - gxm
523: col(2) = row - 1
524: col(3) = row
525: col(4) = row + 1
526: col(5) = row + gxm
527: call MatSetValuesLocal(jac,1,row,5,col,v, &
528: & INSERT_VALUES,ierr)
529: endif
530: 10 continue
531: 20 continue
532: call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
533: call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
534: return
535: end