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: mpiexec -n <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: PetscTruth flg
63: #if defined(PETSC_HAVE_ADIFOR)
64: ISColoring coloring
65: PetscTruth adifor_jacobian,adiformf_jacobian
66: #endif
68: ! Note: Any user-defined Fortran routines (such as FormJacobianLocal)
69: ! MUST be declared as external.
71: external FormInitialGuess
72: external FormFunctionLocal,FormJacobianLocal
73: #if defined(PETSC_HAVE_ADIFOR)
74: external g_FormFunctionLocal,m_FormFunctionLocal
75: #endif
77: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
78: ! Initialize program
79: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
81: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
82: call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
83: call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
85: ! Initialize problem parameters
87: i1 = 1
88: i4 = -4
89: lambda_max = 6.81
90: lambda_min = 0.0
91: lambda = 6.0
92: call PetscOptionsGetReal(PETSC_NULL_CHARACTER,'-par',lambda, &
93: & flg,ierr)
94: if (lambda .ge. lambda_max .or. lambda .le. lambda_min) then
95: if (rank .eq. 0) write(6,*) 'Lambda is out of range'
96: SETERRQ(1,' ',ierr)
97: endif
99: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100: ! Create nonlinear solver context
101: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103: call SNESCreate(PETSC_COMM_WORLD,snes,ierr)
105: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106: ! Create vector data structures; set function evaluation routine
107: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109: ! Create distributed array (DA) to manage parallel grid and vectors
111: ! This really needs only the star-type stencil, but we use the box
112: ! stencil temporarily.
113: call DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR, &
114: & i4,i4,PETSC_DECIDE,PETSC_DECIDE,i1,i1,PETSC_NULL_INTEGER, &
115: & PETSC_NULL_INTEGER,da,ierr)
117: ! Extract global and local vectors from DA; then duplicate for remaining
118: ! vectors that are the same types
120: call DACreateGlobalVector(da,x,ierr)
121: call VecDuplicate(x,r,ierr)
123: ! Get local grid boundaries (for 2-dimensional DA)
125: call DAGetInfo(da,PETSC_NULL_INTEGER,mx,my,PETSC_NULL_INTEGER, &
126: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
127: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
128: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
129: & PETSC_NULL_INTEGER,ierr)
130: call DAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym, &
131: & PETSC_NULL_INTEGER,ierr)
132: call DAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,gxm,gym, &
133: & PETSC_NULL_INTEGER,ierr)
135: ! Here we shift the starting indices up by one so that we can easily
136: ! use the Fortran convention of 1-based indices (rather 0-based indices).
138: xs = xs+1
139: ys = ys+1
140: gxs = gxs+1
141: gys = gys+1
143: ye = ys+ym-1
144: xe = xs+xm-1
145: gye = gys+gym-1
146: gxe = gxs+gxm-1
148: ! Set function evaluation routine and vector
150: call DASetLocalFunction(da,FormFunctionLocal,ierr)
151: call DASetLocalJacobian(da,FormJacobianLocal,ierr)
152: #if defined(PETSC_HAVE_ADIFOR)
153: call DASetLocalAdiforFunction(da, &
154: & g_FormFunctionLocal,ierr)
155: #endif
156: call SNESSetFunction(snes,r,SNESDAFormFunction,da,ierr)
157: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158: ! Create matrix data structure; set Jacobian evaluation routine
159: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161: ! Set Jacobian matrix data structure and default Jacobian evaluation
162: ! routine. User can override with:
163: ! -snes_mf : matrix-free Newton-Krylov method with no preconditioning
164: ! (unless user explicitly sets preconditioner)
165: ! -snes_mf_operator : form preconditioning matrix as set by the user,
166: ! but use matrix-free approx for Jacobian-vector
167: ! products within Newton-Krylov method
168: !
170: call DAGetMatrix(da,MATAIJ,J,ierr)
171: #if defined(PETSC_HAVE_ADIFOR)
172: call PetscOptionsGetTruth(PETSC_NULL_CHARACTER &
173: & ,'-adiformf_jacobian', &
174: & adiformf_jacobian,PETSC_NULL_INTEGER,ierr)
175: if (adiformf_jacobian .eq. 1) then
176: call DASetLocalAdiforMFFunction(da, &
177: & m_FormFunctionLocal,ierr)
178: call MatRegisterDAAD(ierr)
179: call MatCreateDAAD(da,A,ierr)
180: call MatDAADSetSNES(A,snes,ierr)
181: else
182: A = J
183: endif
184: #else
185: A = J
186: #endif
188: call SNESSetJacobian(snes,A,J,SNESDAComputeJacobian, &
189: & da,ierr)
190: #if defined(PETSC_HAVE_ADIFOR)
191: call PetscOptionsGetTruth(PETSC_NULL_CHARACTER &
192: & ,'-adifor_jacobian', &
193: & adifor_jacobian,PETSC_NULL_INTEGER,ierr)
194: if (adifor_jacobian .eq. 1) then
195: call DAGetColoring(da,IS_COLORING_GHOSTED, &
196: & coloring,ierr)
197: call MatSetColoring(J,coloring,ierr)
198: call SNESSetJacobian(snes,A,J,SNESDAComputeJacobianWithAdifor, &
199: & da,ierr)
200: call ISColoringDestroy(coloring,ierr)
201: endif
202: #endif
205: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
206: ! Customize nonlinear solver; set runtime options
207: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
209: ! Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
211: call SNESSetFromOptions(snes,ierr)
212: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
213: ! Evaluate initial guess; then solve nonlinear system.
214: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
216: ! Note: The user should initialize the vector, x, with the initial guess
217: ! for the nonlinear solver prior to calling SNESSolve(). In particular,
218: ! to employ an initial guess of zero, the user should explicitly set
219: ! this vector to zero by calling VecSet().
221: call FormInitialGuess(x,ierr)
222: call SNESSolve(snes,PETSC_NULL_OBJECT,x,ierr)
223: call SNESGetIterationNumber(snes,its,ierr)
224: if (rank .eq. 0) then
225: write(6,100) its
226: endif
227: 100 format('Number of Newton iterations = ',i5)
230: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
231: ! Free work space. All PETSc objects should be destroyed when they
232: ! are no longer needed.
233: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
235: if (A .ne. J) call MatDestroy(A,ierr)
236: call MatDestroy(J,ierr)
237: call VecDestroy(x,ierr)
238: call VecDestroy(r,ierr)
239: call SNESDestroy(snes,ierr)
240: call DADestroy(da,ierr)
241: call PetscFinalize(ierr)
242: end
244: ! ---------------------------------------------------------------------
245: !
246: ! FormInitialGuess - Forms initial approximation.
247: !
248: ! Input Parameters:
249: ! X - vector
250: !
251: ! Output Parameter:
252: ! X - vector
253: !
254: ! Notes:
255: ! This routine serves as a wrapper for the lower-level routine
256: ! "ApplicationInitialGuess", where the actual computations are
257: ! done using the standard Fortran style of treating the local
258: ! vector data as a multidimensional array over the local mesh.
259: ! This routine merely handles ghost point scatters and accesses
260: ! the local vector data via VecGetArray() and VecRestoreArray().
261: !
262: subroutine FormInitialGuess(X,ierr)
263: implicit none
265: #include "ex5f.h"
267: ! Input/output variables:
268: Vec X
269: PetscErrorCode ierr
271: ! Declarations for use with local arrays:
272: PetscScalar lx_v(0:1)
273: PetscOffset lx_i
274: Vec localX
276: 0
278: ! Get a pointer to vector data.
279: ! - For default PETSc vectors, VecGetArray() returns a pointer to
280: ! the data array. Otherwise, the routine is implementation dependent.
281: ! - You MUST call VecRestoreArray() when you no longer need access to
282: ! the array.
283: ! - Note that the Fortran interface to VecGetArray() differs from the
284: ! C version. See the users manual for details.
286: call DMGetLocalVector(da,localX,ierr)
287: call VecGetArray(localX,lx_v,lx_i,ierr)
289: ! Compute initial guess over the locally owned part of the grid
291: call InitialGuessLocal(lx_v(lx_i),ierr)
293: ! Restore vector
295: call VecRestoreArray(localX,lx_v,lx_i,ierr)
297: ! Insert values into global vector
299: call DALocalToGlobal(da,localX,INSERT_VALUES,X,ierr)
300: call DMRestoreLocalVector(da,localX,ierr)
301: return
302: end
304: ! ---------------------------------------------------------------------
305: !
306: ! InitialGuessLocal - Computes initial approximation, called by
307: ! the higher level routine FormInitialGuess().
308: !
309: ! Input Parameter:
310: ! x - local vector data
311: !
312: ! Output Parameters:
313: ! x - local vector data
314: ! ierr - error code
315: !
316: ! Notes:
317: ! This routine uses standard Fortran-style computations over a 2-dim array.
318: !
319: subroutine InitialGuessLocal(x,ierr)
320: implicit none
322: #include "ex5f.h"
324: ! Input/output variables:
325: PetscScalar x(gxs:gxe,gys:gye)
326: PetscErrorCode ierr
328: ! Local variables:
329: PetscInt i,j
330: PetscScalar temp1,temp,hx,hy
331: PetscScalar one
333: ! Set parameters
335: 0
336: one = 1.0
337: hx = one/(dble(mx-1))
338: hy = one/(dble(my-1))
339: temp1 = lambda/(lambda + one)
341: do 20 j=ys,ye
342: temp = dble(min(j-1,my-j))*hy
343: do 10 i=xs,xe
344: if (i .eq. 1 .or. j .eq. 1 &
345: & .or. i .eq. mx .or. j .eq. my) then
346: x(i,j) = 0.0
347: else
348: x(i,j) = temp1 * &
349: & sqrt(min(dble(min(i-1,mx-i)*hx),dble(temp)))
350: endif
351: 10 continue
352: 20 continue
354: return
355: end
357: ! ---------------------------------------------------------------------
358: !
359: ! FormFunctionLocal - Computes nonlinear function, called by
360: ! the higher level routine FormFunction().
361: !
362: ! Input Parameter:
363: ! x - local vector data
364: !
365: ! Output Parameters:
366: ! f - local vector data, f(x)
367: ! ierr - error code
368: !
369: ! Notes:
370: ! This routine uses standard Fortran-style computations over a 2-dim array.
371: !
372: ! Process adifor: FormFunctionLocal
373: !
374: subroutine FormFunctionLocal(info,x,f,dummy,ierr)
376: implicit none
378: #include "ex5f.h"
380: ! Input/output variables:
381: DALocalInfo info(DA_LOCAL_INFO_SIZE)
382: PetscScalar x(gxs:gxe,gys:gye)
383: PetscScalar f(xs:xe,ys:ye)
384: PetscErrorCode ierr
385: PetscObject dummy
387: ! Local variables:
388: PetscScalar two,one,hx,hy
389: PetscScalar hxdhy,hydhx,sc
390: PetscScalar u,uxx,uyy
391: PetscInt i,j
393: xs = info(DA_LOCAL_INFO_XS)+1
394: xe = xs+info(DA_LOCAL_INFO_XM)-1
395: ys = info(DA_LOCAL_INFO_YS)+1
396: ye = ys+info(DA_LOCAL_INFO_YM)-1
397: mx = info(DA_LOCAL_INFO_MX)
398: my = info(DA_LOCAL_INFO_MY)
400: one = 1.0
401: two = 2.0
402: hx = one/dble(mx-1)
403: hy = one/dble(my-1)
404: sc = hx*hy*lambda
405: hxdhy = hx/hy
406: hydhx = hy/hx
408: ! Compute function over the locally owned part of the grid
410: do 20 j=ys,ye
411: do 10 i=xs,xe
412: if (i .eq. 1 .or. j .eq. 1 &
413: & .or. i .eq. mx .or. j .eq. my) then
414: f(i,j) = x(i,j)
415: else
416: u = x(i,j)
417: uxx = hydhx * (two*u &
418: & - x(i-1,j) - x(i+1,j))
419: uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1))
420: f(i,j) = uxx + uyy - sc*exp(u)
421: endif
422: 10 continue
423: 20 continue
425: call PetscLogFlops(11*ym*xm,ierr)
427: return
428: end
430: ! ---------------------------------------------------------------------
431: !
432: ! FormJacobianLocal - Computes Jacobian matrix, called by
433: ! the higher level routine FormJacobian().
434: !
435: ! Input Parameters:
436: ! x - local vector data
437: !
438: ! Output Parameters:
439: ! jac - Jacobian matrix
440: ! jac_prec - optionally different preconditioning matrix (not used here)
441: ! ierr - error code
442: !
443: ! Notes:
444: ! This routine uses standard Fortran-style computations over a 2-dim array.
445: !
446: ! Notes:
447: ! Due to grid point reordering with DAs, we must always work
448: ! with the local grid points, and then transform them to the new
449: ! global numbering with the "ltog" mapping (via DAGetGlobalIndices()).
450: ! We cannot work directly with the global numbers for the original
451: ! uniprocessor grid!
452: !
453: ! Two methods are available for imposing this transformation
454: ! when setting matrix entries:
455: ! (A) MatSetValuesLocal(), using the local ordering (including
456: ! ghost points!)
457: ! - Use DAGetGlobalIndices() to extract the local-to-global map
458: ! - Associate this map with the matrix by calling
459: ! MatSetLocalToGlobalMapping() once
460: ! - Set matrix entries using the local ordering
461: ! by calling MatSetValuesLocal()
462: ! (B) MatSetValues(), using the global ordering
463: ! - Use DAGetGlobalIndices() to extract the local-to-global map
464: ! - Then apply this map explicitly yourself
465: ! - Set matrix entries using the global ordering by calling
466: ! MatSetValues()
467: ! Option (A) seems cleaner/easier in many cases, and is the procedure
468: ! used in this example.
469: !
470: subroutine FormJacobianLocal(info,x,jac,ctx,ierr)
471: implicit none
473: #include "ex5f.h"
475: ! Input/output variables:
476: PetscScalar x(gxs:gxe,gys:gye)
477: Mat jac
478: PetscErrorCode ierr
479: integer ctx
480: DALocalInfo info(DA_LOCAL_INFO_SIZE)
481:
483: ! Local variables:
484: PetscInt row,col(5),i,j,i1,i5
485: PetscScalar two,one,hx,hy,v(5)
486: PetscScalar hxdhy,hydhx,sc
488: ! Set parameters
490: i1 = 1
491: i5 = 5
492: one = 1.0
493: two = 2.0
494: hx = one/dble(mx-1)
495: hy = one/dble(my-1)
496: sc = hx*hy
497: hxdhy = hx/hy
498: hydhx = hy/hx
500: ! Compute entries for the locally owned part of the Jacobian.
501: ! - Currently, all PETSc parallel matrix formats are partitioned by
502: ! contiguous chunks of rows across the processors.
503: ! - Each processor needs to insert only elements that it owns
504: ! locally (but any non-local elements will be sent to the
505: ! appropriate processor during matrix assembly).
506: ! - Here, we set all entries for a particular row at once.
507: ! - We can set matrix entries either using either
508: ! MatSetValuesLocal() or MatSetValues(), as discussed above.
509: ! - Note that MatSetValues() uses 0-based row and column numbers
510: ! in Fortran as well as in C.
512: do 20 j=ys,ye
513: row = (j - gys)*gxm + xs - gxs - 1
514: do 10 i=xs,xe
515: row = row + 1
516: ! boundary points
517: if (i .eq. 1 .or. j .eq. 1 &
518: & .or. i .eq. mx .or. j .eq. my) then
519: ! Some f90 compilers need 4th arg to be of same type in both calls
520: col(1) = row
521: v(1) = one
522: call MatSetValuesLocal(jac,i1,row,i1,col,v, &
523: & INSERT_VALUES,ierr)
524: ! interior grid points
525: else
526: v(1) = -hxdhy
527: v(2) = -hydhx
528: v(3) = two*(hydhx + hxdhy) &
529: & - sc*lambda*exp(x(i,j))
530: v(4) = -hydhx
531: v(5) = -hxdhy
532: col(1) = row - gxm
533: col(2) = row - 1
534: col(3) = row
535: col(4) = row + 1
536: col(5) = row + gxm
537: call MatSetValuesLocal(jac,i1,row,i5,col,v, &
538: & INSERT_VALUES,ierr)
539: endif
540: 10 continue
541: 20 continue
542: call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
543: call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
544: return
545: end