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