1: !
2: ! Description: 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 (DMDAs) to partition the parallel grid.
5: ! The command line options include:
6: ! -par <parameter>, where <parameter> indicates the nonlinearity of the problem
7: ! problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81)
8: !
9: !/*T
10: ! Concepts: SNES^parallel Bratu example
11: ! Concepts: DMDA^using distributed arrays;
12: ! Processors: n
13: !T*/
14: !
15: ! --------------------------------------------------------------------------
16: !
17: ! Solid Fuel Ignition (SFI) problem. This problem is modeled by
18: ! the partial differential equation
19: !
20: ! -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1,
21: !
22: ! with boundary conditions
23: !
24: ! u = 0 for x = 0, x = 1, y = 0, y = 1.
25: !
26: ! A finite difference approximation with the usual 5-point stencil
27: ! is used to discretize the boundary value problem to obtain a nonlinear
28: ! system of equations.
29: !
30: ! The uniprocessor version of this code is snes/examples/tutorials/ex4f.F
31: !
32: ! --------------------------------------------------------------------------
33: ! The following define must be used before including any PETSc include files
34: ! into a module or interface. This is because they can't handle declarations
35: ! in them
36: !
38: module f90module
39: #include <finclude/petscdmdef.h>
40: use petscdmdef
41: type userctx
42: type(DM) da
43: PetscInt xs,xe,xm,gxs,gxe,gxm
44: PetscInt ys,ye,ym,gys,gye,gym
45: PetscInt mx,my
46: PetscMPIInt rank
47: PetscReal lambda
48: end type userctx
50: contains
51: ! ---------------------------------------------------------------------
52: !
53: ! FormFunction - Evaluates nonlinear function, F(x).
54: !
55: ! Input Parameters:
56: ! snes - the SNES context
57: ! X - input vector
58: ! dummy - optional user-defined context, as set by SNESSetFunction()
59: ! (not used here)
60: !
61: ! Output Parameter:
62: ! F - function vector
63: !
64: ! Notes:
65: ! This routine serves as a wrapper for the lower-level routine
66: ! "FormFunctionLocal", where the actual computations are
67: ! done using the standard Fortran style of treating the local
68: ! vector data as a multidimensional array over the local mesh.
69: ! This routine merely handles ghost point scatters and accesses
70: ! the local vector data via VecGetArrayF90() and VecRestoreArrayF90().
71: !
72: subroutine FormFunction(snesIn,X,F,user,ierr)
73: #include <finclude/petscsnesdef.h>
74: use petscsnes
76: ! Input/output variables:
77: type(SNES) snesIn
78: type(Vec) X,F
79: PetscErrorCode ierr
80: type (userctx) user
82: ! Declarations for use with local arrays:
83: PetscScalar,pointer :: lx_v(:),lf_v(:)
84: type(Vec) localX
86: ! Scatter ghost points to local vector, using the 2-step process
87: ! DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
88: ! By placing code between these two statements, computations can
89: ! be done while messages are in transition.
90: call DMGetLocalVector(user%da,localX,ierr)
91: call DMGlobalToLocalBegin(user%da,X,INSERT_VALUES, &
92: & localX,ierr)
93: call DMGlobalToLocalEnd(user%da,X,INSERT_VALUES,localX,ierr)
95: ! Get a pointer to vector data.
96: ! - For default PETSc vectors, VecGetArray90() returns a pointer to
97: ! the data array. Otherwise, the routine is implementation dependent.
98: ! - You MUST call VecRestoreArrayF90() when you no longer need access to
99: ! the array.
100: ! - Note that the interface to VecGetArrayF90() differs from VecGetArray(),
101: ! and is useable from Fortran-90 Only.
103: call VecGetArrayF90(localX,lx_v,ierr)
104: call VecGetArrayF90(F,lf_v,ierr)
106: ! Compute function over the locally owned part of the grid
107: call FormFunctionLocal(lx_v,lf_v,user,ierr)
109: ! Restore vectors
110: call VecRestoreArrayF90(localX,lx_v,ierr)
111: call VecRestoreArrayF90(F,lf_v,ierr)
113: ! Insert values into global vector
115: call DMRestoreLocalVector(user%da,localX,ierr)
116: call PetscLogFlops(11.0d0*user%ym*user%xm,ierr)
118: ! call VecView(X,PETSC_VIEWER_STDOUT_WORLD,ierr)
119: ! call VecView(F,PETSC_VIEWER_STDOUT_WORLD,ierr)
120: return
121: end subroutine formfunction
122: end module f90module
124: module f90moduleinterfaces
125: use f90module
127: Interface SNESSetApplicationContext128: Subroutine SNESSetApplicationContext(snesIn,ctx,ierr)
129: #include <finclude/petscsnesdef.h>
130: use petscsnes
131: use f90module
132: type(SNES) snesIn
133: type(userctx) ctx
134: PetscErrorCode ierr
135: End Subroutine
136: End Interface SNESSetApplicationContext138: Interface SNESGetApplicationContext139: Subroutine SNESGetApplicationContext(snesIn,ctx,ierr)
140: #include <finclude/petscsnesdef.h>
141: use petscsnes
142: use f90module
143: type(SNES) snesIn
144: type(userctx), pointer :: ctx
145: PetscErrorCode ierr
146: End Subroutine
147: End Interface SNESGetApplicationContext148: end module f90moduleinterfaces
150: program main
151: #include <finclude/petscdmdef.h>
152: #include <finclude/petscsnesdef.h>
153: use petscdm
154: use petscsnes
155: use f90module
156: use f90moduleinterfaces
157: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158: ! Variable declarations
159: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
160: !
161: ! Variables:
162: ! mysnes - nonlinear solver
163: ! x, r - solution, residual vectors
164: ! J - Jacobian matrix
165: ! its - iterations for convergence
166: ! Nx, Ny - number of preocessors in x- and y- directions
167: ! matrix_free - flag - 1 indicates matrix-free version
168: !
169: type(SNES) mysnes
170: type(Vec) x,r
171: type(Mat) J
172: PetscErrorCode ierr
173: PetscInt its
174: PetscBool flg,matrix_free
175: PetscInt ione,nfour
176: PetscReal lambda_max,lambda_min
177: type(userctx) user
178: type(userctx), pointer:: puser
180: ! Note: Any user-defined Fortran routines (such as FormJacobian)
181: ! MUST be declared as external.
182: external FormInitialGuess,FormJacobian
184: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185: ! Initialize program
186: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
188: call MPI_Comm_rank(PETSC_COMM_WORLD,user%rank,ierr)
190: ! Initialize problem parameters
191: lambda_max = 6.81
192: lambda_min = 0.0
193: user%lambda = 6.0
194: ione = 1
195: nfour = -4
196: call PetscOptionsGetReal(PETSC_NULL_CHARACTER,'-par', &
197: & user%lambda,flg,ierr)
198: if (user%lambda .ge. lambda_max .or. user%lambda .le. lambda_min) &
199: & then
200: if (user%rank .eq. 0) write(6,*) 'Lambda is out of range'
201: SETERRQ(PETSC_COMM_SELF,1,' ',ierr)
202: endif
204: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
205: ! Create nonlinear solver context
206: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
207: call SNESCreate(PETSC_COMM_WORLD,mysnes,ierr)
209: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
210: ! Create vector data structures; set function evaluation routine
211: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
213: ! Create distributed array (DMDA) to manage parallel grid and vectors
215: ! This really needs only the star-type stencil, but we use the box
216: ! stencil temporarily.
217: call DMDACreate2d(PETSC_COMM_WORLD, &
218: & DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, &
219: & DMDA_STENCIL_BOX,nfour,nfour,PETSC_DECIDE,PETSC_DECIDE, &
220: & ione,ione,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,user%da,ierr)
221: call DMDAGetInfo(user%da,PETSC_NULL_INTEGER,user%mx,user%my, &
222: & PETSC_NULL_INTEGER, &
223: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
224: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
225: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
226: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
227: & PETSC_NULL_INTEGER,ierr)
229: !
230: ! Visualize the distribution of the array across the processors
231: !
232: ! call DMView(user%da,PETSC_VIEWER_DRAW_WORLD,ierr)
234: ! Extract global and local vectors from DMDA; then duplicate for remaining
235: ! vectors that are the same types
236: call DMCreateGlobalVector(user%da,x,ierr)
237: call VecDuplicate(x,r,ierr)
239: ! Get local grid boundaries (for 2-dimensional DMDA)
240: call DMDAGetCorners(user%da,user%xs,user%ys,PETSC_NULL_INTEGER, &
241: & user%xm,user%ym,PETSC_NULL_INTEGER,ierr)
242: call DMDAGetGhostCorners(user%da,user%gxs,user%gys, &
243: & PETSC_NULL_INTEGER,user%gxm,user%gym, &
244: & PETSC_NULL_INTEGER,ierr)
246: ! Here we shift the starting indices up by one so that we can easily
247: ! use the Fortran convention of 1-based indices (rather 0-based indices).
248: user%xs = user%xs+1
249: user%ys = user%ys+1
250: user%gxs = user%gxs+1
251: user%gys = user%gys+1
253: user%ye = user%ys+user%ym-1
254: user%xe = user%xs+user%xm-1
255: user%gye = user%gys+user%gym-1
256: user%gxe = user%gxs+user%gxm-1
258: call SNESSetApplicationContext(mysnes,user,ierr)
260: ! Set function evaluation routine and vector
261: call SNESSetFunction(mysnes,r,FormFunction,user,ierr)
263: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
264: ! Create matrix data structure; set Jacobian evaluation routine
265: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
267: ! Set Jacobian matrix data structure and default Jacobian evaluation
268: ! routine. User can override with:
269: ! -snes_fd : default finite differencing approximation of Jacobian
270: ! -snes_mf : matrix-free Newton-Krylov method with no preconditioning
271: ! (unless user explicitly sets preconditioner)
272: ! -snes_mf_operator : form preconditioning matrix as set by the user,
273: ! but use matrix-free approx for Jacobian-vector
274: ! products within Newton-Krylov method
275: !
276: ! Note: For the parallel case, vectors and matrices MUST be partitioned
277: ! accordingly. When using distributed arrays (DMDAs) to create vectors,
278: ! the DMDAs determine the problem partitioning. We must explicitly
279: ! specify the local matrix dimensions upon its creation for compatibility
280: ! with the vector distribution. Thus, the generic MatCreate() routine
281: ! is NOT sufficient when working with distributed arrays.
282: !
283: ! Note: Here we only approximately preallocate storage space for the
284: ! Jacobian. See the users manual for a discussion of better techniques
285: ! for preallocating matrix memory.
287: call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-snes_mf', &
288: & matrix_free,ierr)
289: if (.not. matrix_free) then
290: call DMSetMatType(user%da,MATAIJ,ierr)
291: call DMCreateMatrix(user%da,J,ierr)
292: call SNESSetJacobian(mysnes,J,J,FormJacobian,user,ierr)
293: endif
295: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
296: ! Customize nonlinear solver; set runtime options
297: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
298: ! Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
299: call SNESSetFromOptions(mysnes,ierr)
301: ! Test Fortran90 wrapper for SNESSet/Get ApplicationContext()
302: call PetscOptionsGetBool(PETSC_NULL_CHARACTER,'-test_appctx', &
303: & flg,PETSC_NULL_CHARACTER,ierr)
304: if (flg) then
305: call SNESGetApplicationContext(mysnes,puser,ierr)
306: endif
308: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
309: ! Evaluate initial guess; then solve nonlinear system.
310: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
311: ! Note: The user should initialize the vector, x, with the initial guess
312: ! for the nonlinear solver prior to calling SNESSolve(). In particular,
313: ! to employ an initial guess of zero, the user should explicitly set
314: ! this vector to zero by calling VecSet().
316: call FormInitialGuess(mysnes,x,ierr)
317: call SNESSolve(mysnes,PETSC_NULL_OBJECT,x,ierr)
318: call SNESGetIterationNumber(mysnes,its,ierr);
319: if (user%rank .eq. 0) then
320: write(6,100) its
321: endif
322: 100 format('Number of SNES iterations = ',i5)
324: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
325: ! Free work space. All PETSc objects should be destroyed when they
326: ! are no longer needed.
327: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
328: if (.not. matrix_free) call MatDestroy(J,ierr)
329: call VecDestroy(x,ierr)
330: call VecDestroy(r,ierr)
331: call SNESDestroy(mysnes,ierr)
332: call DMDestroy(user%da,ierr)
334: call PetscFinalize(ierr)
335: end
337: ! ---------------------------------------------------------------------
338: !
339: ! FormInitialGuess - Forms initial approximation.
340: !
341: ! Input Parameters:
342: ! X - vector
343: !
344: ! Output Parameter:
345: ! X - vector
346: !
347: ! Notes:
348: ! This routine serves as a wrapper for the lower-level routine
349: ! "InitialGuessLocal", where the actual computations are
350: ! done using the standard Fortran style of treating the local
351: ! vector data as a multidimensional array over the local mesh.
352: ! This routine merely handles ghost point scatters and accesses
353: ! the local vector data via VecGetArrayF90() and VecRestoreArrayF90().
354: !
355: subroutine FormInitialGuess(mysnes,X,ierr)
356: #include <finclude/petscsnesdef.h>
357: use petscsnes
358: use f90module
359: use f90moduleinterfaces
360: ! Input/output variables:
361: type(SNES) mysnes
362: type(userctx), pointer:: puser
363: type(Vec) X
364: PetscErrorCode ierr
366: ! Declarations for use with local arrays:
367: PetscScalar,pointer :: lx_v(:)
368: type(Vec) localX
370: 0
371: call SNESGetApplicationContext(mysnes,puser,ierr)
372: ! Get a pointer to vector data.
373: ! - For default PETSc vectors, VecGetArray90() returns a pointer to
374: ! the data array. Otherwise, the routine is implementation dependent.
375: ! - You MUST call VecRestoreArrayF90() when you no longer need access to
376: ! the array.
377: ! - Note that the interface to VecGetArrayF90() differs from VecGetArray(),
378: ! and is useable from Fortran-90 Only.
380: call DMGetLocalVector(puser%da,localX,ierr)
381: call VecGetArrayF90(localX,lx_v,ierr)
383: ! Compute initial guess over the locally owned part of the grid
384: call InitialGuessLocal(puser,lx_v,ierr)
386: ! Restore vector
387: call VecRestoreArrayF90(localX,lx_v,ierr)
389: ! Insert values into global vector
390: call DMLocalToGlobalBegin(puser%da,localX,INSERT_VALUES,X,ierr)
391: call DMLocalToGlobalEnd(puser%da,localX,INSERT_VALUES,X,ierr)
392: call DMRestoreLocalVector(puser%da,localX,ierr)
394: return
395: end
397: ! ---------------------------------------------------------------------
398: !
399: ! InitialGuessLocal - Computes initial approximation, called by
400: ! the higher level routine FormInitialGuess().
401: !
402: ! Input Parameter:
403: ! x - local vector data
404: !
405: ! Output Parameters:
406: ! x - local vector data
407: ! ierr - error code
408: !
409: ! Notes:
410: ! This routine uses standard Fortran-style computations over a 2-dim array.
411: !
412: subroutine InitialGuessLocal(user,x,ierr)
413: #include <finclude/petscsysdef.h>
414: use petscsys
415: use f90module
416: ! Input/output variables:
417: type (userctx) user
418: PetscScalar x(user%gxs:user%gxe, &
419: & user%gys:user%gye)
420: PetscErrorCode ierr
422: ! Local variables:
423: PetscInt i,j
424: PetscScalar temp1,temp,hx,hy
425: PetscScalar one
427: ! Set parameters
429: 0
430: one = 1.0
431: hx = one/(dble(user%mx-1))
432: hy = one/(dble(user%my-1))
433: temp1 = user%lambda/(user%lambda + one)
435: do 20 j=user%ys,user%ye
436: temp = dble(min(j-1,user%my-j))*hy
437: do 10 i=user%xs,user%xe
438: if (i .eq. 1 .or. j .eq. 1 &
439: & .or. i .eq. user%mx .or. j .eq. user%my) then
440: x(i,j) = 0.0
441: else
442: x(i,j) = temp1 * &
443: & sqrt(min(dble(min(i-1,user%mx-i)*hx),dble(temp)))
444: endif
445: 10 continue
446: 20 continue
448: return
449: end
451: ! ---------------------------------------------------------------------
452: !
453: ! FormFunctionLocal - Computes nonlinear function, called by
454: ! the higher level routine FormFunction().
455: !
456: ! Input Parameter:
457: ! x - local vector data
458: !
459: ! Output Parameters:
460: ! f - local vector data, f(x)
461: ! ierr - error code
462: !
463: ! Notes:
464: ! This routine uses standard Fortran-style computations over a 2-dim array.
465: !
466: subroutine FormFunctionLocal(x,f,user,ierr)
467: #include <finclude/petscsysdef.h>
468: use petscsys
469: use f90module
470: ! Input/output variables:
471: type (userctx) user
472: PetscScalar x(user%gxs:user%gxe, &
473: & user%gys:user%gye)
474: PetscScalar f(user%xs:user%xe, &
475: & user%ys:user%ye)
476: PetscErrorCode ierr
478: ! Local variables:
479: PetscScalar two,one,hx,hy,hxdhy,hydhx,sc
480: PetscScalar u,uxx,uyy
481: PetscInt i,j
483: one = 1.0
484: two = 2.0
485: hx = one/dble(user%mx-1)
486: hy = one/dble(user%my-1)
487: sc = hx*hy*user%lambda
488: hxdhy = hx/hy
489: hydhx = hy/hx
491: ! Compute function over the locally owned part of the grid
493: do 20 j=user%ys,user%ye
494: do 10 i=user%xs,user%xe
495: if (i .eq. 1 .or. j .eq. 1 &
496: & .or. i .eq. user%mx .or. j .eq. user%my) then
497: f(i,j) = x(i,j)
498: else
499: u = x(i,j)
500: uxx = hydhx * (two*u &
501: & - x(i-1,j) - x(i+1,j))
502: uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1))
503: f(i,j) = uxx + uyy - sc*exp(u)
504: endif
505: 10 continue
506: 20 continue
507: 0
508: return
509: end
511: ! ---------------------------------------------------------------------
512: !
513: ! FormJacobian - Evaluates Jacobian matrix.
514: !
515: ! Input Parameters:
516: ! snes - the SNES context
517: ! x - input vector
518: ! dummy - optional user-defined context, as set by SNESSetJacobian()
519: ! (not used here)
520: !
521: ! Output Parameters:
522: ! jac - Jacobian matrix
523: ! jac_prec - optionally different preconditioning matrix (not used here)
524: ! flag - flag indicating matrix structure
525: !
526: ! Notes:
527: ! This routine serves as a wrapper for the lower-level routine
528: ! "FormJacobianLocal", where the actual computations are
529: ! done using the standard Fortran style of treating the local
530: ! vector data as a multidimensional array over the local mesh.
531: ! This routine merely accesses the local vector data via
532: ! VecGetArrayF90() and VecRestoreArrayF90().
533: !
534: ! Notes:
535: ! Due to grid point reordering with DMDAs, we must always work
536: ! with the local grid points, and then transform them to the new
537: ! global numbering with the "ltog" mapping
538: ! We cannot work directly with the global numbers for the original
539: ! uniprocessor grid!
540: !
541: ! Two methods are available for imposing this transformation
542: ! when setting matrix entries:
543: ! (A) MatSetValuesLocal(), using the local ordering (including
544: ! ghost points!)
545: ! - Set matrix entries using the local ordering
546: ! by calling MatSetValuesLocal()
547: ! (B) MatSetValues(), using the global ordering
548: ! - Use DMGetLocalToGlobalMapping() then
549: ! ISLocalToGlobalMappingGetIndicesF90() to extract the local-to-global map
550: ! - Then apply this map explicitly yourself
551: ! - Set matrix entries using the global ordering by calling
552: ! MatSetValues()
553: ! Option (A) seems cleaner/easier in many cases, and is the procedure
554: ! used in this example.
555: !
556: subroutine FormJacobian(mysnes,X,jac,jac_prec,user,ierr)
557: #include <finclude/petscsnesdef.h>
558: use petscsnes
559: use f90module
560: ! Input/output variables:
561: type(SNES) mysnes
562: type(Vec) X
563: type(Mat) jac,jac_prec
564: type(userctx) user
565: PetscErrorCode ierr
567: ! Declarations for use with local arrays:
568: PetscScalar,pointer :: lx_v(:)
569: type(Vec) localX
571: ! Scatter ghost points to local vector, using the 2-step process
572: ! DMGlobalToLocalBegin(), DMGlobalToLocalEnd()
573: ! Computations can be done while messages are in transition,
574: ! by placing code between these two statements.
576: call DMGetLocalVector(user%da,localX,ierr)
577: call DMGlobalToLocalBegin(user%da,X,INSERT_VALUES,localX, &
578: & ierr)
579: call DMGlobalToLocalEnd(user%da,X,INSERT_VALUES,localX,ierr)
581: ! Get a pointer to vector data
582: call VecGetArrayF90(localX,lx_v,ierr)
584: ! Compute entries for the locally owned part of the Jacobian preconditioner.
585: call FormJacobianLocal(lx_v,jac_prec,user,ierr)
587: ! Assemble matrix, using the 2-step process:
588: ! MatAssemblyBegin(), MatAssemblyEnd()
589: ! Computations can be done while messages are in transition,
590: ! by placing code between these two statements.
592: call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
593: ! if (jac .ne. jac_prec) then
594: call MatAssemblyBegin(jac_prec,MAT_FINAL_ASSEMBLY,ierr)
595: ! endif
596: call VecRestoreArrayF90(localX,lx_v,ierr)
597: call DMRestoreLocalVector(user%da,localX,ierr)
598: call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
599: ! if (jac .ne. jac_prec) then
600: call MatAssemblyEnd(jac_prec,MAT_FINAL_ASSEMBLY,ierr)
601: ! endif
603: ! Tell the matrix we will never add a new nonzero location to the
604: ! matrix. If we do it will generate an error.
606: call MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE, &
607: & ierr)
609: return
610: end
612: ! ---------------------------------------------------------------------
613: !
614: ! FormJacobianLocal - Computes Jacobian preconditioner matrix,
615: ! called by the higher level routine FormJacobian().
616: !
617: ! Input Parameters:
618: ! x - local vector data
619: !
620: ! Output Parameters:
621: ! jac_prec - Jacobian preconditioner matrix
622: ! ierr - error code
623: !
624: ! Notes:
625: ! This routine uses standard Fortran-style computations over a 2-dim array.
626: !
627: ! Notes:
628: ! Due to grid point reordering with DMDAs, we must always work
629: ! with the local grid points, and then transform them to the new
630: ! global numbering with the "ltog" mapping
631: ! We cannot work directly with the global numbers for the original
632: ! uniprocessor grid!
633: !
634: ! Two methods are available for imposing this transformation
635: ! when setting matrix entries:
636: ! (A) MatSetValuesLocal(), using the local ordering (including
637: ! ghost points!)
638: ! - Set matrix entries using the local ordering
639: ! by calling MatSetValuesLocal()
640: ! (B) MatSetValues(), using the global ordering
641: ! - Set matrix entries using the global ordering by calling
642: ! MatSetValues()
643: ! Option (A) seems cleaner/easier in many cases, and is the procedure
644: ! used in this example.
645: !
646: subroutine FormJacobianLocal(x,jac_prec,user,ierr)
647: #include <finclude/petscmatdef.h>
648: use petscmat
649: use f90module
650: ! Input/output variables:
651: type (userctx) user
652: PetscScalar x(user%gxs:user%gxe, &
653: & user%gys:user%gye)
654: type(Mat) jac_prec
655: PetscErrorCode ierr
657: ! Local variables:
658: PetscInt row,col(5),i,j
659: PetscInt ione,ifive
660: PetscScalar two,one,hx,hy,hxdhy
661: PetscScalar hydhx,sc,v(5)
663: ! Set parameters
664: ione = 1
665: ifive = 5
666: one = 1.0
667: two = 2.0
668: hx = one/dble(user%mx-1)
669: hy = one/dble(user%my-1)
670: sc = hx*hy
671: hxdhy = hx/hy
672: hydhx = hy/hx
674: ! Compute entries for the locally owned part of the Jacobian.
675: ! - Currently, all PETSc parallel matrix formats are partitioned by
676: ! contiguous chunks of rows across the processors.
677: ! - Each processor needs to insert only elements that it owns
678: ! locally (but any non-local elements will be sent to the
679: ! appropriate processor during matrix assembly).
680: ! - Here, we set all entries for a particular row at once.
681: ! - We can set matrix entries either using either
682: ! MatSetValuesLocal() or MatSetValues(), as discussed above.
683: ! - Note that MatSetValues() uses 0-based row and column numbers
684: ! in Fortran as well as in C.
686: do 20 j=user%ys,user%ye
687: row = (j - user%gys)*user%gxm + user%xs - user%gxs - 1
688: do 10 i=user%xs,user%xe
689: row = row + 1
690: ! boundary points
691: if (i .eq. 1 .or. j .eq. 1 &
692: & .or. i .eq. user%mx .or. j .eq. user%my) then
693: col(1) = row
694: v(1) = one
695: call MatSetValuesLocal(jac_prec,ione,row,ione,col,v, &
696: & INSERT_VALUES,ierr)
697: ! interior grid points
698: else
699: v(1) = -hxdhy
700: v(2) = -hydhx
701: v(3) = two*(hydhx + hxdhy) &
702: & - sc*user%lambda*exp(x(i,j))
703: v(4) = -hydhx
704: v(5) = -hxdhy
705: col(1) = row - user%gxm
706: col(2) = row - 1
707: col(3) = row
708: col(4) = row + 1
709: col(5) = row + user%gxm
710: call MatSetValuesLocal(jac_prec,ione,row,ifive,col,v, &
711: & INSERT_VALUES,ierr)
712: endif
713: 10 continue
714: 20 continue
715: return
716: end