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