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: ! Modified from ex5f90.F by Mike McCourt <mccomic@iit.edu>
10: ! for testing Fortran interface on
11: ! SNESLineSearchSet(), SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck()
12: !
13: ! --------------------------------------------------------------------------
14: !
15: ! Solid Fuel Ignition (SFI) problem. This problem is modeled by
16: ! the partial differential equation
17: !
18: ! -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1,
19: !
20: ! with boundary conditions
21: !
22: ! u = 0 for x = 0, x = 1, y = 0, y = 1.
23: !
24: ! A finite difference approximation with the usual 5-point stencil
25: ! is used to discretize the boundary value problem to obtain a nonlinear
26: ! system of equations.
27: !
28: ! The uniprocessor version of this code is snes/examples/tutorials/ex4f.F
29: !
30: ! --------------------------------------------------------------------------
31: ! The following define must be used before including any PETSc include files
32: ! into a module or interface. This is because they can't handle declarations
33: ! in them
34: !
36: module f90module
37: type userctx
38: #include <finclude/petscsysdef.h>
39: #include <finclude/petscvecdef.h>
40: #include <finclude/petscdmdef.h>
41: DM da
42: integer xs,xe,xm,gxs,gxe,gxm
43: integer ys,ye,ym,gys,gye,gym
44: integer mx,my,rank
45: double precision lambda
46: end type userctx
47: contains
48: ! ---------------------------------------------------------------------
49: !
50: ! FormFunction - Evaluates nonlinear function, F(x).
51: !
52: ! Input Parameters:
53: ! snes - the SNES context
54: ! X - input vector
55: ! dummy - optional user-defined context, as set by SNESSetFunction()
56: ! (not used here)
57: !
58: ! Output Parameter:
59: ! F - function vector
60: !
61: ! Notes:
62: ! This routine serves as a wrapper for the lower-level routine
63: ! "FormFunctionLocal", where the actual computations are
64: ! done using the standard Fortran style of treating the local
65: ! vector data as a multidimensional array over the local mesh.
66: ! This routine merely handles ghost point scatters and accesses
67: ! the local vector data via VecGetArrayF90() and VecRestoreArrayF90().
68: !
69: subroutine FormFunction(snes,X,F,user,ierr)
70: implicit none
72: #include <finclude/petscsys.h>
73: #include <finclude/petscvec.h>
74: #include <finclude/petscdmda.h>
75: #include <finclude/petscis.h>
76: #include <finclude/petscmat.h>
77: #include <finclude/petscksp.h>
78: #include <finclude/petscpc.h>
79: #include <finclude/petscsnes.h>
81: #include <finclude/petscvec.h90>
84: ! Input/output variables:
85: SNES snes
86: Vec X,F
87: integer ierr
88: type (userctx) user
90: ! Declarations for use with local arrays:
91: PetscScalar,pointer :: lx_v(:),lf_v(:)
92: Vec localX
94: write(*,*)"Inside FormFunction, user%xm=",user%xm
96: ! Scatter ghost points to local vector, using the 2-step process
97: ! DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
98: ! By placing code between these two statements, computations can
99: ! be done while messages are in transition.
101: call DMGetLocalVector(user%da,localX,ierr)
102: call DMGlobalToLocalBegin(user%da,X,INSERT_VALUES, &
103: & localX,ierr)
104: call DMGlobalToLocalEnd(user%da,X,INSERT_VALUES,localX,ierr)
106: ! Get a pointer to vector data.
107: ! - For default PETSc vectors, VecGetArray90() returns a pointer to
108: ! the data array. Otherwise, the routine is implementation dependent.
109: ! - You MUST call VecRestoreArrayF90() when you no longer need access to
110: ! the array.
111: ! - Note that the interface to VecGetArrayF90() differs from VecGetArray(),
112: ! and is useable from Fortran-90 Only.
114: call VecGetArrayF90(localX,lx_v,ierr)
115: call VecGetArrayF90(F,lf_v,ierr)
117: ! Compute function over the locally owned part of the grid
119: call FormFunctionLocal(lx_v,lf_v,user,ierr)
121: ! Restore vectors
123: call VecRestoreArrayF90(localX,lx_v,ierr)
124: call VecRestoreArrayF90(F,lf_v,ierr)
126: ! Insert values into global vector
128: call DMRestoreLocalVector(user%da,localX,ierr)
129: call PetscLogFlops(11.0d0*user%ym*user%xm,ierr)
131: ! call VecView(X,PETSC_VIEWER_STDOUT_WORLD,ierr)
132: ! call VecView(F,PETSC_VIEWER_STDOUT_WORLD,ierr)
134: return
135: end subroutine formfunction
136: end module f90module
140: program main
141: use f90module
142: implicit none
143: !
144: !
145: #include <finclude/petscsys.h>
146: #include <finclude/petscvec.h>
147: #include <finclude/petscdmda.h>
148: #include <finclude/petscis.h>
149: #include <finclude/petscmat.h>
150: #include <finclude/petscksp.h>
151: #include <finclude/petscpc.h>
152: #include <finclude/petscsnes.h>
153: #include <finclude/petscvec.h90>
154: #include <finclude/petscdmda.h90>
156: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157: ! Variable declarations
158: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
159: !
160: ! Variables:
161: ! snes - nonlinear solver
162: ! x, r - solution, residual vectors
163: ! J - Jacobian matrix
164: ! its - iterations for convergence
165: ! Nx, Ny - number of preocessors in x- and y- directions
166: ! matrix_free - flag - 1 indicates matrix-free version
167: !
168: !
169: SNES snes
170: Vec x,r
171: Mat J
172: integer its,matrix_free,flg,ierr
173: double precision lambda_max,lambda_min
174: type (userctx) user
175: PetscBool test_linesearch,test_check
177: ! Note: Any user-defined Fortran routines (such as FormJacobian)
178: ! MUST be declared as external.
180: external FormInitialGuess,FormJacobian,FormLineSearch
181: external FormPostCheck,FormPreCheck
183: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
184: ! Initialize program
185: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
188: call MPI_Comm_rank(PETSC_COMM_WORLD,user%rank,ierr)
190: ! Initialize problem parameters
192: lambda_max = 6.81
193: lambda_min = 0.0
194: user%lambda = 6.0
195: call PetscOptionsGetReal(PETSC_NULL_CHARACTER,'-par', &
196: & user%lambda,flg,ierr)
197: if (user%lambda .ge. lambda_max .or. user%lambda .le. lambda_min) &
198: & then
199: if (user%rank .eq. 0) write(6,*) 'Lambda is out of range'
200: SETERRQ(PETSC_COMM_SELF,1,' ',ierr)
201: endif
204: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
205: ! Create nonlinear solver context
206: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
208: call SNESCreate(PETSC_COMM_WORLD,snes,ierr)
210: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
211: ! Create vector data structures; set function evaluation routine
212: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
214: ! Create distributed array (DMDA) to manage parallel grid and vectors
216: ! This really needs only the star-type stencil, but we use the box
217: ! stencil temporarily.
218: call DMDACreate2d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE, &
219: & DMDA_BOUNDARY_NONE,DMDA_STENCIL_BOX, &
220: & -4,-4,PETSC_DECIDE,PETSC_DECIDE,1,1, &
221: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,user%da,ierr)
222: call DMDAGetInfo(user%da,PETSC_NULL_INTEGER,user%mx,user%my, &
223: & 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,PETSC_NULL_INTEGER, &
228: & PETSC_NULL_INTEGER,ierr)
229: 230: !
231: ! Visualize the distribution of the array across the processors
232: !
233: ! call DMView(user%da,PETSC_VIEWER_DRAW_WORLD,ierr)
235: ! Extract global and local vectors from DMDA; then duplicate for remaining
236: ! vectors that are the same types
238: call DMCreateGlobalVector(user%da,x,ierr)
239: call VecDuplicate(x,r,ierr)
241: ! Get local grid boundaries (for 2-dimensional DMDA)
243: call DMDAGetCorners(user%da,user%xs,user%ys,PETSC_NULL_INTEGER, &
244: & user%xm,user%ym,PETSC_NULL_INTEGER,ierr)
245: call DMDAGetGhostCorners(user%da,user%gxs,user%gys, &
246: & PETSC_NULL_INTEGER,user%gxm,user%gym, &
247: & PETSC_NULL_INTEGER,ierr)
249: ! Here we shift the starting indices up by one so that we can easily
250: ! use the Fortran convention of 1-based indices (rather 0-based indices).
252: user%xs = user%xs+1
253: user%ys = user%ys+1
254: user%gxs = user%gxs+1
255: user%gys = user%gys+1
257: user%ye = user%ys+user%ym-1
258: user%xe = user%xs+user%xm-1
259: user%gye = user%gys+user%gym-1
260: user%gxe = user%gxs+user%gxm-1
262: ! Set function evaluation routine and vector
264: call SNESSetFunction(snes,r,FormFunction,user,ierr)
266: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
267: ! Create matrix data structure; set Jacobian evaluation routine
268: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
270: ! Set Jacobian matrix data structure and default Jacobian evaluation
271: ! routine. User can override with:
272: ! -snes_fd : default finite differencing approximation of Jacobian
273: ! -snes_mf : matrix-free Newton-Krylov method with no preconditioning
274: ! (unless user explicitly sets preconditioner)
275: ! -snes_mf_operator : form preconditioning matrix as set by the user,
276: ! but use matrix-free approx for Jacobian-vector
277: ! products within Newton-Krylov method
278: !
279: ! Note: For the parallel case, vectors and matrices MUST be partitioned
280: ! accordingly. When using distributed arrays (DMDAs) to create vectors,
281: ! the DMDAs determine the problem partitioning. We must explicitly
282: ! specify the local matrix dimensions upon its creation for compatibility
283: ! with the vector distribution. Thus, the generic MatCreate() routine
284: ! is NOT sufficient when working with distributed arrays.
285: !
286: ! Note: Here we only approximately preallocate storage space for the
287: ! Jacobian. See the users manual for a discussion of better techniques
288: ! for preallocating matrix memory.
290: call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-snes_mf', &
291: & matrix_free,ierr)
292: if (matrix_free .eq. 0) then
293: call DMCreateMatrix(user%da,MATAIJ,J,ierr)
294: call SNESSetJacobian(snes,J,J,FormJacobian,user,ierr)
295: endif
297: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
298: ! Customize nonlinear solver; set runtime options
299: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
301: ! Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
303: call SNESSetFromOptions(snes,ierr)
304: test_linesearch = 0
305: test_check = 0
306: call PetscOptionsGetBool(PETSC_NULL_CHARACTER,'-test_check',
307: & test_check,PETSC_NULL_INTEGER,ierr)
308: if (test_check.eq.1) then
309: call SNESLineSearchSetPreCheck(snes,FormPreCheck,user,ierr)
310: call SNESLineSearchSetPostCheck(snes,FormPostCheck,user,ierr)
311: else
312: call PetscOptionsGetBool(PETSC_NULL_CHARACTER,
313: & '-test_linesearch',test_linesearch,PETSC_NULL_INTEGER,ierr)
314: if (test_linesearch.eq.1) then
315: call SNESLineSearchSet(snes,FormLineSearch,user,ierr)
316: end if
317: end if
318: 319: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
320: ! Evaluate initial guess; then solve nonlinear system.
321: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
323: ! Note: The user should initialize the vector, x, with the initial guess
324: ! for the nonlinear solver prior to calling SNESSolve(). In particular,
325: ! to employ an initial guess of zero, the user should explicitly set
326: ! this vector to zero by calling VecSet().
328: call FormInitialGuess(user,x,ierr)
329: write(*,*)"Before SNESSolve"
330: write(*,*)"user%xm=",user%xm
331: call SNESSolve(snes,PETSC_NULL_OBJECT,x,ierr)
332: call SNESGetIterationNumber(snes,its,ierr);
333: if (user%rank .eq. 0) then
334: write(6,100) its
335: endif
336: 100 format('Number of SNES iterations = ',i5)
338: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
339: ! Free work space. All PETSc objects should be destroyed when they
340: ! are no longer needed.
341: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
343: if (matrix_free .eq. 0) call MatDestroy(J,ierr)
344: call VecDestroy(x,ierr)
345: call VecDestroy(r,ierr)
346: call SNESDestroy(snes,ierr)
347: call DMDestroy(user%da,ierr)
348: call PetscFinalize(ierr)
349: end
351: ! ---------------------------------------------------------------------
352: !
353: ! FormLineSearch - Applies the line search to the step size
354: !
355: subroutine FormLineSearch(snes,user,x,f,g,y,w,fnorm,ynorm,gnorm,
356: & flag,ierr)
358: use f90module
360: #include <finclude/petscsys.h>
361: #include <finclude/petscvec.h>
362: #include <finclude/petscvec.h90>
363: #include <finclude/petscmat.h>
364: #include <finclude/petscmat.h90>
365: #include <finclude/petscksp.h>
366: #include <finclude/petscpc.h>
367: #include <finclude/petscsnes.h>
369: SNES snes
370: type (userctx) user
371: Vec x,f,g,y,w
372: PetscReal fnorm,ynorm,gnorm
373: PetscBool flag
374: PetscErrorCode ierr
376: PetscScalar mone
378: write(*,*)"Inside FormLineSearch, user%xm=",user%xm
379: mone = -1.0d0
380: flag = 0
381: call VecNorm(y,NORM_2,ynorm,ierr)
382: call VecAYPX(y,mone,x,ierr)
383: call SNESComputeFunction(snes,y,g,ierr)
384: call VecNorm(g,NORM_2,gnorm,ierr)
385: return
386: end
388: ! ---------------------------------------------------------------------
389: !
390: ! FormInitialGuess - Forms initial approximation.
391: !
392: ! Input Parameters:
393: ! X - vector
394: !
395: ! Output Parameter:
396: ! X - vector
397: !
398: ! Notes:
399: ! This routine serves as a wrapper for the lower-level routine
400: ! "InitialGuessLocal", where the actual computations are
401: ! done using the standard Fortran style of treating the local
402: ! vector data as a multidimensional array over the local mesh.
403: ! This routine merely handles ghost point scatters and accesses
404: ! the local vector data via VecGetArrayF90() and VecRestoreArrayF90().
405: !
406: subroutine FormInitialGuess(user,X,ierr)
407: use f90module
408: implicit none
410: #include <finclude/petscvec.h90>
411: #include <finclude/petscsys.h>
412: #include <finclude/petscvec.h>
413: #include <finclude/petscdmda.h>
414: #include <finclude/petscis.h>
415: #include <finclude/petscmat.h>
416: #include <finclude/petscksp.h>
417: #include <finclude/petscpc.h>
418: #include <finclude/petscsnes.h>
420: ! Input/output variables:
421: type (userctx) user
422: Vec X
423: integer ierr
424: 425: ! Declarations for use with local arrays:
426: PetscScalar,pointer :: lx_v(:)
427: Vec localX
429: 0
431: ! Get a pointer to vector data.
432: ! - For default PETSc vectors, VecGetArray90() returns a pointer to
433: ! the data array. Otherwise, the routine is implementation dependent.
434: ! - You MUST call VecRestoreArrayF90() when you no longer need access to
435: ! the array.
436: ! - Note that the interface to VecGetArrayF90() differs from VecGetArray(),
437: ! and is useable from Fortran-90 Only.
439: call DMGetLocalVector(user%da,localX,ierr)
440: call VecGetArrayF90(localX,lx_v,ierr)
442: ! Compute initial guess over the locally owned part of the grid
444: call InitialGuessLocal(user,lx_v,ierr)
446: ! Restore vector
448: call VecRestoreArrayF90(localX,lx_v,ierr)
450: ! Insert values into global vector
452: call DMLocalToGlobalBegin(user%da,localX,INSERT_VALUES,X,ierr)
453: call DMLocalToGlobalEnd(user%da,localX,INSERT_VALUES,X,ierr)
454: call DMRestoreLocalVector(user%da,localX,ierr)
456: return
457: end
459: ! ---------------------------------------------------------------------
460: !
461: ! InitialGuessLocal - Computes initial approximation, called by
462: ! the higher level routine FormInitialGuess().
463: !
464: ! Input Parameter:
465: ! x - local vector data
466: !
467: ! Output Parameters:
468: ! x - local vector data
469: ! ierr - error code
470: !
471: ! Notes:
472: ! This routine uses standard Fortran-style computations over a 2-dim array.
473: !
474: subroutine InitialGuessLocal(user,x,ierr)
475: use f90module
476: implicit none
478: #include <finclude/petscsys.h>
479: #include <finclude/petscvec.h>
480: #include <finclude/petscdmda.h>
481: #include <finclude/petscis.h>
482: #include <finclude/petscmat.h>
483: #include <finclude/petscksp.h>
484: #include <finclude/petscpc.h>
485: #include <finclude/petscsnes.h>
487: ! Input/output variables:
488: type (userctx) user
489: PetscScalar x(user%gxs:user%gxe, &
490: & user%gys:user%gye)
491: integer ierr
493: ! Local variables:
494: integer i,j
495: PetscScalar temp1,temp,hx,hy
496: PetscScalar one
498: ! Set parameters
500: 0
501: one = 1.0
502: hx = one/(dble(user%mx-1))
503: hy = one/(dble(user%my-1))
504: temp1 = user%lambda/(user%lambda + one)
506: do 20 j=user%ys,user%ye
507: temp = dble(min(j-1,user%my-j))*hy
508: do 10 i=user%xs,user%xe
509: if (i .eq. 1 .or. j .eq. 1 &
510: & .or. i .eq. user%mx .or. j .eq. user%my) then
511: x(i,j) = 0.0
512: else
513: x(i,j) = temp1 * &
514: & sqrt(min(dble(min(i-1,user%mx-i)*hx),dble(temp)))
515: endif
516: 10 continue
517: 20 continue
519: return
520: end
522: ! ---------------------------------------------------------------------
523: !
524: ! FormFunctionLocal - Computes nonlinear function, called by
525: ! the higher level routine FormFunction().
526: !
527: ! Input Parameter:
528: ! x - local vector data
529: !
530: ! Output Parameters:
531: ! f - local vector data, f(x)
532: ! ierr - error code
533: !
534: ! Notes:
535: ! This routine uses standard Fortran-style computations over a 2-dim array.
536: !
537: subroutine FormFunctionLocal(x,f,user,ierr)
538: use f90module
540: implicit none
542: ! Input/output variables:
543: type (userctx) user
544: PetscScalar x(user%gxs:user%gxe, &
545: & user%gys:user%gye)
546: PetscScalar f(user%xs:user%xe, &
547: & user%ys:user%ye)
548: integer ierr
550: ! Local variables:
551: PetscScalar two,one,hx,hy,hxdhy,hydhx,sc
552: PetscScalar u,uxx,uyy
553: integer i,j
555: one = 1.0
556: two = 2.0
557: hx = one/dble(user%mx-1)
558: hy = one/dble(user%my-1)
559: sc = hx*hy*user%lambda
560: hxdhy = hx/hy
561: hydhx = hy/hx
563: ! Compute function over the locally owned part of the grid
565: do 20 j=user%ys,user%ye
566: do 10 i=user%xs,user%xe
567: if (i .eq. 1 .or. j .eq. 1 &
568: & .or. i .eq. user%mx .or. j .eq. user%my) then
569: f(i,j) = x(i,j)
570: else
571: u = x(i,j)
572: uxx = hydhx * (two*u &
573: & - x(i-1,j) - x(i+1,j))
574: uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1))
575: f(i,j) = uxx + uyy - sc*exp(u)
576: endif
577: 10 continue
578: 20 continue
580: return
581: end
583: ! ---------------------------------------------------------------------
584: !
585: ! FormJacobian - Evaluates Jacobian matrix.
586: !
587: ! Input Parameters:
588: ! snes - the SNES context
589: ! x - input vector
590: ! dummy - optional user-defined context, as set by SNESSetJacobian()
591: ! (not used here)
592: !
593: ! Output Parameters:
594: ! jac - Jacobian matrix
595: ! jac_prec - optionally different preconditioning matrix (not used here)
596: ! flag - flag indicating matrix structure
597: !
598: ! Notes:
599: ! This routine serves as a wrapper for the lower-level routine
600: ! "FormJacobianLocal", where the actual computations are
601: ! done using the standard Fortran style of treating the local
602: ! vector data as a multidimensional array over the local mesh.
603: ! This routine merely accesses the local vector data via
604: ! VecGetArrayF90() and VecRestoreArrayF90().
605: !
606: ! Notes:
607: ! Due to grid point reordering with DMDAs, we must always work
608: ! with the local grid points, and then transform them to the new
609: ! global numbering with the "ltog" mapping (via DMDAGetGlobalIndicesF90()).
610: ! We cannot work directly with the global numbers for the original
611: ! uniprocessor grid!
612: !
613: ! Two methods are available for imposing this transformation
614: ! when setting matrix entries:
615: ! (A) MatSetValuesLocal(), using the local ordering (including
616: ! ghost points!)
617: ! - Use DMDAGetGlobalIndicesF90() to extract the local-to-global map
618: ! - Associate this map with the matrix by calling
619: ! MatSetLocalToGlobalMapping() once
620: ! - Set matrix entries using the local ordering
621: ! by calling MatSetValuesLocal()
622: ! (B) MatSetValues(), using the global ordering
623: ! - Use DMDAGetGlobalIndicesF90() to extract the local-to-global map
624: ! - Then apply this map explicitly yourself
625: ! - Set matrix entries using the global ordering by calling
626: ! MatSetValues()
627: ! Option (A) seems cleaner/easier in many cases, and is the procedure
628: ! used in this example.
629: !
630: subroutine FormJacobian(snes,X,jac,jac_prec,flag,user,ierr)
631: use f90module
632: implicit none
634: #include <finclude/petscsys.h>
635: #include <finclude/petscvec.h>
636: #include <finclude/petscdmda.h>
637: #include <finclude/petscis.h>
638: #include <finclude/petscmat.h>
639: #include <finclude/petscksp.h>
640: #include <finclude/petscpc.h>
641: #include <finclude/petscsnes.h>
643: #include <finclude/petscvec.h90>
645: ! Input/output variables:
646: SNES snes
647: Vec X
648: Mat jac,jac_prec
649: MatStructure flag
650: type(userctx) user
651: integer ierr
653: ! Declarations for use with local arrays:
654: PetscScalar,pointer :: lx_v(:)
655: Vec localX
657: ! Scatter ghost points to local vector, using the 2-step process
658: ! DMGlobalToLocalBegin(), DMGlobalToLocalEnd()
659: ! Computations can be done while messages are in transition,
660: ! by placing code between these two statements.
662: call DMGetLocalVector(user%da,localX,ierr)
663: call DMGlobalToLocalBegin(user%da,X,INSERT_VALUES,localX, &
664: & ierr)
665: call DMGlobalToLocalEnd(user%da,X,INSERT_VALUES,localX,ierr)
667: ! Get a pointer to vector data
669: call VecGetArrayF90(localX,lx_v,ierr)
671: ! Compute entries for the locally owned part of the Jacobian.
673: call FormJacobianLocal(lx_v,jac,jac_prec,user,ierr)
675: ! Assemble matrix, using the 2-step process:
676: ! MatAssemblyBegin(), MatAssemblyEnd()
677: ! Computations can be done while messages are in transition,
678: ! by placing code between these two statements.
680: call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
681: call VecRestoreArrayF90(localX,lx_v,ierr)
682: call DMRestoreLocalVector(user%da,localX,ierr)
683: call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
685: ! Set flag to indicate that the Jacobian matrix retains an identical
686: ! nonzero structure throughout all nonlinear iterations (although the
687: ! values of the entries change). Thus, we can save some work in setting
688: ! up the preconditioner (e.g., no need to redo symbolic factorization for
689: ! ILU/ICC preconditioners).
690: ! - If the nonzero structure of the matrix is different during
691: ! successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
692: ! must be used instead. If you are unsure whether the matrix
693: ! structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
694: ! - Caution: If you specify SAME_NONZERO_PATTERN, PETSc
695: ! believes your assertion and does not check the structure
696: ! of the matrix. If you erroneously claim that the structure
697: ! is the same when it actually is not, the new preconditioner
698: ! will not function correctly. Thus, use this optimization
699: ! feature with caution!
701: flag = SAME_NONZERO_PATTERN
703: ! Tell the matrix we will never add a new nonzero location to the
704: ! matrix. If we do it will generate an error.
706: ! call MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,ierr)
708: return
709: end
711: ! ---------------------------------------------------------------------
712: !
713: ! FormJacobianLocal - Computes Jacobian matrix, called by
714: ! the higher level routine FormJacobian().
715: !
716: ! Input Parameters:
717: ! x - local vector data
718: !
719: ! Output Parameters:
720: ! jac - Jacobian matrix
721: ! jac_prec - optionally different preconditioning matrix (not used here)
722: ! ierr - error code
723: !
724: ! Notes:
725: ! This routine uses standard Fortran-style computations over a 2-dim array.
726: !
727: ! Notes:
728: ! Due to grid point reordering with DMDAs, we must always work
729: ! with the local grid points, and then transform them to the new
730: ! global numbering with the "ltog" mapping (via DMDAGetGlobalIndicesF90()).
731: ! We cannot work directly with the global numbers for the original
732: ! uniprocessor grid!
733: !
734: ! Two methods are available for imposing this transformation
735: ! when setting matrix entries:
736: ! (A) MatSetValuesLocal(), using the local ordering (including
737: ! ghost points!)
738: ! - Use DMDAGetGlobalIndicesF90() to extract the local-to-global map
739: ! - Associate this map with the matrix by calling
740: ! MatSetLocalToGlobalMapping() once
741: ! - Set matrix entries using the local ordering
742: ! by calling MatSetValuesLocal()
743: ! (B) MatSetValues(), using the global ordering
744: ! - Use DMDAGetGlobalIndicesF90() to extract the local-to-global map
745: ! - Then apply this map explicitly yourself
746: ! - Set matrix entries using the global ordering by calling
747: ! MatSetValues()
748: ! Option (A) seems cleaner/easier in many cases, and is the procedure
749: ! used in this example.
750: !
751: subroutine FormJacobianLocal(x,jac,jac_prec,user,ierr)
752: use f90module
753: implicit none
755: #include <finclude/petscsys.h>
756: #include <finclude/petscvec.h>
757: #include <finclude/petscdmda.h>
758: #include <finclude/petscis.h>
759: #include <finclude/petscmat.h>
760: #include <finclude/petscksp.h>
761: #include <finclude/petscpc.h>
762: #include <finclude/petscsnes.h>
764: ! Input/output variables:
765: type (userctx) user
766: PetscScalar x(user%gxs:user%gxe, &
767: & user%gys:user%gye)
768: Mat jac,jac_prec
769: integer ierr
771: ! Local variables:
772: integer row,col(5),i,j
773: PetscScalar two,one,hx,hy,hxdhy
774: PetscScalar hydhx,sc,v(5)
776: ! Set parameters
778: one = 1.0
779: two = 2.0
780: hx = one/dble(user%mx-1)
781: hy = one/dble(user%my-1)
782: sc = hx*hy
783: hxdhy = hx/hy
784: hydhx = hy/hx
786: ! Compute entries for the locally owned part of the Jacobian.
787: ! - Currently, all PETSc parallel matrix formats are partitioned by
788: ! contiguous chunks of rows across the processors.
789: ! - Each processor needs to insert only elements that it owns
790: ! locally (but any non-local elements will be sent to the
791: ! appropriate processor during matrix assembly).
792: ! - Here, we set all entries for a particular row at once.
793: ! - We can set matrix entries either using either
794: ! MatSetValuesLocal() or MatSetValues(), as discussed above.
795: ! - Note that MatSetValues() uses 0-based row and column numbers
796: ! in Fortran as well as in C.
798: do 20 j=user%ys,user%ye
799: row = (j - user%gys)*user%gxm + user%xs - user%gxs - 1
800: do 10 i=user%xs,user%xe
801: row = row + 1
802: ! boundary points
803: if (i .eq. 1 .or. j .eq. 1 &
804: & .or. i .eq. user%mx .or. j .eq. user%my) then
805: col(1) = row
806: v(1) = one
807: call MatSetValuesLocal(jac,1,row,1,col,v, &
808: & INSERT_VALUES,ierr)
809: ! interior grid points
810: else
811: v(1) = -hxdhy
812: v(2) = -hydhx
813: v(3) = two*(hydhx + hxdhy) &
814: & - sc*user%lambda*exp(x(i,j))
815: v(4) = -hydhx
816: v(5) = -hxdhy
817: col(1) = row - user%gxm
818: col(2) = row - 1
819: col(3) = row
820: col(4) = row + 1
821: col(5) = row + user%gxm
822: call MatSetValuesLocal(jac,1,row,5,col,v, &
823: & INSERT_VALUES,ierr)
824: endif
825: 10 continue
826: 20 continue
828: return
829: end
831: subroutine FormPreCheck(snes,X,Y,user,changed_Y,ierr)
832: use f90module
834: SNES snes
835: Vec X,Y
836: type (userctx) user
837: PetscBool changed_Y
838: PetscErrorCode ierr
840: write(*,*)"Inside formPreCheck, user%xm=",user%xm
841: end subroutine formPreCheck
843: subroutine FormPostCheck(snes,X,Y,W,user,changed_Y,changed_W,ierr)
844: use f90module
846: SNES snes
847: Vec X,Y,W
848: type (userctx) user
849: PetscBool changed_Y,changed_W
850: PetscErrorCode ierr
852: write(*,*)"Inside formPostCheck, user%xm=",user%xm
853: end subroutine formPostCheck