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