xref: /petsc/src/snes/tutorials/ex73f90t.F90 (revision 7addb90f52a7936ba144cdab1bb2cc37152af090)
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
42ce78bad3SBarry Smith#include <petsc/finclude/petscdmda.h>
43ce78bad3SBarry Smith#include <petsc/finclude/petscdmcomposite.h>
44c4762a1bSJed Brown#include <petsc/finclude/petscmat.h>
45ce78bad3SBarry Smith      use petscdmda
46ce78bad3SBarry Smith      use petscdmcomposite
47c4762a1bSJed Brown      use petscmat
4877d968b7SBarry Smith      type ex73f90tmodule_type
49c4762a1bSJed Brown        DM::da
50c4762a1bSJed Brown!     temp A block stuff
51c4762a1bSJed Brown        PetscInt mx,my
52c4762a1bSJed Brown        PetscMPIInt rank
53c4762a1bSJed Brown        PetscReal lambda
54c4762a1bSJed Brown!     Mats
55c4762a1bSJed Brown        Mat::Amat,AmatLin,Bmat,CMat,Dmat
56c4762a1bSJed Brown        IS::isPhi,isLambda
5777d968b7SBarry Smith      end type ex73f90tmodule_type
58c4762a1bSJed Brown
5977d968b7SBarry Smith      end module ex73f90tmodule
60c4762a1bSJed Brown
6177d968b7SBarry Smith      module ex73f90tmodule_interfaces
6277d968b7SBarry Smith        use ex73f90tmodule
63c4762a1bSJed Brown
64c4762a1bSJed Brown      Interface SNESSetApplicationContext
65c4762a1bSJed Brown        Subroutine SNESSetApplicationContext(snesIn,ctx,ierr)
66c4762a1bSJed Brown#include <petsc/finclude/petscsnes.h>
67c4762a1bSJed Brown        use petscsnes
6877d968b7SBarry Smith        use ex73f90tmodule
69c4762a1bSJed Brown          SNES::    snesIn
7077d968b7SBarry Smith          type(ex73f90tmodule_type) ctx
71c4762a1bSJed Brown          PetscErrorCode ierr
72c4762a1bSJed Brown        End Subroutine
73c4762a1bSJed Brown      End Interface SNESSetApplicationContext
74c4762a1bSJed Brown
75c4762a1bSJed Brown      Interface SNESGetApplicationContext
76c4762a1bSJed Brown        Subroutine SNESGetApplicationContext(snesIn,ctx,ierr)
77c4762a1bSJed Brown#include <petsc/finclude/petscsnes.h>
78c4762a1bSJed Brown        use petscsnes
7977d968b7SBarry Smith        use ex73f90tmodule
80c4762a1bSJed Brown          SNES::     snesIn
8177d968b7SBarry Smith          type(ex73f90tmodule_type), pointer :: ctx
82c4762a1bSJed Brown          PetscErrorCode ierr
83c4762a1bSJed Brown        End Subroutine
84c4762a1bSJed Brown      End Interface SNESGetApplicationContext
8577d968b7SBarry Smith      end module ex73f90tmodule_interfaces
86c4762a1bSJed Brown
87c00ad2bcSBarry Smith      subroutine MyObjective(snes, x, result, ctx, ierr )
88c00ad2bcSBarry Smith#include <petsc/finclude/petsc.h>
89c00ad2bcSBarry Smith        use petsc
90c00ad2bcSBarry Smith        implicit none
91c00ad2bcSBarry Smith        PetscInt ctx
92c00ad2bcSBarry Smith        Vec            x, f
93c00ad2bcSBarry Smith        SNES           snes
94c00ad2bcSBarry Smith        PetscErrorCode ierr
95c00ad2bcSBarry Smith        PetscScalar    result
96c00ad2bcSBarry Smith        PetscReal      fnorm
97c00ad2bcSBarry Smith
98c00ad2bcSBarry Smith        PetscCall(VecDuplicate(x,f,ierr))
99c00ad2bcSBarry Smith        PetscCall(SNESComputeFunction(snes,x,f,ierr))
100c00ad2bcSBarry Smith        PetscCall(VecNorm(f,NORM_2,fnorm,ierr))
101c00ad2bcSBarry Smith        result = .5*fnorm*fnorm
102c00ad2bcSBarry Smith        PetscCall(VecDestroy(f,ierr))
103c00ad2bcSBarry Smith      end subroutine MyObjective
104c00ad2bcSBarry Smith
105c4762a1bSJed Brown      program main
106ce78bad3SBarry Smith#include <petsc/finclude/petscdmda.h>
107c4762a1bSJed Brown#include <petsc/finclude/petscsnes.h>
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
128ce78bad3SBarry Smith      IS, pointer ::        isglobal(:)
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
137c00ad2bcSBarry Smith      PetscBool        useobjective
138c4762a1bSJed Brown
139c4762a1bSJed Brown!  Note: Any user-defined Fortran routines (such as FormJacobian)
140c4762a1bSJed Brown!  MUST be declared as external.
141c00ad2bcSBarry 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
156c00ad2bcSBarry 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')
159c00ad2bcSBarry 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
1665d83a8b1SBarry 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))
167d8606c27SBarry Smith      PetscCallA(DMSetFromOptions(daphi,ierr))
168d8606c27SBarry Smith      PetscCallA(DMSetUp(daphi,ierr))
169ce78bad3SBarry 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))
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)
298ce78bad3SBarry Smith      PetscCallA(PetscObjectReference(solver%isPhi,ierr))
299c4762a1bSJed Brown      solver%isLambda = isglobal(2)
300ce78bad3SBarry Smith      PetscCallA(PetscObjectReference(solver%isLambda,ierr))
301c4762a1bSJed Brown
302c4762a1bSJed Brown!     cache matrices
303c4762a1bSJed Brown      solver%Amat = Amat
304c4762a1bSJed Brown      solver%Bmat = Bmat
305c4762a1bSJed Brown      solver%Cmat = Cmat
306c4762a1bSJed Brown      solver%Dmat = Dmat
307c4762a1bSJed Brown
308c4762a1bSJed Brown      matArray(1) = Amat
309c4762a1bSJed Brown      matArray(2) = Bmat
310c4762a1bSJed Brown      matArray(3) = Cmat
311c4762a1bSJed Brown      matArray(4) = Dmat
312c4762a1bSJed Brown
313d8606c27SBarry Smith      PetscCallA(MatCreateNest(PETSC_COMM_WORLD,itwo,isglobal,itwo,isglobal,matArray,KKTmat,ierr))
314ce78bad3SBarry Smith      PetscCallA(DMCompositeRestoreGlobalISs(solver%da,isglobal,ierr))
315d8606c27SBarry Smith      PetscCallA(MatSetFromOptions(KKTmat,ierr))
316c4762a1bSJed Brown
317c4762a1bSJed Brown!  Extract global and local vectors from DMDA; then duplicate for remaining
318c4762a1bSJed Brown!     vectors that are the same types
319d8606c27SBarry Smith      PetscCallA(MatCreateVecs(KKTmat,x,PETSC_NULL_VEC,ierr))
320d8606c27SBarry Smith      PetscCallA(VecDuplicate(x,r,ierr))
321c4762a1bSJed Brown
322d8606c27SBarry Smith      PetscCallA(SNESCreate(PETSC_COMM_WORLD,mysnes,ierr))
323c4762a1bSJed Brown
324d8606c27SBarry Smith      PetscCallA(SNESSetDM(mysnes,solver%da,ierr))
325c4762a1bSJed Brown
326d8606c27SBarry Smith      PetscCallA(SNESSetApplicationContext(mysnes,solver,ierr))
327c4762a1bSJed Brown
328d8606c27SBarry Smith      PetscCallA(SNESSetDM(mysnes,solver%da,ierr))
329c4762a1bSJed Brown
330c4762a1bSJed Brown!  Set function evaluation routine and vector
331d8606c27SBarry Smith      PetscCallA(SNESSetFunction(mysnes, r, FormFunction, solver,ierr))
332c00ad2bcSBarry Smith      if (useobjective .eqv. PETSC_TRUE) then
333c00ad2bcSBarry Smith         PetscCallA(SNESSetObjective(mysnes, MyObjective, solver, ierr))
334c00ad2bcSBarry Smith      endif
335d8606c27SBarry Smith      PetscCallA(SNESSetJacobian(mysnes,KKTmat,KKTmat,FormJacobian,solver,ierr))
336c4762a1bSJed Brown
337c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
338c4762a1bSJed Brown!  Customize nonlinear solver; set runtime options
339c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
340c4762a1bSJed Brown!  Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
341d8606c27SBarry Smith      PetscCallA(SNESSetFromOptions(mysnes,ierr))
342c4762a1bSJed Brown
343c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
344c4762a1bSJed Brown!  Evaluate initial guess; then solve nonlinear system.
345c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
346c4762a1bSJed Brown!  Note: The user should initialize the vector, x, with the initial guess
347c4762a1bSJed Brown!  for the nonlinear solver prior to calling SNESSolve().  In particular,
348c4762a1bSJed Brown!  to employ an initial guess of zero, the user should explicitly set
349c4762a1bSJed Brown!  this vector to zero by calling VecSet().
350c4762a1bSJed Brown
351d8606c27SBarry Smith      PetscCallA(FormInitialGuess(mysnes,x,ierr))
352d8606c27SBarry Smith      PetscCallA(SNESSolve(mysnes,PETSC_NULL_VEC,x,ierr))
353d8606c27SBarry Smith      PetscCallA(SNESGetIterationNumber(mysnes,its,ierr))
354c4762a1bSJed Brown      if (solver%rank .eq. 0) then
355c4762a1bSJed Brown         write(6,100) its
356c4762a1bSJed Brown      endif
357c4762a1bSJed Brown  100 format('Number of SNES iterations = ',i5)
358c4762a1bSJed Brown
359c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
360c4762a1bSJed Brown!  Free work space.  All PETSc objects should be destroyed when they
361c4762a1bSJed Brown!  are no longer needed.
362c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
363d8606c27SBarry Smith      PetscCallA(MatDestroy(KKTmat,ierr))
364d8606c27SBarry Smith      PetscCallA(MatDestroy(Amat,ierr))
365d8606c27SBarry Smith      PetscCallA(MatDestroy(Dmat,ierr))
366d8606c27SBarry Smith      PetscCallA(MatDestroy(Bmat,ierr))
367d8606c27SBarry Smith      PetscCallA(MatDestroy(Cmat,ierr))
368d8606c27SBarry Smith      PetscCallA(MatDestroy(solver%AmatLin,ierr))
369d8606c27SBarry Smith      PetscCallA(ISDestroy(solver%isPhi,ierr))
370d8606c27SBarry Smith      PetscCallA(ISDestroy(solver%isLambda,ierr))
371d8606c27SBarry Smith      PetscCallA(VecDestroy(x,ierr))
372d8606c27SBarry Smith      PetscCallA(VecDestroy(x2,ierr))
373d8606c27SBarry Smith      PetscCallA(VecDestroy(x1,ierr))
374d8606c27SBarry Smith      PetscCallA(VecDestroy(x1loc,ierr))
375d8606c27SBarry Smith      PetscCallA(VecDestroy(x2loc,ierr))
376d8606c27SBarry Smith      PetscCallA(VecDestroy(r,ierr))
377d8606c27SBarry Smith      PetscCallA(SNESDestroy(mysnes,ierr))
378d8606c27SBarry Smith      PetscCallA(DMDestroy(solver%da,ierr))
379d8606c27SBarry Smith      PetscCallA(DMDestroy(daphi,ierr))
380d8606c27SBarry Smith      PetscCallA(DMDestroy(dalam,ierr))
381c4762a1bSJed Brown
382d8606c27SBarry Smith      PetscCallA(PetscFinalize(ierr))
383c4762a1bSJed Brown      end
384c4762a1bSJed Brown
385c4762a1bSJed Brown! ---------------------------------------------------------------------
386c4762a1bSJed Brown!
387c4762a1bSJed Brown!  FormInitialGuess - Forms initial approximation.
388c4762a1bSJed Brown!
389c4762a1bSJed Brown!  Input Parameters:
390c4762a1bSJed Brown!  X - vector
391c4762a1bSJed Brown!
392c4762a1bSJed Brown!  Output Parameter:
393c4762a1bSJed Brown!  X - vector
394c4762a1bSJed Brown!
395c4762a1bSJed Brown!  Notes:
396c4762a1bSJed Brown!  This routine serves as a wrapper for the lower-level routine
397c4762a1bSJed Brown!  "InitialGuessLocal", where the actual computations are
398c4762a1bSJed Brown!  done using the standard Fortran style of treating the local
399c4762a1bSJed Brown!  vector data as a multidimensional array over the local mesh.
400c4762a1bSJed Brown!  This routine merely handles ghost point scatters and accesses
401ce78bad3SBarry Smith!  the local vector data via VecGetArray() and VecRestoreArray().
402c4762a1bSJed Brown!
403c4762a1bSJed Brown      subroutine FormInitialGuess(mysnes,Xnest,ierr)
404c4762a1bSJed Brown#include <petsc/finclude/petscsnes.h>
405c4762a1bSJed Brown      use petscsnes
40677d968b7SBarry Smith      use ex73f90tmodule
40777d968b7SBarry Smith      use ex73f90tmodule_interfaces
408c4762a1bSJed Brown      implicit none
409c4762a1bSJed Brown!  Input/output variables:
410c4762a1bSJed Brown      SNES::     mysnes
411c4762a1bSJed Brown      Vec::      Xnest
412c4762a1bSJed Brown      PetscErrorCode ierr
413c4762a1bSJed Brown
414c4762a1bSJed Brown!  Declarations for use with local arrays:
41577d968b7SBarry Smith      type(ex73f90tmodule_type), pointer:: solver
416c4762a1bSJed Brown      Vec::      Xsub(2)
417c4762a1bSJed Brown      PetscInt::  izero,ione,itwo
418c4762a1bSJed Brown
419c4762a1bSJed Brown      izero = 0
420c4762a1bSJed Brown      ione = 1
421c4762a1bSJed Brown      itwo = 2
422c4762a1bSJed Brown      ierr = 0
423d8606c27SBarry Smith      PetscCall(SNESGetApplicationContext(mysnes,solver,ierr))
4245d83a8b1SBarry Smith      PetscCall(DMCompositeGetAccessArray(solver%da,Xnest,itwo,PETSC_NULL_INTEGER_ARRAY,Xsub,ierr))
425c4762a1bSJed Brown
426d8606c27SBarry Smith      PetscCall(InitialGuessLocal(solver,Xsub(1),ierr))
427d8606c27SBarry Smith      PetscCall(VecAssemblyBegin(Xsub(1),ierr))
428d8606c27SBarry Smith      PetscCall(VecAssemblyEnd(Xsub(1),ierr))
429c4762a1bSJed Brown
430c4762a1bSJed Brown!     zero out lambda
431d8606c27SBarry Smith      PetscCall(VecZeroEntries(Xsub(2),ierr))
4325d83a8b1SBarry Smith      PetscCall(DMCompositeRestoreAccessArray(solver%da,Xnest,itwo,PETSC_NULL_INTEGER_ARRAY,Xsub,ierr))
433c4762a1bSJed Brown
434c4762a1bSJed Brown      end subroutine FormInitialGuess
435c4762a1bSJed Brown
436c4762a1bSJed Brown! ---------------------------------------------------------------------
437c4762a1bSJed Brown!
438c4762a1bSJed Brown!  InitialGuessLocal - Computes initial approximation, called by
439c4762a1bSJed Brown!  the higher level routine FormInitialGuess().
440c4762a1bSJed Brown!
441c4762a1bSJed Brown!  Input Parameter:
442c4762a1bSJed Brown!  X1 - local vector data
443c4762a1bSJed Brown!
444c4762a1bSJed Brown!  Output Parameters:
445c4762a1bSJed Brown!  x - local vector data
446c4762a1bSJed Brown!  ierr - error code
447c4762a1bSJed Brown!
448c4762a1bSJed Brown!  Notes:
449c4762a1bSJed Brown!  This routine uses standard Fortran-style computations over a 2-dim array.
450c4762a1bSJed Brown!
451c4762a1bSJed Brown      subroutine InitialGuessLocal(solver,X1,ierr)
452c4762a1bSJed Brown#include <petsc/finclude/petscsys.h>
453c4762a1bSJed Brown      use petscsys
45477d968b7SBarry Smith      use ex73f90tmodule
455c4762a1bSJed Brown      implicit none
456c4762a1bSJed Brown!  Input/output variables:
45777d968b7SBarry Smith      type (ex73f90tmodule_type)         solver
458c4762a1bSJed Brown      Vec::      X1
459c4762a1bSJed Brown      PetscErrorCode ierr
460c4762a1bSJed Brown
461c4762a1bSJed Brown!  Local variables:
462c4762a1bSJed Brown      PetscInt      row,i,j,ione,low,high
463c4762a1bSJed Brown      PetscReal   temp1,temp,hx,hy,v
464c4762a1bSJed Brown      PetscReal   one
465c4762a1bSJed Brown
466c4762a1bSJed Brown!  Set parameters
467c4762a1bSJed Brown      ione = 1
468c4762a1bSJed Brown      ierr   = 0
469c4762a1bSJed Brown      one    = 1.0
470c4762a1bSJed Brown      hx     = one/(solver%mx-1)
471c4762a1bSJed Brown      hy     = one/(solver%my-1)
472c4762a1bSJed Brown      temp1  = solver%lambda/(solver%lambda + one) + one
473c4762a1bSJed Brown
474d8606c27SBarry Smith      PetscCall(VecGetOwnershipRange(X1,low,high,ierr))
475c4762a1bSJed Brown
476c4762a1bSJed Brown      do 20 row=low,high-1
477c4762a1bSJed Brown         j = row/solver%mx
478c4762a1bSJed Brown         i = mod(row,solver%mx)
479c4762a1bSJed Brown         temp = min(j,solver%my-j+1)*hy
480c4762a1bSJed Brown         if (i .eq. 0 .or. j .eq. 0  .or. i .eq. solver%mx-1 .or. j .eq. solver%my-1) then
481c4762a1bSJed Brown            v = 0.0
482c4762a1bSJed Brown         else
483c4762a1bSJed Brown            v = temp1 * sqrt(min(min(i,solver%mx-i+1)*hx,temp))
484c4762a1bSJed Brown         endif
4855d83a8b1SBarry Smith         PetscCall(VecSetValues(X1,ione,[row],[v],INSERT_VALUES,ierr))
486c4762a1bSJed Brown 20   continue
487c4762a1bSJed Brown
488c4762a1bSJed Brown      end subroutine InitialGuessLocal
489c4762a1bSJed Brown
490c4762a1bSJed Brown! ---------------------------------------------------------------------
491c4762a1bSJed Brown!
492c4762a1bSJed Brown!  FormJacobian - Evaluates Jacobian matrix.
493c4762a1bSJed Brown!
494c4762a1bSJed Brown!  Input Parameters:
495c4762a1bSJed Brown!  dummy     - the SNES context
496c4762a1bSJed Brown!  x         - input vector
497c4762a1bSJed Brown!  solver    - solver data
498c4762a1bSJed Brown!
499c4762a1bSJed Brown!  Output Parameters:
500c4762a1bSJed Brown!  jac      - Jacobian matrix
501*7addb90fSBarry Smith!  jac_prec - optionally different matrix used to construct the preconditioner (not used here)
502c4762a1bSJed Brown!
503c4762a1bSJed Brown      subroutine FormJacobian(dummy,X,jac,jac_prec,solver,ierr)
504c4762a1bSJed Brown#include <petsc/finclude/petscsnes.h>
505c4762a1bSJed Brown      use petscsnes
50677d968b7SBarry Smith      use ex73f90tmodule
507c4762a1bSJed Brown      implicit none
508c4762a1bSJed Brown!  Input/output variables:
509c4762a1bSJed Brown      SNES::     dummy
510c4762a1bSJed Brown      Vec::      X
511c4762a1bSJed Brown     Mat::     jac,jac_prec
51277d968b7SBarry Smith      type(ex73f90tmodule_type)  solver
513c4762a1bSJed Brown      PetscErrorCode ierr
514c4762a1bSJed Brown
515c4762a1bSJed Brown!  Declarations for use with local arrays:
516c4762a1bSJed Brown      Vec::      Xsub(1)
517c4762a1bSJed Brown     Mat::     Amat
5182a27bf02SStefano Zampini      PetscInt       ione
519c4762a1bSJed Brown
520c4762a1bSJed Brown      ione = 1
521c4762a1bSJed Brown
5225d83a8b1SBarry Smith      PetscCall(DMCompositeGetAccessArray(solver%da,X,ione,PETSC_NULL_INTEGER_ARRAY,Xsub,ierr))
523c4762a1bSJed Brown
524c4762a1bSJed Brown!     Compute entries for the locally owned part of the Jacobian preconditioner.
525d8606c27SBarry Smith      PetscCall(MatCreateSubMatrix(jac_prec,solver%isPhi,solver%isPhi,MAT_INITIAL_MATRIX,Amat,ierr))
526c4762a1bSJed Brown
527d8606c27SBarry Smith      PetscCall(FormJacobianLocal(Xsub(1),Amat,solver,.true.,ierr))
528d8606c27SBarry Smith      PetscCall(MatDestroy(Amat,ierr)) ! discard our reference
5295d83a8b1SBarry Smith      PetscCall(DMCompositeRestoreAccessArray(solver%da,X,ione,PETSC_NULL_INTEGER_ARRAY,Xsub,ierr))
530c4762a1bSJed Brown
531c4762a1bSJed Brown      ! the rest of the matrix is not touched
532d8606c27SBarry Smith      PetscCall(MatAssemblyBegin(jac_prec,MAT_FINAL_ASSEMBLY,ierr))
533d8606c27SBarry Smith      PetscCall(MatAssemblyEnd(jac_prec,MAT_FINAL_ASSEMBLY,ierr))
534c4762a1bSJed Brown      if (jac .ne. jac_prec) then
535d8606c27SBarry Smith         PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr))
536d8606c27SBarry Smith         PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr))
537c4762a1bSJed Brown      end if
538c4762a1bSJed Brown
539c4762a1bSJed Brown!     Tell the matrix we will never add a new nonzero location to the
540c4762a1bSJed Brown!     matrix. If we do it will generate an error.
541d8606c27SBarry Smith      PetscCall(MatSetOption(jac_prec,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE,ierr))
542c4762a1bSJed Brown
543c4762a1bSJed Brown      end subroutine FormJacobian
544c4762a1bSJed Brown
545c4762a1bSJed Brown! ---------------------------------------------------------------------
546c4762a1bSJed Brown!
547*7addb90fSBarry Smith!  FormJacobianLocal - Computes Jacobian matrix used to compute the preconditioner,
548c4762a1bSJed Brown!  called by the higher level routine FormJacobian().
549c4762a1bSJed Brown!
550c4762a1bSJed Brown!  Input Parameters:
551c4762a1bSJed Brown!  x        - local vector data
552c4762a1bSJed Brown!
553c4762a1bSJed Brown!  Output Parameters:
554*7addb90fSBarry Smith!  jac - Jacobian matrix used to compute the preconditioner
555c4762a1bSJed Brown!  ierr     - error code
556c4762a1bSJed Brown!
557c4762a1bSJed Brown!  Notes:
558c4762a1bSJed Brown!  This routine uses standard Fortran-style computations over a 2-dim array.
559c4762a1bSJed Brown!
560c4762a1bSJed Brown      subroutine FormJacobianLocal(X1,jac,solver,add_nl_term,ierr)
561c4762a1bSJed Brown#include <petsc/finclude/petscmat.h>
56277d968b7SBarry Smith      use ex73f90tmodule
563c4762a1bSJed Brown      implicit none
564c4762a1bSJed Brown!  Input/output variables:
56577d968b7SBarry Smith      type (ex73f90tmodule_type) solver
566c4762a1bSJed Brown      Vec::      X1
567c4762a1bSJed Brown     Mat::     jac
568c4762a1bSJed Brown      logical        add_nl_term
569c4762a1bSJed Brown      PetscErrorCode ierr
570c4762a1bSJed Brown
571c4762a1bSJed Brown!  Local variables:
572c4762a1bSJed Brown      PetscInt    irow,row(1),col(5),i,j
573c4762a1bSJed Brown      PetscInt    ione,ifive,low,high,ii
574c4762a1bSJed Brown      PetscScalar two,one,hx,hy,hy2inv
575c4762a1bSJed Brown      PetscScalar hx2inv,sc,v(5)
576c4762a1bSJed Brown      PetscScalar,pointer :: lx_v(:)
577c4762a1bSJed Brown
578c4762a1bSJed Brown!  Set parameters
579c4762a1bSJed Brown      ione   = 1
580c4762a1bSJed Brown      ifive  = 5
581c4762a1bSJed Brown      one    = 1.0
582c4762a1bSJed Brown      two    = 2.0
583c4762a1bSJed Brown      hx     = one/(solver%mx-1)
584c4762a1bSJed Brown      hy     = one/(solver%my-1)
585c4762a1bSJed Brown      sc     = solver%lambda
586c4762a1bSJed Brown      hx2inv = one/(hx*hx)
587c4762a1bSJed Brown      hy2inv = one/(hy*hy)
588c4762a1bSJed Brown
589d8606c27SBarry Smith      PetscCall(VecGetOwnershipRange(X1,low,high,ierr))
590ce78bad3SBarry Smith      PetscCall(VecGetArrayRead(X1,lx_v,ierr))
591c4762a1bSJed Brown
592c4762a1bSJed Brown      ii = 0
593c4762a1bSJed Brown      do 20 irow=low,high-1
594c4762a1bSJed Brown         j = irow/solver%mx
595c4762a1bSJed Brown         i = mod(irow,solver%mx)
596c4762a1bSJed Brown         ii = ii + 1            ! one based local index
597c4762a1bSJed Brown!     boundary points
598c4762a1bSJed Brown         if (i .eq. 0 .or. j .eq. 0 .or. i .eq. solver%mx-1 .or. j .eq. solver%my-1) then
599c4762a1bSJed Brown            col(1) = irow
600c4762a1bSJed Brown            row(1) = irow
601c4762a1bSJed Brown            v(1)   = one
602d8606c27SBarry Smith            PetscCall(MatSetValues(jac,ione,row,ione,col,v,INSERT_VALUES,ierr))
603c4762a1bSJed Brown!     interior grid points
604c4762a1bSJed Brown         else
605c4762a1bSJed Brown            v(1) = -hy2inv
606c4762a1bSJed Brown            if (j-1==0) v(1) = 0.0
607c4762a1bSJed Brown            v(2) = -hx2inv
608c4762a1bSJed Brown            if (i-1==0) v(2) = 0.0
609c4762a1bSJed Brown            v(3) = two*(hx2inv + hy2inv)
610c4762a1bSJed Brown            if (add_nl_term) v(3) = v(3) - sc*exp(lx_v(ii))
611c4762a1bSJed Brown            v(4) = -hx2inv
612c4762a1bSJed Brown            if (i+1==solver%mx-1) v(4) = 0.0
613c4762a1bSJed Brown            v(5) = -hy2inv
614c4762a1bSJed Brown            if (j+1==solver%my-1) v(5) = 0.0
615c4762a1bSJed Brown            col(1) = irow - solver%mx
616c4762a1bSJed Brown            col(2) = irow - 1
617c4762a1bSJed Brown            col(3) = irow
618c4762a1bSJed Brown            col(4) = irow + 1
619c4762a1bSJed Brown            col(5) = irow + solver%mx
620c4762a1bSJed Brown            row(1) = irow
621d8606c27SBarry Smith            PetscCall(MatSetValues(jac,ione,row,ifive,col,v,INSERT_VALUES,ierr))
622c4762a1bSJed Brown         endif
623c4762a1bSJed Brown 20   continue
624c4762a1bSJed Brown
625ce78bad3SBarry Smith      PetscCall(VecRestoreArrayRead(X1,lx_v,ierr))
626c4762a1bSJed Brown
627c4762a1bSJed Brown      end subroutine FormJacobianLocal
628c4762a1bSJed Brown
629c4762a1bSJed Brown! ---------------------------------------------------------------------
630c4762a1bSJed Brown!
631c4762a1bSJed Brown!  FormFunction - Evaluates nonlinear function, F(x).
632c4762a1bSJed Brown!
633c4762a1bSJed Brown!  Input Parameters:
634c4762a1bSJed Brown!  snes - the SNES context
635c4762a1bSJed Brown!  X - input vector
636c4762a1bSJed Brown!  dummy - optional user-defined context, as set by SNESSetFunction()
637c4762a1bSJed Brown!          (not used here)
638c4762a1bSJed Brown!
639c4762a1bSJed Brown!  Output Parameter:
640c4762a1bSJed Brown!  F - function vector
641c4762a1bSJed Brown!
642c4762a1bSJed Brown      subroutine FormFunction(snesIn,X,F,solver,ierr)
643c4762a1bSJed Brown#include <petsc/finclude/petscsnes.h>
644c4762a1bSJed Brown      use petscsnes
64577d968b7SBarry Smith      use ex73f90tmodule
646c4762a1bSJed Brown      implicit none
647c4762a1bSJed Brown!  Input/output variables:
648c4762a1bSJed Brown      SNES::     snesIn
649c4762a1bSJed Brown     Vec::      X,F
650c4762a1bSJed Brown      PetscErrorCode ierr
65177d968b7SBarry Smith      type (ex73f90tmodule_type) solver
652c4762a1bSJed Brown
653c4762a1bSJed Brown!  Declarations for use with local arrays:
654c4762a1bSJed Brown     Vec::              Xsub(2),Fsub(2)
6552a27bf02SStefano Zampini      PetscInt               itwo
656c4762a1bSJed Brown
657c4762a1bSJed Brown!  Scatter ghost points to local vector, using the 2-step process
658c4762a1bSJed Brown!     DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
659c4762a1bSJed Brown!  By placing code between these two statements, computations can
660c4762a1bSJed Brown!  be done while messages are in transition.
661c4762a1bSJed Brown
662c4762a1bSJed Brown      itwo = 2
6635d83a8b1SBarry Smith      PetscCall(DMCompositeGetAccessArray(solver%da,X,itwo,PETSC_NULL_INTEGER_ARRAY,Xsub,ierr))
6645d83a8b1SBarry Smith      PetscCall(DMCompositeGetAccessArray(solver%da,F,itwo,PETSC_NULL_INTEGER_ARRAY,Fsub,ierr))
665c4762a1bSJed Brown
666d8606c27SBarry Smith      PetscCall(FormFunctionNLTerm( Xsub(1), Fsub(1), solver, ierr))
667d8606c27SBarry Smith      PetscCall(MatMultAdd( solver%AmatLin, Xsub(1), Fsub(1), Fsub(1), ierr))
668c4762a1bSJed Brown
669c4762a1bSJed Brown!     do rest of operator (linear)
670d8606c27SBarry Smith      PetscCall(MatMult(    solver%Cmat, Xsub(1),      Fsub(2), ierr))
671d8606c27SBarry Smith      PetscCall(MatMultAdd( solver%Bmat, Xsub(2), Fsub(1), Fsub(1), ierr))
672d8606c27SBarry Smith      PetscCall(MatMultAdd( solver%Dmat, Xsub(2), Fsub(2), Fsub(2), ierr))
673c4762a1bSJed Brown
6745d83a8b1SBarry Smith      PetscCall(DMCompositeRestoreAccessArray(solver%da,X,itwo,PETSC_NULL_INTEGER_ARRAY,Xsub,ierr))
6755d83a8b1SBarry Smith      PetscCall(DMCompositeRestoreAccessArray(solver%da,F,itwo,PETSC_NULL_INTEGER_ARRAY,Fsub,ierr))
676c4762a1bSJed Brown      end subroutine formfunction
677c4762a1bSJed Brown
678c4762a1bSJed Brown! ---------------------------------------------------------------------
679c4762a1bSJed Brown!
680c4762a1bSJed Brown!  FormFunctionNLTerm - Computes nonlinear function, called by
681c4762a1bSJed Brown!  the higher level routine FormFunction().
682c4762a1bSJed Brown!
683c4762a1bSJed Brown!  Input Parameter:
684c4762a1bSJed Brown!  x - local vector data
685c4762a1bSJed Brown!
686c4762a1bSJed Brown!  Output Parameters:
687c4762a1bSJed Brown!  f - local vector data, f(x)
688c4762a1bSJed Brown!  ierr - error code
689c4762a1bSJed Brown!
690c4762a1bSJed Brown!  Notes:
691c4762a1bSJed Brown!  This routine uses standard Fortran-style computations over a 2-dim array.
692c4762a1bSJed Brown!
693c4762a1bSJed Brown      subroutine FormFunctionNLTerm(X1,F1,solver,ierr)
694c4762a1bSJed Brown#include <petsc/finclude/petscvec.h>
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
710ce78bad3SBarry Smith      PetscCall(VecGetArrayRead(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
729ce78bad3SBarry Smith      PetscCall(VecRestoreArrayRead(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!
746c00ad2bcSBarry Smith!   test:
747c00ad2bcSBarry Smith!      args: -snes_linesearch_type {{l2 cp}separate output} -objective {{false true}shared output}
748c00ad2bcSBarry Smith!
749c00ad2bcSBarry Smith!   test:
750c00ad2bcSBarry Smith!      args: -snes_linesearch_type bt -objective {{false true}separate output}
751c00ad2bcSBarry Smith!
752c4762a1bSJed Brown!TEST*/
753