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) && !defined(PETSC_USE_COMPLEX)
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) && !defined(PETSC_USE_COMPLEX)
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,MATAIJ,J,ierr)
168: #if defined(PETSC_HAVE_ADIFOR) && !defined(PETSC_USE_COMPLEX)
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) && !defined(PETSC_USE_COMPLEX)
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 DAGetColoring(da,IS_COLORING_GHOSTED, &
194: & coloring,ierr)
195: call MatSetColoring(J,coloring,ierr)
196: call SNESSetJacobian(snes,A,J,SNESDAComputeJacobianWithAdifor, &
197: & da,ierr)
198: call ISColoringDestroy(coloring,ierr)
199: endif
200: #endif
203: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204: ! Customize nonlinear solver; set runtime options
205: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
207: ! Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
209: call SNESSetFromOptions(snes,ierr)
211: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
212: ! Evaluate initial guess; then solve nonlinear system.
213: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
215: ! Note: The user should initialize the vector, x, with the initial guess
216: ! for the nonlinear solver prior to calling SNESSolve(). In particular,
217: ! to employ an initial guess of zero, the user should explicitly set
218: ! this vector to zero by calling VecSet().
220: call FormInitialGuess(x,ierr)
221: call SNESSolve(snes,x,ierr)
222: call SNESGetIterationNumber(snes,its,ierr);
223: if (rank .eq. 0) then
224: write(6,100) its
225: endif
226: 100 format('Number of Newton iterations = ',i5)
229: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
230: ! Free work space. All PETSc objects should be destroyed when they
231: ! are no longer needed.
232: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
234: if (A .ne. J) call MatDestroy(A,ierr)
235: call MatDestroy(J,ierr)
236: call VecDestroy(x,ierr)
237: call VecDestroy(r,ierr)
238: call SNESDestroy(snes,ierr)
239: call DADestroy(da,ierr)
240: call PetscFinalize(ierr)
241: end
243: ! ---------------------------------------------------------------------
244: !
245: ! FormInitialGuess - Forms initial approximation.
246: !
247: ! Input Parameters:
248: ! X - vector
249: !
250: ! Output Parameter:
251: ! X - vector
252: !
253: ! Notes:
254: ! This routine serves as a wrapper for the lower-level routine
255: ! "ApplicationInitialGuess", where the actual computations are
256: ! done using the standard Fortran style of treating the local
257: ! vector data as a multidimensional array over the local mesh.
258: ! This routine merely handles ghost point scatters and accesses
259: ! the local vector data via VecGetArray() and VecRestoreArray().
260: !
261: subroutine FormInitialGuess(X,ierr)
262: implicit none
264: #include "ex5f.h"
266: ! Input/output variables:
267: Vec X
268: integer ierr
270: ! Declarations for use with local arrays:
271: PetscScalar lx_v(0:1)
272: PetscOffset lx_i
273: Vec localX
275: 0
277: ! Get a pointer to vector data.
278: ! - For default PETSc vectors, VecGetArray() returns a pointer to
279: ! the data array. Otherwise, the routine is implementation dependent.
280: ! - You MUST call VecRestoreArray() when you no longer need access to
281: ! the array.
282: ! - Note that the Fortran interface to VecGetArray() differs from the
283: ! C version. See the users manual for details.
285: call DAGetLocalVector(da,localX,ierr)
286: call VecGetArray(localX,lx_v,lx_i,ierr)
288: ! Compute initial guess over the locally owned part of the grid
290: call InitialGuessLocal(lx_v(lx_i),ierr)
292: ! Restore vector
294: call VecRestoreArray(localX,lx_v,lx_i,ierr)
296: ! Insert values into global vector
298: call DALocalToGlobal(da,localX,INSERT_VALUES,X,ierr)
299: call DARestoreLocalVector(da,localX,ierr)
300: return
301: end
303: ! ---------------------------------------------------------------------
304: !
305: ! InitialGuessLocal - Computes initial approximation, called by
306: ! the higher level routine FormInitialGuess().
307: !
308: ! Input Parameter:
309: ! x - local vector data
310: !
311: ! Output Parameters:
312: ! x - local vector data
313: ! ierr - error code
314: !
315: ! Notes:
316: ! This routine uses standard Fortran-style computations over a 2-dim array.
317: !
318: subroutine InitialGuessLocal(x,ierr)
319: implicit none
321: #include "ex5f.h"
323: ! Input/output variables:
324: PetscScalar x(gxs:gxe,gys:gye)
325: integer ierr
327: ! Local variables:
328: integer i,j,hxdhy,hydhx
329: PetscScalar temp1,temp,hx,hy,sc,one
331: ! Set parameters
333: 0
334: one = 1.0
335: hx = one/(dble(mx-1))
336: hy = one/(dble(my-1))
337: sc = hx*hy*lambda
338: hxdhy = hx/hy
339: hydhx = hy/hx
340: temp1 = lambda/(lambda + one)
342: do 20 j=ys,ye
343: temp = dble(min(j-1,my-j))*hy
344: do 10 i=xs,xe
345: if (i .eq. 1 .or. j .eq. 1 &
346: & .or. i .eq. mx .or. j .eq. my) then
347: x(i,j) = 0.0
348: else
349: x(i,j) = temp1 * &
350: & sqrt(min(dble(min(i-1,mx-i)*hx),dble(temp)))
351: endif
352: 10 continue
353: 20 continue
355: return
356: end
358: ! ---------------------------------------------------------------------
359: !
360: ! FormFunctionLocal - Computes nonlinear function, called by
361: ! the higher level routine FormFunction().
362: !
363: ! Input Parameter:
364: ! x - local vector data
365: !
366: ! Output Parameters:
367: ! f - local vector data, f(x)
368: ! ierr - error code
369: !
370: ! Notes:
371: ! This routine uses standard Fortran-style computations over a 2-dim array.
372: !
373: ! Process adifor: FormFunctionLocal
374: !
375: subroutine FormFunctionLocal(info,x,f,dummy,ierr)
377: implicit none
379: #include "ex5f.h"
381: ! Input/output variables:
382: DALocalInfo info(DA_LOCAL_INFO_SIZE)
383: PetscScalar x(gxs:gxe,gys:gye)
384: PetscScalar f(xs:xe,ys:ye)
385: integer ierr
386: PetscObject dummy
388: ! Local variables:
389: PetscScalar two,one,hx,hy,hxdhy,hydhx,sc
390: PetscScalar u,uxx,uyy
391: integer i,j
394: xs = info(DA_LOCAL_INFO_XS)+1
395: xe = xs+info(DA_LOCAL_INFO_XM)-1
396: ys = info(DA_LOCAL_INFO_YS)+1
397: ye = ys+info(DA_LOCAL_INFO_YM)-1
398: mx = info(DA_LOCAL_INFO_MX)
399: my = info(DA_LOCAL_INFO_MY)
401: one = 1.0
402: two = 2.0
403: hx = one/dble(mx-1)
404: hy = one/dble(my-1)
405: sc = hx*hy*lambda
406: hxdhy = hx/hy
407: hydhx = hy/hx
409: ! Compute function over the locally owned part of the grid
411: do 20 j=ys,ye
412: do 10 i=xs,xe
413: if (i .eq. 1 .or. j .eq. 1 &
414: & .or. i .eq. mx .or. j .eq. my) then
415: f(i,j) = x(i,j)
416: else
417: u = x(i,j)
418: uxx = hydhx * (two*u &
419: & - x(i-1,j) - x(i+1,j))
420: uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1))
421: f(i,j) = uxx + uyy - sc*exp(u)
422: endif
423: 10 continue
424: 20 continue
426: call PetscLogFlops(11*ym*xm,ierr)
428: return
429: end
431: ! ---------------------------------------------------------------------
432: !
433: ! FormJacobianLocal - Computes Jacobian matrix, called by
434: ! the higher level routine FormJacobian().
435: !
436: ! Input Parameters:
437: ! x - local vector data
438: !
439: ! Output Parameters:
440: ! jac - Jacobian matrix
441: ! jac_prec - optionally different preconditioning matrix (not used here)
442: ! ierr - error code
443: !
444: ! Notes:
445: ! This routine uses standard Fortran-style computations over a 2-dim array.
446: !
447: ! Notes:
448: ! Due to grid point reordering with DAs, we must always work
449: ! with the local grid points, and then transform them to the new
450: ! global numbering with the "ltog" mapping (via DAGetGlobalIndices()).
451: ! We cannot work directly with the global numbers for the original
452: ! uniprocessor grid!
453: !
454: ! Two methods are available for imposing this transformation
455: ! when setting matrix entries:
456: ! (A) MatSetValuesLocal(), using the local ordering (including
457: ! ghost points!)
458: ! - Use DAGetGlobalIndices() to extract the local-to-global map
459: ! - Associate this map with the matrix by calling
460: ! MatSetLocalToGlobalMapping() once
461: ! - Set matrix entries using the local ordering
462: ! by calling MatSetValuesLocal()
463: ! (B) MatSetValues(), using the global ordering
464: ! - Use DAGetGlobalIndices() to extract the local-to-global map
465: ! - Then apply this map explicitly yourself
466: ! - Set matrix entries using the global ordering by calling
467: ! MatSetValues()
468: ! Option (A) seems cleaner/easier in many cases, and is the procedure
469: ! used in this example.
470: !
471: subroutine FormJacobianLocal(info,x,jac,ctx,ierr)
472: implicit none
474: #include "ex5f.h"
476: ! Input/output variables:
477: PetscScalar x(gxs:gxe,gys:gye)
478: Mat jac
479: integer ierr,ctx
480: DALocalInfo info(DA_LOCAL_INFO_SIZE)
481:
483: ! Local variables:
484: integer row,col(5),i,j
485: PetscScalar two,one,hx,hy,hxdhy,hydhx,sc,v(5)
487: ! Set parameters
489: one = 1.0
490: two = 2.0
491: hx = one/dble(mx-1)
492: hy = one/dble(my-1)
493: sc = hx*hy
494: hxdhy = hx/hy
495: hydhx = hy/hx
497: ! Compute entries for the locally owned part of the Jacobian.
498: ! - Currently, all PETSc parallel matrix formats are partitioned by
499: ! contiguous chunks of rows across the processors.
500: ! - Each processor needs to insert only elements that it owns
501: ! locally (but any non-local elements will be sent to the
502: ! appropriate processor during matrix assembly).
503: ! - Here, we set all entries for a particular row at once.
504: ! - We can set matrix entries either using either
505: ! MatSetValuesLocal() or MatSetValues(), as discussed above.
506: ! - Note that MatSetValues() uses 0-based row and column numbers
507: ! in Fortran as well as in C.
509: do 20 j=ys,ye
510: row = (j - gys)*gxm + xs - gxs - 1
511: do 10 i=xs,xe
512: row = row + 1
513: ! boundary points
514: if (i .eq. 1 .or. j .eq. 1 &
515: & .or. i .eq. mx .or. j .eq. my) then
516: ! Some f90 compilers need 4th arg to be of same type in both calls
517: col(1) = row
518: v(1) = one
519: call MatSetValuesLocal(jac,1,row,1,col,v, &
520: & INSERT_VALUES,ierr)
521: ! interior grid points
522: else
523: v(1) = -hxdhy
524: v(2) = -hydhx
525: v(3) = two*(hydhx + hxdhy) &
526: & - sc*lambda*exp(x(i,j))
527: v(4) = -hydhx
528: v(5) = -hxdhy
529: col(1) = row - gxm
530: col(2) = row - 1
531: col(3) = row
532: col(4) = row + 1
533: col(5) = row + gxm
534: call MatSetValuesLocal(jac,1,row,5,col,v, &
535: & INSERT_VALUES,ierr)
536: endif
537: 10 continue
538: 20 continue
539: call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
540: call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
541: return
542: end