xref: /petsc/src/snes/tests/ex1f.F90 (revision e7a95102f46630f317be643b805dc1c3f4655aeb)
1c4762a1bSJed Brown!
2c4762a1bSJed Brown!  Description: This example solves a nonlinear system on 1 processor with SNES.
3c4762a1bSJed Brown!  We solve the  Bratu (SFI - solid fuel ignition) problem in a 2D rectangular
4c4762a1bSJed Brown!  domain.  The command line options include:
5c4762a1bSJed Brown!    -par <parameter>, where <parameter> indicates the nonlinearity of the problem
6c4762a1bSJed Brown!       problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)
7c4762a1bSJed Brown!    -mx <xg>, where <xg> = number of grid points in the x-direction
8c4762a1bSJed Brown!    -my <yg>, where <yg> = number of grid points in the y-direction
9c4762a1bSJed Brown!
10c4762a1bSJed Brown
11c4762a1bSJed Brown!
12c4762a1bSJed Brown!  --------------------------------------------------------------------------
13c4762a1bSJed Brown!
14c4762a1bSJed Brown!  Solid Fuel Ignition (SFI) problem.  This problem is modeled by
15c4762a1bSJed Brown!  the partial differential equation
16c4762a1bSJed Brown!
17c4762a1bSJed Brown!          -Laplacian u - lambda*exp(u) = 0,  0 < x,y < 1,
18c4762a1bSJed Brown!
19c4762a1bSJed Brown!  with boundary conditions
20c4762a1bSJed Brown!
21c4762a1bSJed Brown!           u = 0  for  x = 0, x = 1, y = 0, y = 1.
22c4762a1bSJed Brown!
23c4762a1bSJed Brown!  A finite difference approximation with the usual 5-point stencil
24c4762a1bSJed Brown!  is used to discretize the boundary value problem to obtain a nonlinear
25c4762a1bSJed Brown!  system of equations.
26c4762a1bSJed Brown!
27c4762a1bSJed Brown!  The parallel version of this code is snes/tutorials/ex5f.F
28c4762a1bSJed Brown!
29c4762a1bSJed Brown!  --------------------------------------------------------------------------
30c4762a1bSJed Brown#include <petsc/finclude/petscsnes.h>
31c5e229c2SMartin Diehl#include <petsc/finclude/petscdraw.h>
32*e7a95102SMartin Diehlmodule ex1f_mod
33c4762a1bSJed Brown  use petscsnes
34c4762a1bSJed Brown  implicit none
35*e7a95102SMartin Diehlcontains
36*e7a95102SMartin Diehl  subroutine postcheck(snes, x, y, w, changed_y, changed_w, ctx, ierr)
37c4762a1bSJed Brown    SNES snes
38c4762a1bSJed Brown    PetscReal norm
39c4762a1bSJed Brown    Vec tmp, x, y, w
40c4762a1bSJed Brown    PetscBool changed_w, changed_y
41c4762a1bSJed Brown    PetscErrorCode ierr
42c4762a1bSJed Brown    PetscInt ctx
43c4762a1bSJed Brown    PetscScalar mone
4424fb275aSStefano Zampini    MPI_Comm comm
45c4762a1bSJed Brown
4624fb275aSStefano Zampini    character(len=PETSC_MAX_PATH_LEN) :: outputString
4724fb275aSStefano Zampini
4824fb275aSStefano Zampini    PetscCallA(PetscObjectGetComm(snes, comm, ierr))
49d8606c27SBarry Smith    PetscCallA(VecDuplicate(x, tmp, ierr))
50c4762a1bSJed Brown    mone = -1.0
51d8606c27SBarry Smith    PetscCallA(VecWAXPY(tmp, mone, x, w, ierr))
52d8606c27SBarry Smith    PetscCallA(VecNorm(tmp, NORM_2, norm, ierr))
53d8606c27SBarry Smith    PetscCallA(VecDestroy(tmp, ierr))
5424fb275aSStefano Zampini    write (outputString, *) norm
5524fb275aSStefano Zampini    PetscCallA(PetscPrintf(comm, 'Norm of search step '//trim(outputString)//'\n', ierr))
56c4762a1bSJed Brown  end
57c4762a1bSJed Brown
58*e7a95102SMartin Diehl! ---------------------------------------------------------------------
59*e7a95102SMartin Diehl!
60*e7a95102SMartin Diehl!  FormInitialGuess - Forms initial approximation.
61*e7a95102SMartin Diehl!
62*e7a95102SMartin Diehl!  Input Parameter:
63*e7a95102SMartin Diehl!  X - vector
64*e7a95102SMartin Diehl!
65*e7a95102SMartin Diehl!  Output Parameters:
66*e7a95102SMartin Diehl!  X - vector
67*e7a95102SMartin Diehl!  ierr - error code
68*e7a95102SMartin Diehl!
69*e7a95102SMartin Diehl!  Notes:
70*e7a95102SMartin Diehl!  This routine serves as a wrapper for the lower-level routine
71*e7a95102SMartin Diehl!  "ApplicationInitialGuess", where the actual computations are
72*e7a95102SMartin Diehl!  done using the standard Fortran style of treating the local
73*e7a95102SMartin Diehl!  vector data as a multidimensional array over the local mesh.
74*e7a95102SMartin Diehl!  This routine merely accesses the local vector data via
75*e7a95102SMartin Diehl!  VecGetArray() and VecRestoreArray().
76*e7a95102SMartin Diehl!
77*e7a95102SMartin Diehl  subroutine FormInitialGuess(X, ierr)
78*e7a95102SMartin Diehl
79*e7a95102SMartin Diehl!  Input/output variables:
80*e7a95102SMartin Diehl    Vec X
81*e7a95102SMartin Diehl    PetscErrorCode ierr
82*e7a95102SMartin Diehl
83*e7a95102SMartin Diehl!     Declarations for use with local arrays:
84*e7a95102SMartin Diehl    PetscScalar, pointer :: lx_v(:)
85*e7a95102SMartin Diehl
86*e7a95102SMartin Diehl    ierr = 0
87*e7a95102SMartin Diehl
88*e7a95102SMartin Diehl!  Get a pointer to vector data.
89*e7a95102SMartin Diehl!    - VecGetArray() returns a pointer to the data array.
90*e7a95102SMartin Diehl!    - You MUST call VecRestoreArray() when you no longer need access to
91*e7a95102SMartin Diehl!      the array.
92*e7a95102SMartin Diehl
93*e7a95102SMartin Diehl    PetscCallA(VecGetArray(X, lx_v, ierr))
94*e7a95102SMartin Diehl
95*e7a95102SMartin Diehl!  Compute initial guess
96*e7a95102SMartin Diehl
97*e7a95102SMartin Diehl    PetscCallA(ApplicationInitialGuess(lx_v, ierr))
98*e7a95102SMartin Diehl
99*e7a95102SMartin Diehl!  Restore vector
100*e7a95102SMartin Diehl
101*e7a95102SMartin Diehl    PetscCallA(VecRestoreArray(X, lx_v, ierr))
102*e7a95102SMartin Diehl
103*e7a95102SMartin Diehl  end
104*e7a95102SMartin Diehl
105*e7a95102SMartin Diehl!  ApplicationInitialGuess - Computes initial approximation, called by
106*e7a95102SMartin Diehl!  the higher level routine FormInitialGuess().
107*e7a95102SMartin Diehl!
108*e7a95102SMartin Diehl!  Input Parameter:
109*e7a95102SMartin Diehl!  x - local vector data
110*e7a95102SMartin Diehl!
111*e7a95102SMartin Diehl!  Output Parameters:
112*e7a95102SMartin Diehl!  f - local vector data, f(x)
113*e7a95102SMartin Diehl!  ierr - error code
114*e7a95102SMartin Diehl!
115*e7a95102SMartin Diehl!  Notes:
116*e7a95102SMartin Diehl!  This routine uses standard Fortran-style computations over a 2-dim array.
117*e7a95102SMartin Diehl!
118*e7a95102SMartin Diehl  subroutine ApplicationInitialGuess(x, ierr)
119*e7a95102SMartin Diehl
120*e7a95102SMartin Diehl!  Common blocks:
121*e7a95102SMartin Diehl    PetscReal lambda
122*e7a95102SMartin Diehl    PetscInt mx, my
123*e7a95102SMartin Diehl    PetscBool fd_coloring
124*e7a95102SMartin Diehl    common/params/lambda, mx, my, fd_coloring
125*e7a95102SMartin Diehl
126*e7a95102SMartin Diehl!  Input/output variables:
127*e7a95102SMartin Diehl    PetscScalar x(mx, my)
128*e7a95102SMartin Diehl    PetscErrorCode ierr
129*e7a95102SMartin Diehl
130*e7a95102SMartin Diehl!  Local variables:
131*e7a95102SMartin Diehl    PetscInt i, j
132*e7a95102SMartin Diehl    PetscReal temp1, temp, hx, hy, one
133*e7a95102SMartin Diehl
134*e7a95102SMartin Diehl!  Set parameters
135*e7a95102SMartin Diehl
136*e7a95102SMartin Diehl    ierr = 0
137*e7a95102SMartin Diehl    one = 1.0
138*e7a95102SMartin Diehl    hx = one/(mx - 1)
139*e7a95102SMartin Diehl    hy = one/(my - 1)
140*e7a95102SMartin Diehl    temp1 = lambda/(lambda + one)
141*e7a95102SMartin Diehl
142*e7a95102SMartin Diehl    do j = 1, my
143*e7a95102SMartin Diehl      temp = min(j - 1, my - j)*hy
144*e7a95102SMartin Diehl      do i = 1, mx
145*e7a95102SMartin Diehl        if (i == 1 .or. j == 1 .or. i == mx .or. j == my) then
146*e7a95102SMartin Diehl          x(i, j) = 0.0
147*e7a95102SMartin Diehl        else
148*e7a95102SMartin Diehl          x(i, j) = temp1*sqrt(min(min(i - 1, mx - i)*hx, temp))
149*e7a95102SMartin Diehl        end if
150*e7a95102SMartin Diehl      end do
151*e7a95102SMartin Diehl    end do
152*e7a95102SMartin Diehl
153*e7a95102SMartin Diehl  end
154*e7a95102SMartin Diehl
155*e7a95102SMartin Diehl! ---------------------------------------------------------------------
156*e7a95102SMartin Diehl!
157*e7a95102SMartin Diehl!  FormFunction - Evaluates nonlinear function, F(x).
158*e7a95102SMartin Diehl!
159*e7a95102SMartin Diehl!  Input Parameters:
160*e7a95102SMartin Diehl!  snes  - the SNES context
161*e7a95102SMartin Diehl!  X     - input vector
162*e7a95102SMartin Diehl!  dummy - optional user-defined context, as set by SNESSetFunction()
163*e7a95102SMartin Diehl!          (not used here)
164*e7a95102SMartin Diehl!
165*e7a95102SMartin Diehl!  Output Parameter:
166*e7a95102SMartin Diehl!  F     - vector with newly computed function
167*e7a95102SMartin Diehl!
168*e7a95102SMartin Diehl!  Notes:
169*e7a95102SMartin Diehl!  This routine serves as a wrapper for the lower-level routine
170*e7a95102SMartin Diehl!  "ApplicationFunction", where the actual computations are
171*e7a95102SMartin Diehl!  done using the standard Fortran style of treating the local
172*e7a95102SMartin Diehl!  vector data as a multidimensional array over the local mesh.
173*e7a95102SMartin Diehl!  This routine merely accesses the local vector data via
174*e7a95102SMartin Diehl!  VecGetArray() and VecRestoreArray().
175*e7a95102SMartin Diehl!
176*e7a95102SMartin Diehl  subroutine FormFunction(snes, X, F, fdcoloring, ierr)
177*e7a95102SMartin Diehl
178*e7a95102SMartin Diehl!  Input/output variables:
179*e7a95102SMartin Diehl    SNES snes
180*e7a95102SMartin Diehl    Vec X, F
181*e7a95102SMartin Diehl    PetscErrorCode ierr
182*e7a95102SMartin Diehl    MatFDColoring fdcoloring
183*e7a95102SMartin Diehl
184*e7a95102SMartin Diehl!  Common blocks:
185*e7a95102SMartin Diehl    PetscReal lambda
186*e7a95102SMartin Diehl    PetscInt mx, my
187*e7a95102SMartin Diehl    PetscBool fd_coloring
188*e7a95102SMartin Diehl    common/params/lambda, mx, my, fd_coloring
189*e7a95102SMartin Diehl
190*e7a95102SMartin Diehl!  Declarations for use with local arrays:
191*e7a95102SMartin Diehl    PetscScalar, pointer :: lx_v(:), lf_v(:)
192*e7a95102SMartin Diehl    PetscInt, pointer :: indices(:)
193*e7a95102SMartin Diehl
194*e7a95102SMartin Diehl!  Get pointers to vector data.
195*e7a95102SMartin Diehl!    - VecGetArray() returns a pointer to the data array.
196*e7a95102SMartin Diehl!    - You MUST call VecRestoreArray() when you no longer need access to
197*e7a95102SMartin Diehl!      the array.
198*e7a95102SMartin Diehl
199*e7a95102SMartin Diehl    PetscCallA(VecGetArrayRead(X, lx_v, ierr))
200*e7a95102SMartin Diehl    PetscCallA(VecGetArray(F, lf_v, ierr))
201*e7a95102SMartin Diehl
202*e7a95102SMartin Diehl!  Compute function
203*e7a95102SMartin Diehl
204*e7a95102SMartin Diehl    PetscCallA(ApplicationFunction(lx_v, lf_v, ierr))
205*e7a95102SMartin Diehl
206*e7a95102SMartin Diehl!  Restore vectors
207*e7a95102SMartin Diehl
208*e7a95102SMartin Diehl    PetscCallA(VecRestoreArrayRead(X, lx_v, ierr))
209*e7a95102SMartin Diehl    PetscCallA(VecRestoreArray(F, lf_v, ierr))
210*e7a95102SMartin Diehl
211*e7a95102SMartin Diehl    PetscCallA(PetscLogFlops(11.0d0*mx*my, ierr))
212*e7a95102SMartin Diehl!
213*e7a95102SMartin Diehl!     fdcoloring is in the common block and used here ONLY to test the
214*e7a95102SMartin Diehl!     calls to MatFDColoringGetPerturbedColumns() and  MatFDColoringRestorePerturbedColumns()
215*e7a95102SMartin Diehl!
216*e7a95102SMartin Diehl    if (fd_coloring) then
217*e7a95102SMartin Diehl      PetscCallA(MatFDColoringGetPerturbedColumns(fdcoloring, PETSC_NULL_INTEGER, indices, ierr))
218*e7a95102SMartin Diehl      print *, 'Indices from GetPerturbedColumns'
219*e7a95102SMartin Diehl      write (*, 1000) indices
220*e7a95102SMartin Diehl1000  format(50i4)
221*e7a95102SMartin Diehl      PetscCallA(MatFDColoringRestorePerturbedColumns(fdcoloring, PETSC_NULL_INTEGER, indices, ierr))
222*e7a95102SMartin Diehl    end if
223*e7a95102SMartin Diehl  end
224*e7a95102SMartin Diehl
225*e7a95102SMartin Diehl! ---------------------------------------------------------------------
226*e7a95102SMartin Diehl!
227*e7a95102SMartin Diehl!  ApplicationFunction - Computes nonlinear function, called by
228*e7a95102SMartin Diehl!  the higher level routine FormFunction().
229*e7a95102SMartin Diehl!
230*e7a95102SMartin Diehl!  Input Parameter:
231*e7a95102SMartin Diehl!  x    - local vector data
232*e7a95102SMartin Diehl!
233*e7a95102SMartin Diehl!  Output Parameters:
234*e7a95102SMartin Diehl!  f    - local vector data, f(x)
235*e7a95102SMartin Diehl!  ierr - error code
236*e7a95102SMartin Diehl!
237*e7a95102SMartin Diehl!  Notes:
238*e7a95102SMartin Diehl!  This routine uses standard Fortran-style computations over a 2-dim array.
239*e7a95102SMartin Diehl!
240*e7a95102SMartin Diehl  subroutine ApplicationFunction(x, f, ierr)
241*e7a95102SMartin Diehl
242*e7a95102SMartin Diehl!  Common blocks:
243*e7a95102SMartin Diehl    PetscReal lambda
244*e7a95102SMartin Diehl    PetscInt mx, my
245*e7a95102SMartin Diehl    PetscBool fd_coloring
246*e7a95102SMartin Diehl    common/params/lambda, mx, my, fd_coloring
247*e7a95102SMartin Diehl
248*e7a95102SMartin Diehl!  Input/output variables:
249*e7a95102SMartin Diehl    PetscScalar x(mx, my), f(mx, my)
250*e7a95102SMartin Diehl    PetscErrorCode ierr
251*e7a95102SMartin Diehl
252*e7a95102SMartin Diehl!  Local variables:
253*e7a95102SMartin Diehl    PetscScalar two, one, hx, hy
254*e7a95102SMartin Diehl    PetscScalar hxdhy, hydhx, sc
255*e7a95102SMartin Diehl    PetscScalar u, uxx, uyy
256*e7a95102SMartin Diehl    PetscInt i, j
257*e7a95102SMartin Diehl
258*e7a95102SMartin Diehl    ierr = 0
259*e7a95102SMartin Diehl    one = 1.0
260*e7a95102SMartin Diehl    two = 2.0
261*e7a95102SMartin Diehl    hx = one/(mx - 1)
262*e7a95102SMartin Diehl    hy = one/(my - 1)
263*e7a95102SMartin Diehl    sc = hx*hy*lambda
264*e7a95102SMartin Diehl    hxdhy = hx/hy
265*e7a95102SMartin Diehl    hydhx = hy/hx
266*e7a95102SMartin Diehl
267*e7a95102SMartin Diehl!  Compute function
268*e7a95102SMartin Diehl
269*e7a95102SMartin Diehl    do j = 1, my
270*e7a95102SMartin Diehl      do i = 1, mx
271*e7a95102SMartin Diehl        if (i == 1 .or. j == 1 .or. i == mx .or. j == my) then
272*e7a95102SMartin Diehl          f(i, j) = x(i, j)
273*e7a95102SMartin Diehl        else
274*e7a95102SMartin Diehl          u = x(i, j)
275*e7a95102SMartin Diehl          uxx = hydhx*(two*u - x(i - 1, j) - x(i + 1, j))
276*e7a95102SMartin Diehl          uyy = hxdhy*(two*u - x(i, j - 1) - x(i, j + 1))
277*e7a95102SMartin Diehl          f(i, j) = uxx + uyy - sc*exp(u)
278*e7a95102SMartin Diehl        end if
279*e7a95102SMartin Diehl      end do
280*e7a95102SMartin Diehl    end do
281*e7a95102SMartin Diehl
282*e7a95102SMartin Diehl  end
283*e7a95102SMartin Diehl
284*e7a95102SMartin Diehl! ---------------------------------------------------------------------
285*e7a95102SMartin Diehl!
286*e7a95102SMartin Diehl!  FormJacobian - Evaluates Jacobian matrix.
287*e7a95102SMartin Diehl!
288*e7a95102SMartin Diehl!  Input Parameters:
289*e7a95102SMartin Diehl!  snes    - the SNES context
290*e7a95102SMartin Diehl!  x       - input vector
291*e7a95102SMartin Diehl!  dummy   - optional user-defined context, as set by SNESSetJacobian()
292*e7a95102SMartin Diehl!            (not used here)
293*e7a95102SMartin Diehl!
294*e7a95102SMartin Diehl!  Output Parameters:
295*e7a95102SMartin Diehl!  jac      - Jacobian matrix
296*e7a95102SMartin Diehl!  jac_prec - optionally different matrix used to construct the preconditioner (not used here)
297*e7a95102SMartin Diehl!
298*e7a95102SMartin Diehl!  Notes:
299*e7a95102SMartin Diehl!  This routine serves as a wrapper for the lower-level routine
300*e7a95102SMartin Diehl!  "ApplicationJacobian", where the actual computations are
301*e7a95102SMartin Diehl!  done using the standard Fortran style of treating the local
302*e7a95102SMartin Diehl!  vector data as a multidimensional array over the local mesh.
303*e7a95102SMartin Diehl!  This routine merely accesses the local vector data via
304*e7a95102SMartin Diehl!  VecGetArray() and VecRestoreArray().
305*e7a95102SMartin Diehl!
306*e7a95102SMartin Diehl  subroutine FormJacobian(snes, X, jac, jac_prec, dummy, ierr)
307*e7a95102SMartin Diehl
308*e7a95102SMartin Diehl!  Input/output variables:
309*e7a95102SMartin Diehl    SNES snes
310*e7a95102SMartin Diehl    Vec X
311*e7a95102SMartin Diehl    Mat jac, jac_prec
312*e7a95102SMartin Diehl    PetscErrorCode ierr
313*e7a95102SMartin Diehl    integer dummy
314*e7a95102SMartin Diehl
315*e7a95102SMartin Diehl!  Common blocks:
316*e7a95102SMartin Diehl    PetscReal lambda
317*e7a95102SMartin Diehl    PetscInt mx, my
318*e7a95102SMartin Diehl    PetscBool fd_coloring
319*e7a95102SMartin Diehl    common/params/lambda, mx, my, fd_coloring
320*e7a95102SMartin Diehl
321*e7a95102SMartin Diehl!  Declarations for use with local array:
322*e7a95102SMartin Diehl    PetscScalar, pointer :: lx_v(:)
323*e7a95102SMartin Diehl
324*e7a95102SMartin Diehl!  Get a pointer to vector data
325*e7a95102SMartin Diehl
326*e7a95102SMartin Diehl    PetscCallA(VecGetArrayRead(X, lx_v, ierr))
327*e7a95102SMartin Diehl
328*e7a95102SMartin Diehl!  Compute Jacobian entries
329*e7a95102SMartin Diehl
330*e7a95102SMartin Diehl    PetscCallA(ApplicationJacobian(lx_v, jac, jac_prec, ierr))
331*e7a95102SMartin Diehl
332*e7a95102SMartin Diehl!  Restore vector
333*e7a95102SMartin Diehl
334*e7a95102SMartin Diehl    PetscCallA(VecRestoreArrayRead(X, lx_v, ierr))
335*e7a95102SMartin Diehl
336*e7a95102SMartin Diehl!  Assemble matrix
337*e7a95102SMartin Diehl
338*e7a95102SMartin Diehl    PetscCallA(MatAssemblyBegin(jac_prec, MAT_FINAL_ASSEMBLY, ierr))
339*e7a95102SMartin Diehl    PetscCallA(MatAssemblyEnd(jac_prec, MAT_FINAL_ASSEMBLY, ierr))
340*e7a95102SMartin Diehl
341*e7a95102SMartin Diehl  end
342*e7a95102SMartin Diehl
343*e7a95102SMartin Diehl! ---------------------------------------------------------------------
344*e7a95102SMartin Diehl!
345*e7a95102SMartin Diehl!  ApplicationJacobian - Computes Jacobian matrix, called by
346*e7a95102SMartin Diehl!  the higher level routine FormJacobian().
347*e7a95102SMartin Diehl!
348*e7a95102SMartin Diehl!  Input Parameters:
349*e7a95102SMartin Diehl!  x        - local vector data
350*e7a95102SMartin Diehl!
351*e7a95102SMartin Diehl!  Output Parameters:
352*e7a95102SMartin Diehl!  jac      - Jacobian matrix
353*e7a95102SMartin Diehl!  jac_prec - optionally different matrix used to construct the preconditioner (not used here)
354*e7a95102SMartin Diehl!  ierr     - error code
355*e7a95102SMartin Diehl!
356*e7a95102SMartin Diehl!  Notes:
357*e7a95102SMartin Diehl!  This routine uses standard Fortran-style computations over a 2-dim array.
358*e7a95102SMartin Diehl!
359*e7a95102SMartin Diehl  subroutine ApplicationJacobian(x, jac, jac_prec, ierr)
360*e7a95102SMartin Diehl!  Common blocks:
361*e7a95102SMartin Diehl    PetscReal lambda
362*e7a95102SMartin Diehl    PetscInt mx, my
363*e7a95102SMartin Diehl    PetscBool fd_coloring
364*e7a95102SMartin Diehl    common/params/lambda, mx, my, fd_coloring
365*e7a95102SMartin Diehl
366*e7a95102SMartin Diehl!  Input/output variables:
367*e7a95102SMartin Diehl    PetscScalar x(mx, my)
368*e7a95102SMartin Diehl    Mat jac, jac_prec
369*e7a95102SMartin Diehl    PetscErrorCode ierr
370*e7a95102SMartin Diehl
371*e7a95102SMartin Diehl!  Local variables:
372*e7a95102SMartin Diehl    PetscInt i, j, row(1), col(5), i1, i5
373*e7a95102SMartin Diehl    PetscScalar two, one, hx, hy
374*e7a95102SMartin Diehl    PetscScalar hxdhy, hydhx, sc, v(5)
375*e7a95102SMartin Diehl
376*e7a95102SMartin Diehl!  Set parameters
377*e7a95102SMartin Diehl
378*e7a95102SMartin Diehl    i1 = 1
379*e7a95102SMartin Diehl    i5 = 5
380*e7a95102SMartin Diehl    one = 1.0
381*e7a95102SMartin Diehl    two = 2.0
382*e7a95102SMartin Diehl    hx = one/(mx - 1)
383*e7a95102SMartin Diehl    hy = one/(my - 1)
384*e7a95102SMartin Diehl    sc = hx*hy
385*e7a95102SMartin Diehl    hxdhy = hx/hy
386*e7a95102SMartin Diehl    hydhx = hy/hx
387*e7a95102SMartin Diehl
388*e7a95102SMartin Diehl!  Compute entries of the Jacobian matrix
389*e7a95102SMartin Diehl!   - Here, we set all entries for a particular row at once.
390*e7a95102SMartin Diehl!   - Note that MatSetValues() uses 0-based row and column numbers
391*e7a95102SMartin Diehl!     in Fortran as well as in C.
392*e7a95102SMartin Diehl
393*e7a95102SMartin Diehl    do j = 1, my
394*e7a95102SMartin Diehl      row(1) = (j - 1)*mx - 1
395*e7a95102SMartin Diehl      do i = 1, mx
396*e7a95102SMartin Diehl        row(1) = row(1) + 1
397*e7a95102SMartin Diehl!           boundary points
398*e7a95102SMartin Diehl        if (i == 1 .or. j == 1 .or. i == mx .or. j == my) then
399*e7a95102SMartin Diehl          PetscCallA(MatSetValues(jac_prec, i1, row, i1, row, [one], INSERT_VALUES, ierr))
400*e7a95102SMartin Diehl!           interior grid points
401*e7a95102SMartin Diehl        else
402*e7a95102SMartin Diehl          v(1) = -hxdhy
403*e7a95102SMartin Diehl          v(2) = -hydhx
404*e7a95102SMartin Diehl          v(3) = two*(hydhx + hxdhy) - sc*lambda*exp(x(i, j))
405*e7a95102SMartin Diehl          v(4) = -hydhx
406*e7a95102SMartin Diehl          v(5) = -hxdhy
407*e7a95102SMartin Diehl          col(1) = row(1) - mx
408*e7a95102SMartin Diehl          col(2) = row(1) - 1
409*e7a95102SMartin Diehl          col(3) = row(1)
410*e7a95102SMartin Diehl          col(4) = row(1) + 1
411*e7a95102SMartin Diehl          col(5) = row(1) + mx
412*e7a95102SMartin Diehl          PetscCallA(MatSetValues(jac_prec, i1, row, i5, col, v, INSERT_VALUES, ierr))
413*e7a95102SMartin Diehl        end if
414*e7a95102SMartin Diehl      end do
415*e7a95102SMartin Diehl    end do
416*e7a95102SMartin Diehl
417*e7a95102SMartin Diehl  end
418*e7a95102SMartin Diehl
419*e7a95102SMartin Diehlend module ex1f_mod
420c4762a1bSJed Brownprogram main
421ce78bad3SBarry Smith  use petscdraw
422c4762a1bSJed Brown  use petscsnes
423*e7a95102SMartin Diehl  use ex1f_mod
424c4762a1bSJed Brown  implicit none
425c4762a1bSJed Brown!
426c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
427c4762a1bSJed Brown!                   Variable declarations
428c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
429c4762a1bSJed Brown!
430c4762a1bSJed Brown!  Variables:
431c4762a1bSJed Brown!     snes        - nonlinear solver
432c4762a1bSJed Brown!     x, r        - solution, residual vectors
433c4762a1bSJed Brown!     J           - Jacobian matrix
434c4762a1bSJed Brown!     its         - iterations for convergence
435c4762a1bSJed Brown!     matrix_free - flag - 1 indicates matrix-free version
436c4762a1bSJed Brown!     lambda      - nonlinearity parameter
437c4762a1bSJed Brown!     draw        - drawing context
438c4762a1bSJed Brown!
439c4762a1bSJed Brown  SNES snes
440c4762a1bSJed Brown  MatColoring mc
441c4762a1bSJed Brown  Vec x, r
442c4762a1bSJed Brown  PetscDraw draw
443c4762a1bSJed Brown  Mat J
444c4762a1bSJed Brown  PetscBool matrix_free, flg, fd_coloring
445c4762a1bSJed Brown  PetscErrorCode ierr
446c4762a1bSJed Brown  PetscInt its, N, mx, my, i5
447c4762a1bSJed Brown  PetscMPIInt size, rank
448c4762a1bSJed Brown  PetscReal lambda_max, lambda_min, lambda
449c4762a1bSJed Brown  MatFDColoring fdcoloring
450c4762a1bSJed Brown  ISColoring iscoloring
451c4762a1bSJed Brown  PetscBool pc
452ce78bad3SBarry Smith  integer4 imx, imy
45324fb275aSStefano Zampini  character(len=PETSC_MAX_PATH_LEN) :: outputString
45442ce371bSBarry Smith  PetscScalar, pointer :: lx_v(:)
455ce78bad3SBarry Smith  integer4 xl, yl, width, height
456c4762a1bSJed Brown
457c4762a1bSJed Brown!  Store parameters in common block
458c4762a1bSJed Brown
459c4762a1bSJed Brown  common/params/lambda, mx, my, fd_coloring
460c4762a1bSJed Brown
461c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
462c4762a1bSJed Brown!  Initialize program
463c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
464c4762a1bSJed Brown
465d8606c27SBarry Smith  PetscCallA(PetscInitialize(ierr))
466d8606c27SBarry Smith  PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD, size, ierr))
467d8606c27SBarry Smith  PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr))
468c4762a1bSJed Brown
4694820e4eaSBarry Smith  PetscCheckA(size == 1, PETSC_COMM_SELF, PETSC_ERR_WRONG_MPI_SIZE, 'This is a uniprocessor example only')
470c4762a1bSJed Brown
471c4762a1bSJed Brown!  Initialize problem parameters
472c4762a1bSJed Brown  i5 = 5
473c4762a1bSJed Brown  lambda_max = 6.81
474c4762a1bSJed Brown  lambda_min = 0.0
475c4762a1bSJed Brown  lambda = 6.0
476c4762a1bSJed Brown  mx = 4
477c4762a1bSJed Brown  my = 4
478d8606c27SBarry Smith  PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-mx', mx, flg, ierr))
479d8606c27SBarry Smith  PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-my', my, flg, ierr))
480d8606c27SBarry Smith  PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-par', lambda, flg, ierr))
4814820e4eaSBarry Smith  PetscCheckA(lambda < lambda_max .and. lambda > lambda_min, PETSC_COMM_SELF, PETSC_ERR_USER, 'Lambda out of range ')
482c4762a1bSJed Brown  N = mx*my
483c4762a1bSJed Brown  pc = PETSC_FALSE
484d8606c27SBarry Smith  PetscCallA(PetscOptionsGetBool(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-pc', pc, PETSC_NULL_BOOL, ierr))
485c4762a1bSJed Brown
486c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
487c4762a1bSJed Brown!  Create nonlinear solver context
488c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
489c4762a1bSJed Brown
490d8606c27SBarry Smith  PetscCallA(SNESCreate(PETSC_COMM_WORLD, snes, ierr))
491c4762a1bSJed Brown
492c4762a1bSJed Brown  if (pc .eqv. PETSC_TRUE) then
493d8606c27SBarry Smith    PetscCallA(SNESSetType(snes, SNESNEWTONTR, ierr))
494d8606c27SBarry Smith    PetscCallA(SNESNewtonTRSetPostCheck(snes, postcheck, snes, ierr))
495c4762a1bSJed Brown  end if
496c4762a1bSJed Brown
497c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
498c4762a1bSJed Brown!  Create vector data structures; set function evaluation routine
499c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
500c4762a1bSJed Brown
501d8606c27SBarry Smith  PetscCallA(VecCreate(PETSC_COMM_WORLD, x, ierr))
502d8606c27SBarry Smith  PetscCallA(VecSetSizes(x, PETSC_DECIDE, N, ierr))
503d8606c27SBarry Smith  PetscCallA(VecSetFromOptions(x, ierr))
504d8606c27SBarry Smith  PetscCallA(VecDuplicate(x, r, ierr))
505c4762a1bSJed Brown
506c4762a1bSJed Brown!  Set function evaluation routine and vector.  Whenever the nonlinear
507c4762a1bSJed Brown!  solver needs to evaluate the nonlinear function, it will call this
508c4762a1bSJed Brown!  routine.
509c4762a1bSJed Brown!   - Note that the final routine argument is the user-defined
510c4762a1bSJed Brown!     context that provides application-specific data for the
511c4762a1bSJed Brown!     function evaluation routine.
512c4762a1bSJed Brown
513d8606c27SBarry Smith  PetscCallA(SNESSetFunction(snes, r, FormFunction, fdcoloring, ierr))
514c4762a1bSJed Brown
515c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
516c4762a1bSJed Brown!  Create matrix data structure; set Jacobian evaluation routine
517c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
518c4762a1bSJed Brown
519c4762a1bSJed Brown!  Create matrix. Here we only approximately preallocate storage space
520c4762a1bSJed Brown!  for the Jacobian.  See the users manual for a discussion of better
521c4762a1bSJed Brown!  techniques for preallocating matrix memory.
522c4762a1bSJed Brown
523d8606c27SBarry Smith  PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-snes_mf', matrix_free, ierr))
524c4762a1bSJed Brown  if (.not. matrix_free) then
5255d83a8b1SBarry Smith    PetscCallA(MatCreateSeqAIJ(PETSC_COMM_WORLD, N, N, i5, PETSC_NULL_INTEGER_ARRAY, J, ierr))
526c4762a1bSJed Brown  end if
527c4762a1bSJed Brown
528c4762a1bSJed Brown!
529c4762a1bSJed Brown!     This option will cause the Jacobian to be computed via finite differences
530c4762a1bSJed Brown!    efficiently using a coloring of the columns of the matrix.
531c4762a1bSJed Brown!
532c4762a1bSJed Brown  fd_coloring = .false.
533d8606c27SBarry Smith  PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-snes_fd_coloring', fd_coloring, ierr))
534c4762a1bSJed Brown  if (fd_coloring) then
535c4762a1bSJed Brown
536c4762a1bSJed Brown!
537c4762a1bSJed Brown!      This initializes the nonzero structure of the Jacobian. This is artificial
538c4762a1bSJed Brown!      because clearly if we had a routine to compute the Jacobian we won't need
539c4762a1bSJed Brown!      to use finite differences.
540c4762a1bSJed Brown!
541d8606c27SBarry Smith    PetscCallA(FormJacobian(snes, x, J, J, 0, ierr))
542c4762a1bSJed Brown!
543c4762a1bSJed Brown!       Color the matrix, i.e. determine groups of columns that share no common
544a5b23f4aSJose E. Roman!      rows. These columns in the Jacobian can all be computed simultaneously.
545c4762a1bSJed Brown!
546d8606c27SBarry Smith    PetscCallA(MatColoringCreate(J, mc, ierr))
547d8606c27SBarry Smith    PetscCallA(MatColoringSetType(mc, MATCOLORINGNATURAL, ierr))
548d8606c27SBarry Smith    PetscCallA(MatColoringSetFromOptions(mc, ierr))
549d8606c27SBarry Smith    PetscCallA(MatColoringApply(mc, iscoloring, ierr))
550d8606c27SBarry Smith    PetscCallA(MatColoringDestroy(mc, ierr))
551c4762a1bSJed Brown!
552c4762a1bSJed Brown!       Create the data structure that SNESComputeJacobianDefaultColor() uses
553c4762a1bSJed Brown!       to compute the actual Jacobians via finite differences.
554c4762a1bSJed Brown!
555d8606c27SBarry Smith    PetscCallA(MatFDColoringCreate(J, iscoloring, fdcoloring, ierr))
556d8606c27SBarry Smith    PetscCallA(MatFDColoringSetFunction(fdcoloring, FormFunction, fdcoloring, ierr))
557d8606c27SBarry Smith    PetscCallA(MatFDColoringSetFromOptions(fdcoloring, ierr))
558d8606c27SBarry Smith    PetscCallA(MatFDColoringSetUp(J, iscoloring, fdcoloring, ierr))
559c4762a1bSJed Brown!
560c4762a1bSJed Brown!        Tell SNES to use the routine SNESComputeJacobianDefaultColor()
561c4762a1bSJed Brown!      to compute Jacobians.
562c4762a1bSJed Brown!
563d8606c27SBarry Smith    PetscCallA(SNESSetJacobian(snes, J, J, SNESComputeJacobianDefaultColor, fdcoloring, ierr))
564f51a5268SBarry Smith    PetscCallA(SNESGetJacobian(snes, J, PETSC_NULL_MAT, PETSC_NULL_FUNCTION, PETSC_NULL_INTEGER, ierr))
565f51a5268SBarry Smith    PetscCallA(SNESGetJacobian(snes, J, PETSC_NULL_MAT, PETSC_NULL_FUNCTION, fdcoloring, ierr))
566d8606c27SBarry Smith    PetscCallA(ISColoringDestroy(iscoloring, ierr))
567c4762a1bSJed Brown
568c4762a1bSJed Brown  else if (.not. matrix_free) then
569c4762a1bSJed Brown
570c4762a1bSJed Brown!  Set Jacobian matrix data structure and default Jacobian evaluation
571c4762a1bSJed Brown!  routine.  Whenever the nonlinear solver needs to compute the
572c4762a1bSJed Brown!  Jacobian matrix, it will call this routine.
573c4762a1bSJed Brown!   - Note that the final routine argument is the user-defined
574c4762a1bSJed Brown!     context that provides application-specific data for the
575c4762a1bSJed Brown!     Jacobian evaluation routine.
576c4762a1bSJed Brown!   - The user can override with:
577c4762a1bSJed Brown!      -snes_fd : default finite differencing approximation of Jacobian
578c4762a1bSJed Brown!      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
579c4762a1bSJed Brown!                 (unless user explicitly sets preconditioner)
5807addb90fSBarry Smith!      -snes_mf_operator : form matrix from which to construct the preconditioner as set by the user,
581c4762a1bSJed Brown!                          but use matrix-free approx for Jacobian-vector
582c4762a1bSJed Brown!                          products within Newton-Krylov method
583c4762a1bSJed Brown!
584d8606c27SBarry Smith    PetscCallA(SNESSetJacobian(snes, J, J, FormJacobian, 0, ierr))
585c4762a1bSJed Brown  end if
586c4762a1bSJed Brown
587c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
588c4762a1bSJed Brown!  Customize nonlinear solver; set runtime options
589c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
590c4762a1bSJed Brown
591c4762a1bSJed Brown!  Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
592c4762a1bSJed Brown
593d8606c27SBarry Smith  PetscCallA(SNESSetFromOptions(snes, ierr))
594c4762a1bSJed Brown
595c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
596c4762a1bSJed Brown!  Evaluate initial guess; then solve nonlinear system.
597c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
598c4762a1bSJed Brown
599c4762a1bSJed Brown!  Note: The user should initialize the vector, x, with the initial guess
600c4762a1bSJed Brown!  for the nonlinear solver prior to calling SNESSolve().  In particular,
601c4762a1bSJed Brown!  to employ an initial guess of zero, the user should explicitly set
602c4762a1bSJed Brown!  this vector to zero by calling VecSet().
603c4762a1bSJed Brown
604d8606c27SBarry Smith  PetscCallA(FormInitialGuess(x, ierr))
605d8606c27SBarry Smith  PetscCallA(SNESSolve(snes, PETSC_NULL_VEC, x, ierr))
606d8606c27SBarry Smith  PetscCallA(SNESGetIterationNumber(snes, its, ierr))
60724fb275aSStefano Zampini  write (outputString, *) its
60824fb275aSStefano Zampini  PetscCallA(PetscPrintf(PETSC_COMM_WORLD, 'Number of SNES iterations = '//trim(outputString)//'\n', ierr))
609c4762a1bSJed Brown
610c4762a1bSJed Brown!  PetscDraw contour plot of solution
611c4762a1bSJed Brown
612ce78bad3SBarry Smith  xl = 300
613ce78bad3SBarry Smith  yl = 0
614ce78bad3SBarry Smith  width = 300
615ce78bad3SBarry Smith  height = 300
616ce78bad3SBarry Smith  PetscCallA(PetscDrawCreate(PETSC_COMM_WORLD, PETSC_NULL_CHARACTER, 'Solution', xl, yl, width, height, draw, ierr))
617d8606c27SBarry Smith  PetscCallA(PetscDrawSetFromOptions(draw, ierr))
618c4762a1bSJed Brown
619ce78bad3SBarry Smith  PetscCallA(VecGetArrayRead(x, lx_v, ierr))
620ce78bad3SBarry Smith  imx = int(mx, kind=kind(imx))
621ce78bad3SBarry Smith  imy = int(my, kind=kind(imy))
622ce78bad3SBarry Smith  PetscCallA(PetscDrawTensorContour(draw, imx, imy, PETSC_NULL_REAL_ARRAY, PETSC_NULL_REAL_ARRAY, lx_v, ierr))
623ce78bad3SBarry Smith  PetscCallA(VecRestoreArrayRead(x, lx_v, ierr))
624c4762a1bSJed Brown
625c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
626c4762a1bSJed Brown!  Free work space.  All PETSc objects should be destroyed when they
627c4762a1bSJed Brown!  are no longer needed.
628c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
629c4762a1bSJed Brown
630d8606c27SBarry Smith  if (.not. matrix_free) PetscCallA(MatDestroy(J, ierr))
631d8606c27SBarry Smith  if (fd_coloring) PetscCallA(MatFDColoringDestroy(fdcoloring, ierr))
632c4762a1bSJed Brown
633d8606c27SBarry Smith  PetscCallA(VecDestroy(x, ierr))
634d8606c27SBarry Smith  PetscCallA(VecDestroy(r, ierr))
635d8606c27SBarry Smith  PetscCallA(SNESDestroy(snes, ierr))
636d8606c27SBarry Smith  PetscCallA(PetscDrawDestroy(draw, ierr))
637d8606c27SBarry Smith  PetscCallA(PetscFinalize(ierr))
638c4762a1bSJed Brownend
639c4762a1bSJed Brown!
640c4762a1bSJed Brown!/*TEST
641c4762a1bSJed Brown!
642c4762a1bSJed Brown!   build:
643ce78bad3SBarry Smith!      requires: !single !complex
644c4762a1bSJed Brown!
645c4762a1bSJed Brown!   test:
646c4762a1bSJed Brown!      args: -snes_monitor_short -nox -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always
647c4762a1bSJed Brown!
648c4762a1bSJed Brown!   test:
649c4762a1bSJed Brown!      suffix: 2
650c4762a1bSJed Brown!      args: -snes_monitor_short -nox -snes_fd -ksp_gmres_cgs_refinement_type refine_always
651c4762a1bSJed Brown!
652c4762a1bSJed Brown!   test:
653c4762a1bSJed Brown!      suffix: 3
654c4762a1bSJed Brown!      args: -snes_monitor_short -nox -snes_fd_coloring -mat_coloring_type sl -ksp_gmres_cgs_refinement_type refine_always
655c4762a1bSJed Brown!      filter: sort -b
656c4762a1bSJed Brown!      filter_output: sort -b
657c4762a1bSJed Brown!
658c4762a1bSJed Brown!   test:
659c4762a1bSJed Brown!     suffix: 4
660c4762a1bSJed Brown!     args: -pc -par 6.807 -nox
661c4762a1bSJed Brown!
662c4762a1bSJed Brown!TEST*/
663