xref: /petsc/src/snes/tutorials/ex73f90t.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!  This system (A) is augmented with constraints:
10c4762a1bSJed Brown!
11c4762a1bSJed Brown!    A -B   *  phi  =  rho
12c4762a1bSJed Brown!   -C  I      lam  = 0
13c4762a1bSJed Brown!
14c4762a1bSJed Brown!  where I is the identity, A is the "normal" Poisson equation, B is the "distributor" of the
15c4762a1bSJed Brown!  total flux (the first block equation is the flux surface averaging equation).  The second
16c4762a1bSJed Brown!  equation  lambda = C * x enforces the surface flux auxiliary equation.  B and C have all
17c4762a1bSJed Brown!  positive entries, areas in C and fraction of area in B.
18c4762a1bSJed Brown!
19c4762a1bSJed Brown
20c4762a1bSJed Brown!
21c4762a1bSJed Brown!  --------------------------------------------------------------------------
22c4762a1bSJed Brown!
23c4762a1bSJed Brown!  Solid Fuel Ignition (SFI) problem.  This problem is modeled by
24c4762a1bSJed Brown!  the partial differential equation
25c4762a1bSJed Brown!
26c4762a1bSJed Brown!          -Laplacian u - lambda*exp(u) = 0,  0 < x,y < 1,
27c4762a1bSJed Brown!
28c4762a1bSJed Brown!  with boundary conditions
29c4762a1bSJed Brown!
30c4762a1bSJed Brown!           u = 0  for  x = 0, x = 1, y = 0, y = 1.
31c4762a1bSJed Brown!
32c4762a1bSJed Brown!  A finite difference approximation with the usual 5-point stencil
33c4762a1bSJed Brown!  is used to discretize the boundary value problem to obtain a nonlinear
34c4762a1bSJed Brown!  system of equations.
35c4762a1bSJed Brown!
36c4762a1bSJed Brown!  --------------------------------------------------------------------------
37c4762a1bSJed Brown!  The following define must be used before including any PETSc include files
38c4762a1bSJed Brown!  into a module or interface. This is because they can't handle declarations
39c4762a1bSJed Brown!  in them
40c4762a1bSJed Brown!
41c5e229c2SMartin Diehl#include <petsc/finclude/petsc.h>
4277d968b7SBarry Smithmodule ex73f90tmodule
43ce78bad3SBarry Smith  use petscdmda
44ce78bad3SBarry Smith  use petscdmcomposite
45c4762a1bSJed Brown  use petscmat
46*e7a95102SMartin Diehl  use petscsys
47*e7a95102SMartin Diehl  use petscsnes
48*e7a95102SMartin Diehl  implicit none
4977d968b7SBarry Smith  type ex73f90tmodule_type
50c4762a1bSJed Brown    DM::da
51c4762a1bSJed Brown!     temp A block stuff
52c4762a1bSJed Brown    PetscInt mx, my
53c4762a1bSJed Brown    PetscMPIInt rank
54c4762a1bSJed Brown    PetscReal lambda
55c4762a1bSJed Brown!     Mats
56c4762a1bSJed Brown    Mat::Amat, AmatLin, Bmat, CMat, Dmat
57c4762a1bSJed Brown    IS::isPhi, isLambda
5877d968b7SBarry Smith  end type ex73f90tmodule_type
59c4762a1bSJed Brown
60*e7a95102SMartin Diehlcontains
61c00ad2bcSBarry Smith  subroutine MyObjective(snes, x, result, ctx, ierr)
62c00ad2bcSBarry Smith    PetscInt ctx
63c00ad2bcSBarry Smith    Vec x, f
64c00ad2bcSBarry Smith    SNES snes
65c00ad2bcSBarry Smith    PetscErrorCode ierr
66c00ad2bcSBarry Smith    PetscScalar result
67c00ad2bcSBarry Smith    PetscReal fnorm
68c00ad2bcSBarry Smith
69c00ad2bcSBarry Smith    PetscCall(VecDuplicate(x, f, ierr))
70c00ad2bcSBarry Smith    PetscCall(SNESComputeFunction(snes, x, f, ierr))
71c00ad2bcSBarry Smith    PetscCall(VecNorm(f, NORM_2, fnorm, ierr))
72c00ad2bcSBarry Smith    result = .5*fnorm*fnorm
73c00ad2bcSBarry Smith    PetscCall(VecDestroy(f, ierr))
74c00ad2bcSBarry Smith  end subroutine MyObjective
75c00ad2bcSBarry Smith
76*e7a95102SMartin Diehl! ---------------------------------------------------------------------
77*e7a95102SMartin Diehl!
78*e7a95102SMartin Diehl!  FormInitialGuess - Forms initial approximation.
79*e7a95102SMartin Diehl!
80*e7a95102SMartin Diehl!  Input Parameters:
81*e7a95102SMartin Diehl!  X - vector
82*e7a95102SMartin Diehl!
83*e7a95102SMartin Diehl!  Output Parameter:
84*e7a95102SMartin Diehl!  X - vector
85*e7a95102SMartin Diehl!
86*e7a95102SMartin Diehl!  Notes:
87*e7a95102SMartin Diehl!  This routine serves as a wrapper for the lower-level routine
88*e7a95102SMartin Diehl!  "InitialGuessLocal", where the actual computations are
89*e7a95102SMartin Diehl!  done using the standard Fortran style of treating the local
90*e7a95102SMartin Diehl!  vector data as a multidimensional array over the local mesh.
91*e7a95102SMartin Diehl!  This routine merely handles ghost point scatters and accesses
92*e7a95102SMartin Diehl!  the local vector data via VecGetArray() and VecRestoreArray().
93*e7a95102SMartin Diehl!
94*e7a95102SMartin Diehl  subroutine FormInitialGuess(mysnes, Xnest, ierr)
95*e7a95102SMartin Diehl!  Input/output variables:
96*e7a95102SMartin Diehl    SNES::     mysnes
97*e7a95102SMartin Diehl    Vec::      Xnest
98*e7a95102SMartin Diehl    PetscErrorCode ierr
99*e7a95102SMartin Diehl
100*e7a95102SMartin Diehl!  Declarations for use with local arrays:
101*e7a95102SMartin Diehl    type(ex73f90tmodule_type), pointer:: solver
102*e7a95102SMartin Diehl    Vec::      Xsub(2)
103*e7a95102SMartin Diehl    PetscInt::  izero, ione, itwo
104*e7a95102SMartin Diehl
105*e7a95102SMartin Diehl    izero = 0
106*e7a95102SMartin Diehl    ione = 1
107*e7a95102SMartin Diehl    itwo = 2
108*e7a95102SMartin Diehl    ierr = 0
109*e7a95102SMartin Diehl    PetscCall(SNESGetApplicationContext(mysnes, solver, ierr))
110*e7a95102SMartin Diehl    PetscCall(DMCompositeGetAccessArray(solver%da, Xnest, itwo, PETSC_NULL_INTEGER_ARRAY, Xsub, ierr))
111*e7a95102SMartin Diehl
112*e7a95102SMartin Diehl    PetscCall(InitialGuessLocal(solver, Xsub(1), ierr))
113*e7a95102SMartin Diehl    PetscCall(VecAssemblyBegin(Xsub(1), ierr))
114*e7a95102SMartin Diehl    PetscCall(VecAssemblyEnd(Xsub(1), ierr))
115*e7a95102SMartin Diehl
116*e7a95102SMartin Diehl!     zero out lambda
117*e7a95102SMartin Diehl    PetscCall(VecZeroEntries(Xsub(2), ierr))
118*e7a95102SMartin Diehl    PetscCall(DMCompositeRestoreAccessArray(solver%da, Xnest, itwo, PETSC_NULL_INTEGER_ARRAY, Xsub, ierr))
119*e7a95102SMartin Diehl
120*e7a95102SMartin Diehl  end subroutine FormInitialGuess
121*e7a95102SMartin Diehl
122*e7a95102SMartin Diehl! ---------------------------------------------------------------------
123*e7a95102SMartin Diehl!
124*e7a95102SMartin Diehl!  InitialGuessLocal - Computes initial approximation, called by
125*e7a95102SMartin Diehl!  the higher level routine FormInitialGuess().
126*e7a95102SMartin Diehl!
127*e7a95102SMartin Diehl!  Input Parameter:
128*e7a95102SMartin Diehl!  X1 - local vector data
129*e7a95102SMartin Diehl!
130*e7a95102SMartin Diehl!  Output Parameters:
131*e7a95102SMartin Diehl!  x - local vector data
132*e7a95102SMartin Diehl!  ierr - error code
133*e7a95102SMartin Diehl!
134*e7a95102SMartin Diehl!  Notes:
135*e7a95102SMartin Diehl!  This routine uses standard Fortran-style computations over a 2-dim array.
136*e7a95102SMartin Diehl!
137*e7a95102SMartin Diehl  subroutine InitialGuessLocal(solver, X1, ierr)
138*e7a95102SMartin Diehl!  Input/output variables:
139*e7a95102SMartin Diehl    type(ex73f90tmodule_type) solver
140*e7a95102SMartin Diehl    Vec::      X1
141*e7a95102SMartin Diehl    PetscErrorCode ierr
142*e7a95102SMartin Diehl
143*e7a95102SMartin Diehl!  Local variables:
144*e7a95102SMartin Diehl    PetscInt row, i, j, ione, low, high
145*e7a95102SMartin Diehl    PetscReal temp1, temp, hx, hy, v
146*e7a95102SMartin Diehl    PetscReal one
147*e7a95102SMartin Diehl
148*e7a95102SMartin Diehl!  Set parameters
149*e7a95102SMartin Diehl    ione = 1
150*e7a95102SMartin Diehl    ierr = 0
151*e7a95102SMartin Diehl    one = 1.0
152*e7a95102SMartin Diehl    hx = one/(solver%mx - 1)
153*e7a95102SMartin Diehl    hy = one/(solver%my - 1)
154*e7a95102SMartin Diehl    temp1 = solver%lambda/(solver%lambda + one) + one
155*e7a95102SMartin Diehl
156*e7a95102SMartin Diehl    PetscCall(VecGetOwnershipRange(X1, low, high, ierr))
157*e7a95102SMartin Diehl
158*e7a95102SMartin Diehl    do row = low, high - 1
159*e7a95102SMartin Diehl      j = row/solver%mx
160*e7a95102SMartin Diehl      i = mod(row, solver%mx)
161*e7a95102SMartin Diehl      temp = min(j, solver%my - j + 1)*hy
162*e7a95102SMartin Diehl      if (i == 0 .or. j == 0 .or. i == solver%mx - 1 .or. j == solver%my - 1) then
163*e7a95102SMartin Diehl        v = 0.0
164*e7a95102SMartin Diehl      else
165*e7a95102SMartin Diehl        v = temp1*sqrt(min(min(i, solver%mx - i + 1)*hx, temp))
166*e7a95102SMartin Diehl      end if
167*e7a95102SMartin Diehl      PetscCall(VecSetValues(X1, ione, [row], [v], INSERT_VALUES, ierr))
168*e7a95102SMartin Diehl    end do
169*e7a95102SMartin Diehl
170*e7a95102SMartin Diehl  end subroutine InitialGuessLocal
171*e7a95102SMartin Diehl
172*e7a95102SMartin Diehl! ---------------------------------------------------------------------
173*e7a95102SMartin Diehl!
174*e7a95102SMartin Diehl!  FormJacobian - Evaluates Jacobian matrix.
175*e7a95102SMartin Diehl!
176*e7a95102SMartin Diehl!  Input Parameters:
177*e7a95102SMartin Diehl!  dummy     - the SNES context
178*e7a95102SMartin Diehl!  x         - input vector
179*e7a95102SMartin Diehl!  solver    - solver data
180*e7a95102SMartin Diehl!
181*e7a95102SMartin Diehl!  Output Parameters:
182*e7a95102SMartin Diehl!  jac      - Jacobian matrix
183*e7a95102SMartin Diehl!  jac_prec - optionally different matrix used to construct the preconditioner (not used here)
184*e7a95102SMartin Diehl!
185*e7a95102SMartin Diehl  subroutine FormJacobian(dummy, X, jac, jac_prec, solver, ierr)
186*e7a95102SMartin Diehl!  Input/output variables:
187*e7a95102SMartin Diehl    SNES::     dummy
188*e7a95102SMartin Diehl    Vec::      X
189*e7a95102SMartin Diehl    Mat::     jac, jac_prec
190*e7a95102SMartin Diehl    type(ex73f90tmodule_type) solver
191*e7a95102SMartin Diehl    PetscErrorCode ierr
192*e7a95102SMartin Diehl
193*e7a95102SMartin Diehl!  Declarations for use with local arrays:
194*e7a95102SMartin Diehl    Vec::      Xsub(1)
195*e7a95102SMartin Diehl    Mat::     Amat
196*e7a95102SMartin Diehl    PetscInt ione
197*e7a95102SMartin Diehl
198*e7a95102SMartin Diehl    ione = 1
199*e7a95102SMartin Diehl
200*e7a95102SMartin Diehl    PetscCall(DMCompositeGetAccessArray(solver%da, X, ione, PETSC_NULL_INTEGER_ARRAY, Xsub, ierr))
201*e7a95102SMartin Diehl
202*e7a95102SMartin Diehl!     Compute entries for the locally owned part of the Jacobian preconditioner.
203*e7a95102SMartin Diehl    PetscCall(MatCreateSubMatrix(jac_prec, solver%isPhi, solver%isPhi, MAT_INITIAL_MATRIX, Amat, ierr))
204*e7a95102SMartin Diehl
205*e7a95102SMartin Diehl    PetscCall(FormJacobianLocal(Xsub(1), Amat, solver, .true., ierr))
206*e7a95102SMartin Diehl    PetscCall(MatDestroy(Amat, ierr)) ! discard our reference
207*e7a95102SMartin Diehl    PetscCall(DMCompositeRestoreAccessArray(solver%da, X, ione, PETSC_NULL_INTEGER_ARRAY, Xsub, ierr))
208*e7a95102SMartin Diehl
209*e7a95102SMartin Diehl    ! the rest of the matrix is not touched
210*e7a95102SMartin Diehl    PetscCall(MatAssemblyBegin(jac_prec, MAT_FINAL_ASSEMBLY, ierr))
211*e7a95102SMartin Diehl    PetscCall(MatAssemblyEnd(jac_prec, MAT_FINAL_ASSEMBLY, ierr))
212*e7a95102SMartin Diehl    if (jac /= jac_prec) then
213*e7a95102SMartin Diehl      PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY, ierr))
214*e7a95102SMartin Diehl      PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY, ierr))
215*e7a95102SMartin Diehl    end if
216*e7a95102SMartin Diehl
217*e7a95102SMartin Diehl!     Tell the matrix we will never add a new nonzero location to the
218*e7a95102SMartin Diehl!     matrix. If we do it will generate an error.
219*e7a95102SMartin Diehl    PetscCall(MatSetOption(jac_prec, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE, ierr))
220*e7a95102SMartin Diehl
221*e7a95102SMartin Diehl  end subroutine FormJacobian
222*e7a95102SMartin Diehl
223*e7a95102SMartin Diehl! ---------------------------------------------------------------------
224*e7a95102SMartin Diehl!
225*e7a95102SMartin Diehl!  FormJacobianLocal - Computes Jacobian matrix used to compute the preconditioner,
226*e7a95102SMartin Diehl!  called by the higher level routine FormJacobian().
227*e7a95102SMartin Diehl!
228*e7a95102SMartin Diehl!  Input Parameters:
229*e7a95102SMartin Diehl!  x        - local vector data
230*e7a95102SMartin Diehl!
231*e7a95102SMartin Diehl!  Output Parameters:
232*e7a95102SMartin Diehl!  jac - Jacobian matrix used to compute the preconditioner
233*e7a95102SMartin Diehl!  ierr     - error code
234*e7a95102SMartin Diehl!
235*e7a95102SMartin Diehl!  Notes:
236*e7a95102SMartin Diehl!  This routine uses standard Fortran-style computations over a 2-dim array.
237*e7a95102SMartin Diehl!
238*e7a95102SMartin Diehl  subroutine FormJacobianLocal(X1, jac, solver, add_nl_term, ierr)
239*e7a95102SMartin Diehl!  Input/output variables:
240*e7a95102SMartin Diehl    type(ex73f90tmodule_type) solver
241*e7a95102SMartin Diehl    Vec::      X1
242*e7a95102SMartin Diehl    Mat::     jac
243*e7a95102SMartin Diehl    logical add_nl_term
244*e7a95102SMartin Diehl    PetscErrorCode ierr
245*e7a95102SMartin Diehl
246*e7a95102SMartin Diehl!  Local variables:
247*e7a95102SMartin Diehl    PetscInt irow, row(1), col(5), i, j
248*e7a95102SMartin Diehl    PetscInt ione, ifive, low, high, ii
249*e7a95102SMartin Diehl    PetscScalar two, one, hx, hy, hy2inv
250*e7a95102SMartin Diehl    PetscScalar hx2inv, sc, v(5)
251*e7a95102SMartin Diehl    PetscScalar, pointer :: lx_v(:)
252*e7a95102SMartin Diehl
253*e7a95102SMartin Diehl!  Set parameters
254*e7a95102SMartin Diehl    ione = 1
255*e7a95102SMartin Diehl    ifive = 5
256*e7a95102SMartin Diehl    one = 1.0
257*e7a95102SMartin Diehl    two = 2.0
258*e7a95102SMartin Diehl    hx = one/(solver%mx - 1)
259*e7a95102SMartin Diehl    hy = one/(solver%my - 1)
260*e7a95102SMartin Diehl    sc = solver%lambda
261*e7a95102SMartin Diehl    hx2inv = one/(hx*hx)
262*e7a95102SMartin Diehl    hy2inv = one/(hy*hy)
263*e7a95102SMartin Diehl
264*e7a95102SMartin Diehl    PetscCall(VecGetOwnershipRange(X1, low, high, ierr))
265*e7a95102SMartin Diehl    PetscCall(VecGetArrayRead(X1, lx_v, ierr))
266*e7a95102SMartin Diehl
267*e7a95102SMartin Diehl    ii = 0
268*e7a95102SMartin Diehl    do irow = low, high - 1
269*e7a95102SMartin Diehl      j = irow/solver%mx
270*e7a95102SMartin Diehl      i = mod(irow, solver%mx)
271*e7a95102SMartin Diehl      ii = ii + 1            ! one based local index
272*e7a95102SMartin Diehl!     boundary points
273*e7a95102SMartin Diehl      if (i == 0 .or. j == 0 .or. i == solver%mx - 1 .or. j == solver%my - 1) then
274*e7a95102SMartin Diehl        col(1) = irow
275*e7a95102SMartin Diehl        row(1) = irow
276*e7a95102SMartin Diehl        v(1) = one
277*e7a95102SMartin Diehl        PetscCall(MatSetValues(jac, ione, row, ione, col, v, INSERT_VALUES, ierr))
278*e7a95102SMartin Diehl!     interior grid points
279*e7a95102SMartin Diehl      else
280*e7a95102SMartin Diehl        v(1) = -hy2inv
281*e7a95102SMartin Diehl        if (j - 1 == 0) v(1) = 0.0
282*e7a95102SMartin Diehl        v(2) = -hx2inv
283*e7a95102SMartin Diehl        if (i - 1 == 0) v(2) = 0.0
284*e7a95102SMartin Diehl        v(3) = two*(hx2inv + hy2inv)
285*e7a95102SMartin Diehl        if (add_nl_term) v(3) = v(3) - sc*exp(lx_v(ii))
286*e7a95102SMartin Diehl        v(4) = -hx2inv
287*e7a95102SMartin Diehl        if (i + 1 == solver%mx - 1) v(4) = 0.0
288*e7a95102SMartin Diehl        v(5) = -hy2inv
289*e7a95102SMartin Diehl        if (j + 1 == solver%my - 1) v(5) = 0.0
290*e7a95102SMartin Diehl        col(1) = irow - solver%mx
291*e7a95102SMartin Diehl        col(2) = irow - 1
292*e7a95102SMartin Diehl        col(3) = irow
293*e7a95102SMartin Diehl        col(4) = irow + 1
294*e7a95102SMartin Diehl        col(5) = irow + solver%mx
295*e7a95102SMartin Diehl        row(1) = irow
296*e7a95102SMartin Diehl        PetscCall(MatSetValues(jac, ione, row, ifive, col, v, INSERT_VALUES, ierr))
297*e7a95102SMartin Diehl      end if
298*e7a95102SMartin Diehl    end do
299*e7a95102SMartin Diehl
300*e7a95102SMartin Diehl    PetscCall(VecRestoreArrayRead(X1, lx_v, ierr))
301*e7a95102SMartin Diehl
302*e7a95102SMartin Diehl  end subroutine FormJacobianLocal
303*e7a95102SMartin Diehl
304*e7a95102SMartin Diehl! ---------------------------------------------------------------------
305*e7a95102SMartin Diehl!
306*e7a95102SMartin Diehl!  FormFunction - Evaluates nonlinear function, F(x).
307*e7a95102SMartin Diehl!
308*e7a95102SMartin Diehl!  Input Parameters:
309*e7a95102SMartin Diehl!  snes - the SNES context
310*e7a95102SMartin Diehl!  X - input vector
311*e7a95102SMartin Diehl!  dummy - optional user-defined context, as set by SNESSetFunction()
312*e7a95102SMartin Diehl!          (not used here)
313*e7a95102SMartin Diehl!
314*e7a95102SMartin Diehl!  Output Parameter:
315*e7a95102SMartin Diehl!  F - function vector
316*e7a95102SMartin Diehl!
317*e7a95102SMartin Diehl  subroutine FormFunction(snesIn, X, F, solver, ierr)
318*e7a95102SMartin Diehl!  Input/output variables:
319*e7a95102SMartin Diehl    SNES::     snesIn
320*e7a95102SMartin Diehl    Vec::      X, F
321*e7a95102SMartin Diehl    PetscErrorCode ierr
322*e7a95102SMartin Diehl    type(ex73f90tmodule_type) solver
323*e7a95102SMartin Diehl
324*e7a95102SMartin Diehl!  Declarations for use with local arrays:
325*e7a95102SMartin Diehl    Vec::              Xsub(2), Fsub(2)
326*e7a95102SMartin Diehl    PetscInt itwo
327*e7a95102SMartin Diehl
328*e7a95102SMartin Diehl!  Scatter ghost points to local vector, using the 2-step process
329*e7a95102SMartin Diehl!     DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
330*e7a95102SMartin Diehl!  By placing code between these two statements, computations can
331*e7a95102SMartin Diehl!  be done while messages are in transition.
332*e7a95102SMartin Diehl
333*e7a95102SMartin Diehl    itwo = 2
334*e7a95102SMartin Diehl    PetscCall(DMCompositeGetAccessArray(solver%da, X, itwo, PETSC_NULL_INTEGER_ARRAY, Xsub, ierr))
335*e7a95102SMartin Diehl    PetscCall(DMCompositeGetAccessArray(solver%da, F, itwo, PETSC_NULL_INTEGER_ARRAY, Fsub, ierr))
336*e7a95102SMartin Diehl
337*e7a95102SMartin Diehl    PetscCall(FormFunctionNLTerm(Xsub(1), Fsub(1), solver, ierr))
338*e7a95102SMartin Diehl    PetscCall(MatMultAdd(solver%AmatLin, Xsub(1), Fsub(1), Fsub(1), ierr))
339*e7a95102SMartin Diehl
340*e7a95102SMartin Diehl!     do rest of operator (linear)
341*e7a95102SMartin Diehl    PetscCall(MatMult(solver%Cmat, Xsub(1), Fsub(2), ierr))
342*e7a95102SMartin Diehl    PetscCall(MatMultAdd(solver%Bmat, Xsub(2), Fsub(1), Fsub(1), ierr))
343*e7a95102SMartin Diehl    PetscCall(MatMultAdd(solver%Dmat, Xsub(2), Fsub(2), Fsub(2), ierr))
344*e7a95102SMartin Diehl
345*e7a95102SMartin Diehl    PetscCall(DMCompositeRestoreAccessArray(solver%da, X, itwo, PETSC_NULL_INTEGER_ARRAY, Xsub, ierr))
346*e7a95102SMartin Diehl    PetscCall(DMCompositeRestoreAccessArray(solver%da, F, itwo, PETSC_NULL_INTEGER_ARRAY, Fsub, ierr))
347*e7a95102SMartin Diehl  end subroutine formfunction
348*e7a95102SMartin Diehl
349*e7a95102SMartin Diehl! ---------------------------------------------------------------------
350*e7a95102SMartin Diehl!
351*e7a95102SMartin Diehl!  FormFunctionNLTerm - Computes nonlinear function, called by
352*e7a95102SMartin Diehl!  the higher level routine FormFunction().
353*e7a95102SMartin Diehl!
354*e7a95102SMartin Diehl!  Input Parameter:
355*e7a95102SMartin Diehl!  x - local vector data
356*e7a95102SMartin Diehl!
357*e7a95102SMartin Diehl!  Output Parameters:
358*e7a95102SMartin Diehl!  f - local vector data, f(x)
359*e7a95102SMartin Diehl!  ierr - error code
360*e7a95102SMartin Diehl!
361*e7a95102SMartin Diehl!  Notes:
362*e7a95102SMartin Diehl!  This routine uses standard Fortran-style computations over a 2-dim array.
363*e7a95102SMartin Diehl!
364*e7a95102SMartin Diehl  subroutine FormFunctionNLTerm(X1, F1, solver, ierr)
365*e7a95102SMartin Diehl!  Input/output variables:
366*e7a95102SMartin Diehl    type(ex73f90tmodule_type) solver
367*e7a95102SMartin Diehl    Vec::      X1, F1
368*e7a95102SMartin Diehl    PetscErrorCode ierr
369*e7a95102SMartin Diehl!  Local variables:
370*e7a95102SMartin Diehl    PetscScalar sc
371*e7a95102SMartin Diehl    PetscScalar u, v(1)
372*e7a95102SMartin Diehl    PetscInt i, j, low, high, ii, ione, irow, row(1)
373*e7a95102SMartin Diehl    PetscScalar, pointer :: lx_v(:)
374*e7a95102SMartin Diehl
375*e7a95102SMartin Diehl    sc = solver%lambda
376*e7a95102SMartin Diehl    ione = 1
377*e7a95102SMartin Diehl
378*e7a95102SMartin Diehl    PetscCall(VecGetArrayRead(X1, lx_v, ierr))
379*e7a95102SMartin Diehl    PetscCall(VecGetOwnershipRange(X1, low, high, ierr))
380*e7a95102SMartin Diehl
381*e7a95102SMartin Diehl!     Compute function over the locally owned part of the grid
382*e7a95102SMartin Diehl    ii = 0
383*e7a95102SMartin Diehl    do irow = low, high - 1
384*e7a95102SMartin Diehl      j = irow/solver%mx
385*e7a95102SMartin Diehl      i = mod(irow, solver%mx)
386*e7a95102SMartin Diehl      ii = ii + 1            ! one based local index
387*e7a95102SMartin Diehl      row(1) = irow
388*e7a95102SMartin Diehl      if (i == 0 .or. j == 0 .or. i == solver%mx - 1 .or. j == solver%my - 1) then
389*e7a95102SMartin Diehl        v(1) = 0.0
390*e7a95102SMartin Diehl      else
391*e7a95102SMartin Diehl        u = lx_v(ii)
392*e7a95102SMartin Diehl        v(1) = -sc*exp(u)
393*e7a95102SMartin Diehl      end if
394*e7a95102SMartin Diehl      PetscCall(VecSetValues(F1, ione, row, v, INSERT_VALUES, ierr))
395*e7a95102SMartin Diehl    end do
396*e7a95102SMartin Diehl
397*e7a95102SMartin Diehl    PetscCall(VecRestoreArrayRead(X1, lx_v, ierr))
398*e7a95102SMartin Diehl
399*e7a95102SMartin Diehl    PetscCall(VecAssemblyBegin(F1, ierr))
400*e7a95102SMartin Diehl    PetscCall(VecAssemblyEnd(F1, ierr))
401*e7a95102SMartin Diehl
402*e7a95102SMartin Diehl    ierr = 0
403*e7a95102SMartin Diehl  end subroutine FormFunctionNLTerm
404*e7a95102SMartin Diehl
405*e7a95102SMartin Diehlend module ex73f90tmodule
406*e7a95102SMartin Diehl
407c4762a1bSJed Brownprogram main
408c4762a1bSJed Brown  use petscsnes
40977d968b7SBarry Smith  use ex73f90tmodule
410c4762a1bSJed Brown  implicit none
411c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
412c4762a1bSJed Brown!                   Variable declarations
413c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
414c4762a1bSJed Brown!
415c4762a1bSJed Brown!  Variables:
416c4762a1bSJed Brown!     mysnes      - nonlinear solver
417c4762a1bSJed Brown!     x, r        - solution, residual vectors
418c4762a1bSJed Brown!     J           - Jacobian matrix
419c4762a1bSJed Brown!     its         - iterations for convergence
420c4762a1bSJed Brown!     Nx, Ny      - number of preocessors in x- and y- directions
421c4762a1bSJed Brown!
422c4762a1bSJed Brown  SNES::       mysnes
423c4762a1bSJed Brown  Vec::        x, r, x2, x1, x1loc, x2loc
424c4762a1bSJed Brown  Mat::       Amat, Bmat, Cmat, Dmat, KKTMat, matArray(4)
425c4762a1bSJed Brown!      Mat::       tmat
426c4762a1bSJed Brown  DM::       daphi, dalam
427ce78bad3SBarry Smith  IS, pointer ::        isglobal(:)
428c4762a1bSJed Brown  PetscErrorCode ierr
429c4762a1bSJed Brown  PetscInt its, N1, N2, i, j, irow, row(1)
430c4762a1bSJed Brown  PetscInt col(1), low, high, lamlow, lamhigh
431c4762a1bSJed Brown  PetscBool flg
432c4762a1bSJed Brown  PetscInt ione, nfour, itwo, nloc, nloclam
433c4762a1bSJed Brown  PetscReal lambda_max, lambda_min
43477d968b7SBarry Smith  type(ex73f90tmodule_type) solver
435c4762a1bSJed Brown  PetscScalar bval(1), cval(1), one
436c00ad2bcSBarry Smith  PetscBool useobjective
437c4762a1bSJed Brown
438c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
439c4762a1bSJed Brown!  Initialize program
440c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
441d8606c27SBarry Smith  PetscCallA(PetscInitialize(ierr))
442d8606c27SBarry Smith  PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, solver%rank, ierr))
443c4762a1bSJed Brown
444c4762a1bSJed Brown!  Initialize problem parameters
445c4762a1bSJed Brown  lambda_max = 6.81_PETSC_REAL_KIND
446c4762a1bSJed Brown  lambda_min = 0.0
447c4762a1bSJed Brown  solver%lambda = 6.0
448c4762a1bSJed Brown  ione = 1
449c4762a1bSJed Brown  nfour = 4
450c4762a1bSJed Brown  itwo = 2
451c00ad2bcSBarry Smith  useobjective = PETSC_FALSE
452d8606c27SBarry Smith  PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-par', solver%lambda, flg, ierr))
4534820e4eaSBarry Smith  PetscCheckA(solver%lambda <= lambda_max .and. solver%lambda >= lambda_min, PETSC_COMM_SELF, PETSC_ERR_USER, 'Lambda provided with -par is out of range')
454c00ad2bcSBarry Smith  PetscCallA(PetscOptionsGetBool(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-objective', useobjective, flg, ierr))
455c4762a1bSJed Brown
456c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
457c4762a1bSJed Brown!  Create vector data structures; set function evaluation routine
458c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
459c4762a1bSJed Brown
460c4762a1bSJed Brown!     just get size
4615d83a8b1SBarry 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, daphi, ierr))
462d8606c27SBarry Smith  PetscCallA(DMSetFromOptions(daphi, ierr))
463d8606c27SBarry Smith  PetscCallA(DMSetUp(daphi, ierr))
464ce78bad3SBarry Smith  PetscCallA(DMDAGetInfo(daphi, PETSC_NULL_INTEGER, solver%mx, solver%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))
465c4762a1bSJed Brown  N1 = solver%my*solver%mx
466c4762a1bSJed Brown  N2 = solver%my
467c4762a1bSJed Brown  flg = .false.
468d8606c27SBarry Smith  PetscCallA(PetscOptionsGetBool(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-no_constraints', flg, flg, ierr))
469c4762a1bSJed Brown  if (flg) then
470c4762a1bSJed Brown    N2 = 0
471c4762a1bSJed Brown  end if
472c4762a1bSJed Brown
473d8606c27SBarry Smith  PetscCallA(DMDestroy(daphi, ierr))
474c4762a1bSJed Brown
475c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
476c4762a1bSJed Brown!  Create matrix data structure; set Jacobian evaluation routine
477c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
478d8606c27SBarry Smith  PetscCallA(DMShellCreate(PETSC_COMM_WORLD, daphi, ierr))
479d8606c27SBarry Smith  PetscCallA(DMSetOptionsPrefix(daphi, 'phi_', ierr))
480d8606c27SBarry Smith  PetscCallA(DMSetFromOptions(daphi, ierr))
481c4762a1bSJed Brown
482d8606c27SBarry Smith  PetscCallA(VecCreate(PETSC_COMM_WORLD, x1, ierr))
483d8606c27SBarry Smith  PetscCallA(VecSetSizes(x1, PETSC_DECIDE, N1, ierr))
484d8606c27SBarry Smith  PetscCallA(VecSetFromOptions(x1, ierr))
485c4762a1bSJed Brown
486d8606c27SBarry Smith  PetscCallA(VecGetOwnershipRange(x1, low, high, ierr))
487c4762a1bSJed Brown  nloc = high - low
488c4762a1bSJed Brown
489d8606c27SBarry Smith  PetscCallA(MatCreate(PETSC_COMM_WORLD, Amat, ierr))
490d8606c27SBarry Smith  PetscCallA(MatSetSizes(Amat, PETSC_DECIDE, PETSC_DECIDE, N1, N1, ierr))
491d8606c27SBarry Smith  PetscCallA(MatSetUp(Amat, ierr))
492c4762a1bSJed Brown
493d8606c27SBarry Smith  PetscCallA(MatCreate(PETSC_COMM_WORLD, solver%AmatLin, ierr))
494d8606c27SBarry Smith  PetscCallA(MatSetSizes(solver%AmatLin, PETSC_DECIDE, PETSC_DECIDE, N1, N1, ierr))
495d8606c27SBarry Smith  PetscCallA(MatSetUp(solver%AmatLin, ierr))
496c4762a1bSJed Brown
497d8606c27SBarry Smith  PetscCallA(FormJacobianLocal(x1, solver%AmatLin, solver, .false., ierr))
498d8606c27SBarry Smith  PetscCallA(MatAssemblyBegin(solver%AmatLin, MAT_FINAL_ASSEMBLY, ierr))
499d8606c27SBarry Smith  PetscCallA(MatAssemblyEnd(solver%AmatLin, MAT_FINAL_ASSEMBLY, ierr))
500c4762a1bSJed Brown
501d8606c27SBarry Smith  PetscCallA(DMShellSetGlobalVector(daphi, x1, ierr))
502d8606c27SBarry Smith  PetscCallA(DMShellSetMatrix(daphi, Amat, ierr))
503c4762a1bSJed Brown
504d8606c27SBarry Smith  PetscCallA(VecCreate(PETSC_COMM_SELF, x1loc, ierr))
505d8606c27SBarry Smith  PetscCallA(VecSetSizes(x1loc, nloc, nloc, ierr))
506d8606c27SBarry Smith  PetscCallA(VecSetFromOptions(x1loc, ierr))
507d8606c27SBarry Smith  PetscCallA(DMShellSetLocalVector(daphi, x1loc, ierr))
508c4762a1bSJed Brown
509c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
510c4762a1bSJed Brown!  Create B, C, & D matrices
511c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
512d8606c27SBarry Smith  PetscCallA(MatCreate(PETSC_COMM_WORLD, Cmat, ierr))
513d8606c27SBarry Smith  PetscCallA(MatSetSizes(Cmat, PETSC_DECIDE, PETSC_DECIDE, N2, N1, ierr))
514d8606c27SBarry Smith  PetscCallA(MatSetUp(Cmat, ierr))
515c4762a1bSJed Brown!      create data for C and B
516d8606c27SBarry Smith  PetscCallA(MatCreate(PETSC_COMM_WORLD, Bmat, ierr))
517d8606c27SBarry Smith  PetscCallA(MatSetSizes(Bmat, PETSC_DECIDE, PETSC_DECIDE, N1, N2, ierr))
518d8606c27SBarry Smith  PetscCallA(MatSetUp(Bmat, ierr))
519c4762a1bSJed Brown!     create data for D
520d8606c27SBarry Smith  PetscCallA(MatCreate(PETSC_COMM_WORLD, Dmat, ierr))
521d8606c27SBarry Smith  PetscCallA(MatSetSizes(Dmat, PETSC_DECIDE, PETSC_DECIDE, N2, N2, ierr))
522d8606c27SBarry Smith  PetscCallA(MatSetUp(Dmat, ierr))
523c4762a1bSJed Brown
524d8606c27SBarry Smith  PetscCallA(VecCreate(PETSC_COMM_WORLD, x2, ierr))
525d8606c27SBarry Smith  PetscCallA(VecSetSizes(x2, PETSC_DECIDE, N2, ierr))
526d8606c27SBarry Smith  PetscCallA(VecSetFromOptions(x2, ierr))
527c4762a1bSJed Brown
528d8606c27SBarry Smith  PetscCallA(VecGetOwnershipRange(x2, lamlow, lamhigh, ierr))
529c4762a1bSJed Brown  nloclam = lamhigh - lamlow
530c4762a1bSJed Brown
531c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
532c4762a1bSJed Brown!  Set fake B and C
533c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
534c4762a1bSJed Brown  one = 1.0
5354820e4eaSBarry Smith  if (N2 > 0) then
536c4762a1bSJed Brown    bval(1) = -one/(solver%mx - 2)
537c4762a1bSJed Brown!     cval = -one/(solver%my*solver%mx)
538c4762a1bSJed Brown    cval(1) = -one
539ceeae6b5SMartin Diehl    do irow = low, high - 1
540c4762a1bSJed Brown      j = irow/solver%mx   ! row in domain
541c4762a1bSJed Brown      i = mod(irow, solver%mx)
542c4762a1bSJed Brown      row(1) = irow
543c4762a1bSJed Brown      col(1) = j
5444820e4eaSBarry Smith      if (i == 0 .or. j == 0 .or. i == solver%mx - 1 .or. j == solver%my - 1) then
545c4762a1bSJed Brown        !     no op
546c4762a1bSJed Brown      else
547d8606c27SBarry Smith        PetscCallA(MatSetValues(Bmat, ione, row, ione, col, bval, INSERT_VALUES, ierr))
548c4762a1bSJed Brown      end if
549c4762a1bSJed Brown      row(1) = j
550d8606c27SBarry Smith      PetscCallA(MatSetValues(Cmat, ione, row, ione, row, cval, INSERT_VALUES, ierr))
551ceeae6b5SMartin Diehl    end do
552c4762a1bSJed Brown  end if
553d8606c27SBarry Smith  PetscCallA(MatAssemblyBegin(Bmat, MAT_FINAL_ASSEMBLY, ierr))
554d8606c27SBarry Smith  PetscCallA(MatAssemblyEnd(Bmat, MAT_FINAL_ASSEMBLY, ierr))
555d8606c27SBarry Smith  PetscCallA(MatAssemblyBegin(Cmat, MAT_FINAL_ASSEMBLY, ierr))
556d8606c27SBarry Smith  PetscCallA(MatAssemblyEnd(Cmat, MAT_FINAL_ASSEMBLY, ierr))
557c4762a1bSJed Brown
558c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
559da81f932SPierre Jolivet!  Set D (identity)
560c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
561ceeae6b5SMartin Diehl  do j = lamlow, lamhigh - 1
562c4762a1bSJed Brown    row(1) = j
563c4762a1bSJed Brown    cval(1) = one
564d8606c27SBarry Smith    PetscCallA(MatSetValues(Dmat, ione, row, ione, row, cval, INSERT_VALUES, ierr))
565ceeae6b5SMartin Diehl  end do
566d8606c27SBarry Smith  PetscCallA(MatAssemblyBegin(Dmat, MAT_FINAL_ASSEMBLY, ierr))
567d8606c27SBarry Smith  PetscCallA(MatAssemblyEnd(Dmat, MAT_FINAL_ASSEMBLY, ierr))
568c4762a1bSJed Brown
569c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
570c4762a1bSJed Brown!  DM for lambda (dalam) : temp driver for A block, setup A block solver data
571c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
572d8606c27SBarry Smith  PetscCallA(DMShellCreate(PETSC_COMM_WORLD, dalam, ierr))
573d8606c27SBarry Smith  PetscCallA(DMShellSetGlobalVector(dalam, x2, ierr))
574d8606c27SBarry Smith  PetscCallA(DMShellSetMatrix(dalam, Dmat, ierr))
575c4762a1bSJed Brown
576d8606c27SBarry Smith  PetscCallA(VecCreate(PETSC_COMM_SELF, x2loc, ierr))
577d8606c27SBarry Smith  PetscCallA(VecSetSizes(x2loc, nloclam, nloclam, ierr))
578d8606c27SBarry Smith  PetscCallA(VecSetFromOptions(x2loc, ierr))
579d8606c27SBarry Smith  PetscCallA(DMShellSetLocalVector(dalam, x2loc, ierr))
580c4762a1bSJed Brown
581d8606c27SBarry Smith  PetscCallA(DMSetOptionsPrefix(dalam, 'lambda_', ierr))
582d8606c27SBarry Smith  PetscCallA(DMSetFromOptions(dalam, ierr))
583c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
584c4762a1bSJed Brown!  Create field split DA
585c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
586d8606c27SBarry Smith  PetscCallA(DMCompositeCreate(PETSC_COMM_WORLD, solver%da, ierr))
587d8606c27SBarry Smith  PetscCallA(DMCompositeAddDM(solver%da, daphi, ierr))
588d8606c27SBarry Smith  PetscCallA(DMCompositeAddDM(solver%da, dalam, ierr))
589d8606c27SBarry Smith  PetscCallA(DMSetFromOptions(solver%da, ierr))
590d8606c27SBarry Smith  PetscCallA(DMSetUp(solver%da, ierr))
591d8606c27SBarry Smith  PetscCallA(DMCompositeGetGlobalISs(solver%da, isglobal, ierr))
592c4762a1bSJed Brown  solver%isPhi = isglobal(1)
593ce78bad3SBarry Smith  PetscCallA(PetscObjectReference(solver%isPhi, ierr))
594c4762a1bSJed Brown  solver%isLambda = isglobal(2)
595ce78bad3SBarry Smith  PetscCallA(PetscObjectReference(solver%isLambda, ierr))
596c4762a1bSJed Brown
597c4762a1bSJed Brown!     cache matrices
598c4762a1bSJed Brown  solver%Amat = Amat
599c4762a1bSJed Brown  solver%Bmat = Bmat
600c4762a1bSJed Brown  solver%Cmat = Cmat
601c4762a1bSJed Brown  solver%Dmat = Dmat
602c4762a1bSJed Brown
603c4762a1bSJed Brown  matArray(1) = Amat
604c4762a1bSJed Brown  matArray(2) = Bmat
605c4762a1bSJed Brown  matArray(3) = Cmat
606c4762a1bSJed Brown  matArray(4) = Dmat
607c4762a1bSJed Brown
608d8606c27SBarry Smith  PetscCallA(MatCreateNest(PETSC_COMM_WORLD, itwo, isglobal, itwo, isglobal, matArray, KKTmat, ierr))
609ce78bad3SBarry Smith  PetscCallA(DMCompositeRestoreGlobalISs(solver%da, isglobal, ierr))
610d8606c27SBarry Smith  PetscCallA(MatSetFromOptions(KKTmat, ierr))
611c4762a1bSJed Brown
612c4762a1bSJed Brown!  Extract global and local vectors from DMDA; then duplicate for remaining
613c4762a1bSJed Brown!     vectors that are the same types
614d8606c27SBarry Smith  PetscCallA(MatCreateVecs(KKTmat, x, PETSC_NULL_VEC, ierr))
615d8606c27SBarry Smith  PetscCallA(VecDuplicate(x, r, ierr))
616c4762a1bSJed Brown
617d8606c27SBarry Smith  PetscCallA(SNESCreate(PETSC_COMM_WORLD, mysnes, ierr))
618c4762a1bSJed Brown
619d8606c27SBarry Smith  PetscCallA(SNESSetDM(mysnes, solver%da, ierr))
620c4762a1bSJed Brown
621d8606c27SBarry Smith  PetscCallA(SNESSetApplicationContext(mysnes, solver, ierr))
622c4762a1bSJed Brown
623d8606c27SBarry Smith  PetscCallA(SNESSetDM(mysnes, solver%da, ierr))
624c4762a1bSJed Brown
625c4762a1bSJed Brown!  Set function evaluation routine and vector
626d8606c27SBarry Smith  PetscCallA(SNESSetFunction(mysnes, r, FormFunction, solver, ierr))
627c00ad2bcSBarry Smith  if (useobjective .eqv. PETSC_TRUE) then
628c00ad2bcSBarry Smith    PetscCallA(SNESSetObjective(mysnes, MyObjective, solver, ierr))
629c00ad2bcSBarry Smith  end if
630d8606c27SBarry Smith  PetscCallA(SNESSetJacobian(mysnes, KKTmat, KKTmat, FormJacobian, solver, ierr))
631c4762a1bSJed Brown
632c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
633c4762a1bSJed Brown!  Customize nonlinear solver; set runtime options
634c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
635c4762a1bSJed Brown!  Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
636d8606c27SBarry Smith  PetscCallA(SNESSetFromOptions(mysnes, ierr))
637c4762a1bSJed Brown
638c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
639c4762a1bSJed Brown!  Evaluate initial guess; then solve nonlinear system.
640c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
641c4762a1bSJed Brown!  Note: The user should initialize the vector, x, with the initial guess
642c4762a1bSJed Brown!  for the nonlinear solver prior to calling SNESSolve().  In particular,
643c4762a1bSJed Brown!  to employ an initial guess of zero, the user should explicitly set
644c4762a1bSJed Brown!  this vector to zero by calling VecSet().
645c4762a1bSJed Brown
646d8606c27SBarry Smith  PetscCallA(FormInitialGuess(mysnes, x, ierr))
647d8606c27SBarry Smith  PetscCallA(SNESSolve(mysnes, PETSC_NULL_VEC, x, ierr))
648d8606c27SBarry Smith  PetscCallA(SNESGetIterationNumber(mysnes, its, ierr))
6494820e4eaSBarry Smith  if (solver%rank == 0) then
650c4762a1bSJed Brown    write (6, 100) its
651c4762a1bSJed Brown  end if
652c4762a1bSJed Brown100 format('Number of SNES iterations = ', i5)
653c4762a1bSJed Brown
654c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
655c4762a1bSJed Brown!  Free work space.  All PETSc objects should be destroyed when they
656c4762a1bSJed Brown!  are no longer needed.
657c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
658d8606c27SBarry Smith  PetscCallA(MatDestroy(KKTmat, ierr))
659d8606c27SBarry Smith  PetscCallA(MatDestroy(Amat, ierr))
660d8606c27SBarry Smith  PetscCallA(MatDestroy(Dmat, ierr))
661d8606c27SBarry Smith  PetscCallA(MatDestroy(Bmat, ierr))
662d8606c27SBarry Smith  PetscCallA(MatDestroy(Cmat, ierr))
663d8606c27SBarry Smith  PetscCallA(MatDestroy(solver%AmatLin, ierr))
664d8606c27SBarry Smith  PetscCallA(ISDestroy(solver%isPhi, ierr))
665d8606c27SBarry Smith  PetscCallA(ISDestroy(solver%isLambda, ierr))
666d8606c27SBarry Smith  PetscCallA(VecDestroy(x, ierr))
667d8606c27SBarry Smith  PetscCallA(VecDestroy(x2, ierr))
668d8606c27SBarry Smith  PetscCallA(VecDestroy(x1, ierr))
669d8606c27SBarry Smith  PetscCallA(VecDestroy(x1loc, ierr))
670d8606c27SBarry Smith  PetscCallA(VecDestroy(x2loc, ierr))
671d8606c27SBarry Smith  PetscCallA(VecDestroy(r, ierr))
672d8606c27SBarry Smith  PetscCallA(SNESDestroy(mysnes, ierr))
673d8606c27SBarry Smith  PetscCallA(DMDestroy(solver%da, ierr))
674d8606c27SBarry Smith  PetscCallA(DMDestroy(daphi, ierr))
675d8606c27SBarry Smith  PetscCallA(DMDestroy(dalam, ierr))
676c4762a1bSJed Brown
677d8606c27SBarry Smith  PetscCallA(PetscFinalize(ierr))
678c4762a1bSJed Brownend
679c4762a1bSJed Brown
680c4762a1bSJed Brown!/*TEST
681c4762a1bSJed Brown!
682c4762a1bSJed Brown!   build:
6839e489221SSatish Balay!      requires: !single !complex
684c4762a1bSJed Brown!
685c4762a1bSJed Brown!   test:
686c4762a1bSJed Brown!      nsize: 4
68773f7197eSJed Brown!      args: -par 5.0 -da_grid_x 10 -da_grid_y 10 -snes_monitor_short -snes_linesearch_type basic -snes_converged_reason -ksp_type fgmres -ksp_norm_type unpreconditioned -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type upper -ksp_monitor_short -fieldsplit_lambda_ksp_type preonly -fieldsplit_lambda_pc_type jacobi -fieldsplit_phi_pc_type gamg -fieldsplit_phi_pc_gamg_esteig_ksp_type cg -fieldsplit_phi_pc_gamg_esteig_ksp_max_it 10 -fieldsplit_phi_pc_gamg_agg_nsmooths 1 -fieldsplit_phi_pc_gamg_threshold 0.
688c4762a1bSJed Brown!
689c00ad2bcSBarry Smith!   test:
690a99ef635SJonas Heinzmann!      args: -snes_linesearch_type {{secant cp}separate output} -objective {{false true}shared output}
691c00ad2bcSBarry Smith!
692c00ad2bcSBarry Smith!   test:
693c00ad2bcSBarry Smith!      args: -snes_linesearch_type bt -objective {{false true}separate output}
694c00ad2bcSBarry Smith!
695c4762a1bSJed Brown!TEST*/
696