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