Actual source code: ex5f.F
1: !
2: ! Description: This example solves a nonlinear system in parallel with SNES.
3: ! We solve the Bratu (SFI - solid fuel ignition) problem in a 2D rectangular
4: ! domain, using distributed arrays (DAs) to partition the parallel grid.
5: ! The command line options include:
6: ! -par <param>, where <param> indicates the nonlinearity of the problem
7: ! problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81)
8: !
9: ! Program usage: mpirun -np <procs> ex5f [-help] [all PETSc options]
10: !
11: !/*T
12: ! Concepts: SNES^parallel Bratu example
13: ! Concepts: DA^using distributed arrays;
14: ! Processors: n
15: !T*/
16: !
17: ! --------------------------------------------------------------------------
18: !
19: ! Solid Fuel Ignition (SFI) problem. This problem is modeled by
20: ! the partial differential equation
21: !
22: ! -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1,
23: !
24: ! with boundary conditions
25: !
26: ! u = 0 for x = 0, x = 1, y = 0, y = 1.
27: !
28: ! A finite difference approximation with the usual 5-point stencil
29: ! is used to discretize the boundary value problem to obtain a nonlinear
30: ! system of equations.
31: !
32: ! --------------------------------------------------------------------------
34: program main
35: implicit none
36: !
37: ! We place common blocks, variable declarations, and other include files
38: ! needed for this code in the single file ex5f.h. We then need to include
39: ! only this file throughout the various routines in this program. See
40: ! additional comments in the file ex5f.h.
41: !
42: #include "ex5f.h"
44: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
45: ! Variable declarations
46: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47: !
48: ! Variables:
49: ! snes - nonlinear solver
50: ! x, r - solution, residual vectors
51: ! J - Jacobian matrix
52: ! its - iterations for convergence
53: !
54: ! See additional variable declarations in the file ex5f.h
55: !
56: SNES snes
57: Vec x,r
58: Mat J,A
59: PetscInt its,i1,i4
60: PetscErrorCode ierr
61: double precision lambda_max,lambda_min
62: ISColoring coloring
63: PetscTruth adifor_jacobian,adiformf_jacobian,flg
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: i1 = 1
85: i4 = -4
86: lambda_max = 6.81
87: lambda_min = 0.0
88: lambda = 6.0
89: call PetscOptionsGetReal(PETSC_NULL_CHARACTER,'-par',lambda, &
90: & flg,ierr)
91: if (lambda .ge. lambda_max .or. lambda .le. lambda_min) then
92: if (rank .eq. 0) write(6,*) 'Lambda is out of range'
93: SETERRQ(1,' ',ierr)
94: endif
96: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97: ! Create nonlinear solver context
98: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100: call SNESCreate(PETSC_COMM_WORLD,snes,ierr)
102: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103: ! Create vector data structures; set function evaluation routine
104: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106: ! Create distributed array (DA) to manage parallel grid and vectors
108: ! This really needs only the star-type stencil, but we use the box
109: ! stencil temporarily.
110: call DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR, &
111: & i4,i4,PETSC_DECIDE,PETSC_DECIDE,i1,i1,PETSC_NULL_INTEGER, &
112: & PETSC_NULL_INTEGER,da,ierr)
114: ! Extract global and local vectors from DA; then duplicate for remaining
115: ! vectors that are the same types
117: call DACreateGlobalVector(da,x,ierr)
118: call VecDuplicate(x,r,ierr)
120: ! Get local grid boundaries (for 2-dimensional DA)
122: call DAGetInfo(da,PETSC_NULL_INTEGER,mx,my,PETSC_NULL_INTEGER, &
123: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
124: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
125: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
126: & PETSC_NULL_INTEGER,ierr)
127: call DAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym, &
128: & PETSC_NULL_INTEGER,ierr)
129: call DAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,gxm,gym, &
130: & PETSC_NULL_INTEGER,ierr)
132: ! Here we shift the starting indices up by one so that we can easily
133: ! use the Fortran convention of 1-based indices (rather 0-based indices).
135: xs = xs+1
136: ys = ys+1
137: gxs = gxs+1
138: gys = gys+1
140: ye = ys+ym-1
141: xe = xs+xm-1
142: gye = gys+gym-1
143: gxe = gxs+gxm-1
145: ! Set function evaluation routine and vector
147: call DASetLocalFunction(da,FormFunctionLocal,ierr)
148: call DASetLocalJacobian(da,FormJacobianLocal,ierr)
149: #if defined(PETSC_HAVE_ADIFOR) && !defined(PETSC_USE_COMPLEX)
150: call DASetLocalAdiforFunction(da, &
151: & g_FormFunctionLocal,ierr)
152: #endif
153: 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,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)
187: #if defined(PETSC_HAVE_ADIFOR) && !defined(PETSC_USE_COMPLEX)
188: call PetscOptionsGetLogical(PETSC_NULL_CHARACTER &
189: & ,'-adifor_jacobian', &
190: & adifor_jacobian,PETSC_NULL_INTEGER,ierr)
191: if (adifor_jacobian .eq. 1) then
192: call DAGetColoring(da,IS_COLORING_GHOSTED, &
193: & coloring,ierr)
194: call MatSetColoring(J,coloring,ierr)
195: call SNESSetJacobian(snes,A,J,SNESDAComputeJacobianWithAdifor, &
196: & da,ierr)
197: call ISColoringDestroy(coloring,ierr)
198: endif
199: #endif
202: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
203: ! Customize nonlinear solver; set runtime options
204: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
206: ! Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
208: call SNESSetFromOptions(snes,ierr)
209: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
210: ! Evaluate initial guess; then solve nonlinear system.
211: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
213: ! Note: The user should initialize the vector, x, with the initial guess
214: ! for the nonlinear solver prior to calling SNESSolve(). In particular,
215: ! to employ an initial guess of zero, the user should explicitly set
216: ! this vector to zero by calling VecSet().
218: call FormInitialGuess(x,ierr)
219: call SNESSolve(snes,x,ierr)
220: call SNESGetIterationNumber(snes,its,ierr);
221: if (rank .eq. 0) then
222: write(6,100) its
223: endif
224: 100 format('Number of Newton iterations = ',i5)
227: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
228: ! Free work space. All PETSc objects should be destroyed when they
229: ! are no longer needed.
230: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
232: if (A .ne. J) call MatDestroy(A,ierr)
233: call MatDestroy(J,ierr)
234: call VecDestroy(x,ierr)
235: call VecDestroy(r,ierr)
236: call SNESDestroy(snes,ierr)
237: call DADestroy(da,ierr)
238: call PetscFinalize(ierr)
239: end
241: ! ---------------------------------------------------------------------
242: !
243: ! FormInitialGuess - Forms initial approximation.
244: !
245: ! Input Parameters:
246: ! X - vector
247: !
248: ! Output Parameter:
249: ! X - vector
250: !
251: ! Notes:
252: ! This routine serves as a wrapper for the lower-level routine
253: ! "ApplicationInitialGuess", where the actual computations are
254: ! done using the standard Fortran style of treating the local
255: ! vector data as a multidimensional array over the local mesh.
256: ! This routine merely handles ghost point scatters and accesses
257: ! the local vector data via VecGetArray() and VecRestoreArray().
258: !
259: subroutine FormInitialGuess(X,ierr)
260: implicit none
262: #include "ex5f.h"
264: ! Input/output variables:
265: Vec X
266: PetscErrorCode ierr
268: ! Declarations for use with local arrays:
269: PetscScalar lx_v(0:1)
270: PetscOffset lx_i
271: Vec localX
273: 0
275: ! Get a pointer to vector data.
276: ! - For default PETSc vectors, VecGetArray() returns a pointer to
277: ! the data array. Otherwise, the routine is implementation dependent.
278: ! - You MUST call VecRestoreArray() when you no longer need access to
279: ! the array.
280: ! - Note that the Fortran interface to VecGetArray() differs from the
281: ! C version. See the users manual for details.
283: call DAGetLocalVector(da,localX,ierr)
284: call VecGetArray(localX,lx_v,lx_i,ierr)
286: ! Compute initial guess over the locally owned part of the grid
288: call InitialGuessLocal(lx_v(lx_i),ierr)
290: ! Restore vector
292: call VecRestoreArray(localX,lx_v,lx_i,ierr)
294: ! Insert values into global vector
296: call DALocalToGlobal(da,localX,INSERT_VALUES,X,ierr)
297: call DARestoreLocalVector(da,localX,ierr)
298: return
299: end
301: ! ---------------------------------------------------------------------
302: !
303: ! InitialGuessLocal - Computes initial approximation, called by
304: ! the higher level routine FormInitialGuess().
305: !
306: ! Input Parameter:
307: ! x - local vector data
308: !
309: ! Output Parameters:
310: ! x - local vector data
311: ! ierr - error code
312: !
313: ! Notes:
314: ! This routine uses standard Fortran-style computations over a 2-dim array.
315: !
316: subroutine InitialGuessLocal(x,ierr)
317: implicit none
319: #include "ex5f.h"
321: ! Input/output variables:
322: PetscScalar x(gxs:gxe,gys:gye)
323: PetscErrorCode ierr
325: ! Local variables:
326: PetscInt i,j
327: PetscScalar temp1,temp,hx,hy,sc,one,hxdhy,hydhx
329: ! Set parameters
331: 0
332: one = 1.0
333: hx = one/(dble(mx-1))
334: hy = one/(dble(my-1))
335: sc = hx*hy*lambda
336: hxdhy = hx/hy
337: hydhx = hy/hx
338: temp1 = lambda/(lambda + one)
340: do 20 j=ys,ye
341: temp = dble(min(j-1,my-j))*hy
342: do 10 i=xs,xe
343: if (i .eq. 1 .or. j .eq. 1 &
344: & .or. i .eq. mx .or. j .eq. my) then
345: x(i,j) = 0.0
346: else
347: x(i,j) = temp1 * &
348: & sqrt(min(dble(min(i-1,mx-i)*hx),dble(temp)))
349: endif
350: 10 continue
351: 20 continue
353: return
354: end
356: ! ---------------------------------------------------------------------
357: !
358: ! FormFunctionLocal - Computes nonlinear function, called by
359: ! the higher level routine FormFunction().
360: !
361: ! Input Parameter:
362: ! x - local vector data
363: !
364: ! Output Parameters:
365: ! f - local vector data, f(x)
366: ! ierr - error code
367: !
368: ! Notes:
369: ! This routine uses standard Fortran-style computations over a 2-dim array.
370: !
371: ! Process adifor: FormFunctionLocal
372: !
373: subroutine FormFunctionLocal(info,x,f,dummy,ierr)
375: implicit none
377: #include "ex5f.h"
379: ! Input/output variables:
380: DALocalInfo info(DA_LOCAL_INFO_SIZE)
381: PetscScalar x(gxs:gxe,gys:gye)
382: PetscScalar f(xs:xe,ys:ye)
383: PetscErrorCode ierr
384: PetscObject dummy
386: ! Local variables:
387: PetscScalar two,one,hx,hy,hxdhy,hydhx,sc
388: PetscScalar u,uxx,uyy
389: PetscInt i,j
391: xs = info(DA_LOCAL_INFO_XS)+1
392: xe = xs+info(DA_LOCAL_INFO_XM)-1
393: ys = info(DA_LOCAL_INFO_YS)+1
394: ye = ys+info(DA_LOCAL_INFO_YM)-1
395: mx = info(DA_LOCAL_INFO_MX)
396: my = info(DA_LOCAL_INFO_MY)
398: one = 1.0
399: two = 2.0
400: hx = one/dble(mx-1)
401: hy = one/dble(my-1)
402: sc = hx*hy*lambda
403: hxdhy = hx/hy
404: hydhx = hy/hx
406: ! Compute function over the locally owned part of the grid
408: do 20 j=ys,ye
409: do 10 i=xs,xe
410: if (i .eq. 1 .or. j .eq. 1 &
411: & .or. i .eq. mx .or. j .eq. my) then
412: f(i,j) = x(i,j)
413: else
414: u = x(i,j)
415: uxx = hydhx * (two*u &
416: & - x(i-1,j) - x(i+1,j))
417: uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1))
418: f(i,j) = uxx + uyy - sc*exp(u)
419: endif
420: 10 continue
421: 20 continue
423: call PetscLogFlops(11*ym*xm,ierr)
425: return
426: end
428: ! ---------------------------------------------------------------------
429: !
430: ! FormJacobianLocal - Computes Jacobian matrix, called by
431: ! the higher level routine FormJacobian().
432: !
433: ! Input Parameters:
434: ! x - local vector data
435: !
436: ! Output Parameters:
437: ! jac - Jacobian matrix
438: ! jac_prec - optionally different preconditioning matrix (not used here)
439: ! ierr - error code
440: !
441: ! Notes:
442: ! This routine uses standard Fortran-style computations over a 2-dim array.
443: !
444: ! Notes:
445: ! Due to grid point reordering with DAs, we must always work
446: ! with the local grid points, and then transform them to the new
447: ! global numbering with the "ltog" mapping (via DAGetGlobalIndices()).
448: ! We cannot work directly with the global numbers for the original
449: ! uniprocessor grid!
450: !
451: ! Two methods are available for imposing this transformation
452: ! when setting matrix entries:
453: ! (A) MatSetValuesLocal(), using the local ordering (including
454: ! ghost points!)
455: ! - Use DAGetGlobalIndices() to extract the local-to-global map
456: ! - Associate this map with the matrix by calling
457: ! MatSetLocalToGlobalMapping() once
458: ! - Set matrix entries using the local ordering
459: ! by calling MatSetValuesLocal()
460: ! (B) MatSetValues(), using the global ordering
461: ! - Use DAGetGlobalIndices() to extract the local-to-global map
462: ! - Then apply this map explicitly yourself
463: ! - Set matrix entries using the global ordering by calling
464: ! MatSetValues()
465: ! Option (A) seems cleaner/easier in many cases, and is the procedure
466: ! used in this example.
467: !
468: subroutine FormJacobianLocal(info,x,jac,ctx,ierr)
469: implicit none
471: #include "ex5f.h"
473: ! Input/output variables:
474: PetscScalar x(gxs:gxe,gys:gye)
475: Mat jac
476: PetscErrorCode ierr
477: integer ctx
478: DALocalInfo info(DA_LOCAL_INFO_SIZE)
479:
481: ! Local variables:
482: PetscInt row,col(5),i,j,i1,i5
483: PetscScalar two,one,hx,hy,hxdhy,hydhx,sc,v(5)
485: ! Set parameters
487: i1 = 1
488: i5 = 5
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,i1,row,i1,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,i1,row,i5,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