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