xref: /petsc/src/snes/tutorials/ex5f90t.F90 (revision dcb3e68992f1c4897946af7e8406e2b4165e50f2)
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!
10c4762a1bSJed Brown!  --------------------------------------------------------------------------
11c4762a1bSJed Brown!
12c4762a1bSJed Brown!  Solid Fuel Ignition (SFI) problem.  This problem is modeled by
13c4762a1bSJed Brown!  the partial differential equation
14c4762a1bSJed Brown!
15c4762a1bSJed Brown!          -Laplacian u - lambda*exp(u) = 0,  0 < x,y < 1,
16c4762a1bSJed Brown!
17c4762a1bSJed Brown!  with boundary conditions
18c4762a1bSJed Brown!
19c4762a1bSJed Brown!           u = 0  for  x = 0, x = 1, y = 0, y = 1.
20c4762a1bSJed Brown!
21c4762a1bSJed Brown!  A finite difference approximation with the usual 5-point stencil
22c4762a1bSJed Brown!  is used to discretize the boundary value problem to obtain a nonlinear
23c4762a1bSJed Brown!  system of equations.
24c4762a1bSJed Brown!
25c4762a1bSJed Brown!  The uniprocessor version of this code is snes/tutorials/ex4f.F
26c4762a1bSJed Brown!
27c4762a1bSJed Brown!  --------------------------------------------------------------------------
28c4762a1bSJed Brown!  The following define must be used before including any PETSc include files
29c4762a1bSJed Brown!  into a module or interface. This is because they can't handle declarations
30c4762a1bSJed Brown!  in them
31c4762a1bSJed Brown!
32c4762a1bSJed Brown
3377d968b7SBarry Smith      module ex5f90tmodule
34c4762a1bSJed Brown#include <petsc/finclude/petscdm.h>
35c4762a1bSJed Brown      use petscdmdef
36c4762a1bSJed Brown      type userctx
37c4762a1bSJed Brown        type(tDM) da
38c4762a1bSJed Brown        PetscInt xs,xe,xm,gxs,gxe,gxm
39c4762a1bSJed Brown        PetscInt ys,ye,ym,gys,gye,gym
40c4762a1bSJed Brown        PetscInt mx,my
41c4762a1bSJed Brown        PetscMPIInt rank
42c4762a1bSJed Brown        PetscReal lambda
43c4762a1bSJed Brown      end type userctx
44c4762a1bSJed Brown
45c4762a1bSJed Brown      contains
46c4762a1bSJed Brown! ---------------------------------------------------------------------
47c4762a1bSJed Brown!
48c4762a1bSJed Brown!  FormFunction - Evaluates nonlinear function, F(x).
49c4762a1bSJed Brown!
50c4762a1bSJed Brown!  Input Parameters:
51c4762a1bSJed Brown!  snes - the SNES context
52c4762a1bSJed Brown!  X - input vector
53c4762a1bSJed Brown!  dummy - optional user-defined context, as set by SNESSetFunction()
54c4762a1bSJed Brown!          (not used here)
55c4762a1bSJed Brown!
56c4762a1bSJed Brown!  Output Parameter:
57c4762a1bSJed Brown!  F - function vector
58c4762a1bSJed Brown!
59c4762a1bSJed Brown!  Notes:
60c4762a1bSJed Brown!  This routine serves as a wrapper for the lower-level routine
61c4762a1bSJed Brown!  "FormFunctionLocal", where the actual computations are
62c4762a1bSJed Brown!  done using the standard Fortran style of treating the local
63c4762a1bSJed Brown!  vector data as a multidimensional array over the local mesh.
64c4762a1bSJed Brown!  This routine merely handles ghost point scatters and accesses
65c4762a1bSJed Brown!  the local vector data via VecGetArrayF90() and VecRestoreArrayF90().
66c4762a1bSJed Brown!
67c4762a1bSJed Brown      subroutine FormFunction(snesIn,X,F,user,ierr)
68c4762a1bSJed Brown#include <petsc/finclude/petscsnes.h>
69c4762a1bSJed Brown      use petscsnes
70c4762a1bSJed Brown
71c4762a1bSJed Brown!  Input/output variables:
72c4762a1bSJed Brown      type(tSNES)     snesIn
73c4762a1bSJed Brown      type(tVec)      X,F
74c4762a1bSJed Brown      PetscErrorCode ierr
75c4762a1bSJed Brown      type (userctx) user
76c4762a1bSJed Brown
77c4762a1bSJed Brown!  Declarations for use with local arrays:
78c4762a1bSJed Brown      PetscScalar,pointer :: lx_v(:),lf_v(:)
79c4762a1bSJed Brown      type(tVec)              localX
80c4762a1bSJed Brown
81c4762a1bSJed Brown!  Scatter ghost points to local vector, using the 2-step process
82c4762a1bSJed Brown!     DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
83c4762a1bSJed Brown!  By placing code between these two statements, computations can
84c4762a1bSJed Brown!  be done while messages are in transition.
85d8606c27SBarry Smith      PetscCall(DMGetLocalVector(user%da,localX,ierr))
86d8606c27SBarry Smith      PetscCall(DMGlobalToLocalBegin(user%da,X,INSERT_VALUES,localX,ierr))
87d8606c27SBarry Smith      PetscCall(DMGlobalToLocalEnd(user%da,X,INSERT_VALUES,localX,ierr))
88c4762a1bSJed Brown
89c4762a1bSJed Brown!  Get a pointer to vector data.
9042ce371bSBarry Smith!    - VecGetArray90() returns a pointer to the data array.
91c4762a1bSJed Brown!    - You MUST call VecRestoreArrayF90() when you no longer need access to
92c4762a1bSJed Brown!      the array.
93c4762a1bSJed Brown
94d8606c27SBarry Smith      PetscCall(VecGetArrayF90(localX,lx_v,ierr))
95d8606c27SBarry Smith      PetscCall(VecGetArrayF90(F,lf_v,ierr))
96c4762a1bSJed Brown
97c4762a1bSJed Brown!  Compute function over the locally owned part of the grid
98d8606c27SBarry Smith      PetscCall(FormFunctionLocal(lx_v,lf_v,user,ierr))
99c4762a1bSJed Brown
100c4762a1bSJed Brown!  Restore vectors
101d8606c27SBarry Smith      PetscCall(VecRestoreArrayF90(localX,lx_v,ierr))
102d8606c27SBarry Smith      PetscCall(VecRestoreArrayF90(F,lf_v,ierr))
103c4762a1bSJed Brown
104c4762a1bSJed Brown!  Insert values into global vector
105c4762a1bSJed Brown
106d8606c27SBarry Smith      PetscCall(DMRestoreLocalVector(user%da,localX,ierr))
107d8606c27SBarry Smith      PetscCall(PetscLogFlops(11.0d0*user%ym*user%xm,ierr))
108c4762a1bSJed Brown
109d8606c27SBarry Smith!      PetscCall(VecView(X,PETSC_VIEWER_STDOUT_WORLD,ierr))
110d8606c27SBarry Smith!      PetscCall(VecView(F,PETSC_VIEWER_STDOUT_WORLD,ierr))
111c4762a1bSJed Brown      return
112c4762a1bSJed Brown      end subroutine formfunction
11377d968b7SBarry Smith      end module ex5f90tmodule
114c4762a1bSJed Brown
115c4762a1bSJed Brown      module f90moduleinterfacest
11677d968b7SBarry Smith        use ex5f90tmodule
117c4762a1bSJed Brown
118c4762a1bSJed Brown      Interface SNESSetApplicationContext
119c4762a1bSJed Brown        Subroutine SNESSetApplicationContext(snesIn,ctx,ierr)
120c4762a1bSJed Brown#include <petsc/finclude/petscsnes.h>
121c4762a1bSJed Brown        use petscsnes
12277d968b7SBarry Smith        use ex5f90tmodule
123c4762a1bSJed Brown          type(tSNES)    snesIn
124c4762a1bSJed Brown          type(userctx) ctx
125c4762a1bSJed Brown          PetscErrorCode ierr
126c4762a1bSJed Brown        End Subroutine
127c4762a1bSJed Brown      End Interface SNESSetApplicationContext
128c4762a1bSJed Brown
129c4762a1bSJed Brown      Interface SNESGetApplicationContext
130c4762a1bSJed Brown        Subroutine SNESGetApplicationContext(snesIn,ctx,ierr)
131c4762a1bSJed Brown#include <petsc/finclude/petscsnes.h>
132c4762a1bSJed Brown        use petscsnes
13377d968b7SBarry Smith        use ex5f90tmodule
134c4762a1bSJed Brown          type(tSNES)     snesIn
135c4762a1bSJed Brown          type(userctx), pointer :: ctx
136c4762a1bSJed Brown          PetscErrorCode ierr
137c4762a1bSJed Brown        End Subroutine
138c4762a1bSJed Brown      End Interface SNESGetApplicationContext
139c4762a1bSJed Brown      end module f90moduleinterfacest
140c4762a1bSJed Brown
141c4762a1bSJed Brown      program main
142c4762a1bSJed Brown#include <petsc/finclude/petscdm.h>
143c4762a1bSJed Brown#include <petsc/finclude/petscsnes.h>
144c4762a1bSJed Brown      use petscdmda
145c4762a1bSJed Brown      use petscdm
146c4762a1bSJed Brown      use petscsnes
14777d968b7SBarry Smith      use ex5f90tmodule
148c4762a1bSJed Brown      use f90moduleinterfacest
149c4762a1bSJed Brown      implicit none
150c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151c4762a1bSJed Brown!                   Variable declarations
152c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153c4762a1bSJed Brown!
154c4762a1bSJed Brown!  Variables:
155c4762a1bSJed Brown!     mysnes      - nonlinear solver
156c4762a1bSJed Brown!     x, r        - solution, residual vectors
157c4762a1bSJed Brown!     J           - Jacobian matrix
158c4762a1bSJed Brown!     its         - iterations for convergence
159c4762a1bSJed Brown!     Nx, Ny      - number of preocessors in x- and y- directions
160c4762a1bSJed Brown!     matrix_free - flag - 1 indicates matrix-free version
161c4762a1bSJed Brown!
162c4762a1bSJed Brown      type(tSNES)       mysnes
163c4762a1bSJed Brown      type(tVec)        x,r
164c4762a1bSJed Brown      type(tMat)        J
165c4762a1bSJed Brown      PetscErrorCode   ierr
166c4762a1bSJed Brown      PetscInt         its
167c4762a1bSJed Brown      PetscBool        flg,matrix_free,set
168c4762a1bSJed Brown      PetscInt         ione,nfour
169c4762a1bSJed Brown      PetscReal lambda_max,lambda_min
170c4762a1bSJed Brown      type(userctx)    user
171c4762a1bSJed Brown      type(userctx), pointer:: puser
172c4762a1bSJed Brown      type(tPetscOptions) :: options
173c4762a1bSJed Brown
174c4762a1bSJed Brown!  Note: Any user-defined Fortran routines (such as FormJacobian)
175c4762a1bSJed Brown!  MUST be declared as external.
176c4762a1bSJed Brown      external FormInitialGuess,FormJacobian
177c4762a1bSJed Brown
178c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
179c4762a1bSJed Brown!  Initialize program
180c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
181d8606c27SBarry Smith      PetscCallA(PetscInitialize(ierr))
182d8606c27SBarry Smith      PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,user%rank,ierr))
183c4762a1bSJed Brown
184c4762a1bSJed Brown!  Initialize problem parameters
185c4762a1bSJed Brown      options%v = 0
186c4762a1bSJed Brown      lambda_max  = 6.81
187c4762a1bSJed Brown      lambda_min  = 0.0
188c4762a1bSJed Brown      user%lambda = 6.0
189c4762a1bSJed Brown      ione = 1
190c4762a1bSJed Brown      nfour = 4
191d8606c27SBarry Smith      PetscCallA(PetscOptionsGetReal(options,PETSC_NULL_CHARACTER,'-par',user%lambda,flg,ierr))
192*dcb3e689SBarry Smith      PetscCheckA(user%lambda .lt. lambda_max .and. user%lambda .gt. lambda_min,PETSC_COMM_SELF,PETSC_ERR_USER,'Lambda provided with -par is out of range')
193c4762a1bSJed Brown
194c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
195c4762a1bSJed Brown!  Create nonlinear solver context
196c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
197d8606c27SBarry Smith      PetscCallA(SNESCreate(PETSC_COMM_WORLD,mysnes,ierr))
198c4762a1bSJed Brown
199c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
200c4762a1bSJed Brown!  Create vector data structures; set function evaluation routine
201c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
202c4762a1bSJed Brown
203c4762a1bSJed Brown!  Create distributed array (DMDA) to manage parallel grid and vectors
204c4762a1bSJed Brown
205c4762a1bSJed Brown! This really needs only the star-type stencil, but we use the box
206c4762a1bSJed Brown! stencil temporarily.
207d8606c27SBarry 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,user%da,ierr))
208d8606c27SBarry Smith      PetscCallA(DMSetFromOptions(user%da,ierr))
209d8606c27SBarry Smith      PetscCallA(DMSetUp(user%da,ierr))
210d8606c27SBarry Smith      PetscCallA(DMDAGetInfo(user%da,PETSC_NULL_INTEGER,user%mx,user%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))
211c4762a1bSJed Brown
212c4762a1bSJed Brown!
213c4762a1bSJed Brown!   Visualize the distribution of the array across the processors
214c4762a1bSJed Brown!
215d8606c27SBarry Smith!     PetscCallA(DMView(user%da,PETSC_VIEWER_DRAW_WORLD,ierr))
216c4762a1bSJed Brown
217c4762a1bSJed Brown!  Extract global and local vectors from DMDA; then duplicate for remaining
218c4762a1bSJed Brown!  vectors that are the same types
219d8606c27SBarry Smith      PetscCallA(DMCreateGlobalVector(user%da,x,ierr))
220d8606c27SBarry Smith      PetscCallA(VecDuplicate(x,r,ierr))
221c4762a1bSJed Brown
222c4762a1bSJed Brown!  Get local grid boundaries (for 2-dimensional DMDA)
223d8606c27SBarry Smith      PetscCallA(DMDAGetCorners(user%da,user%xs,user%ys,PETSC_NULL_INTEGER,user%xm,user%ym,PETSC_NULL_INTEGER,ierr))
224d8606c27SBarry Smith      PetscCallA(DMDAGetGhostCorners(user%da,user%gxs,user%gys,PETSC_NULL_INTEGER,user%gxm,user%gym,PETSC_NULL_INTEGER,ierr))
225c4762a1bSJed Brown
226c4762a1bSJed Brown!  Here we shift the starting indices up by one so that we can easily
227c4762a1bSJed Brown!  use the Fortran convention of 1-based indices (rather 0-based indices).
228c4762a1bSJed Brown      user%xs  = user%xs+1
229c4762a1bSJed Brown      user%ys  = user%ys+1
230c4762a1bSJed Brown      user%gxs = user%gxs+1
231c4762a1bSJed Brown      user%gys = user%gys+1
232c4762a1bSJed Brown
233c4762a1bSJed Brown      user%ye  = user%ys+user%ym-1
234c4762a1bSJed Brown      user%xe  = user%xs+user%xm-1
235c4762a1bSJed Brown      user%gye = user%gys+user%gym-1
236c4762a1bSJed Brown      user%gxe = user%gxs+user%gxm-1
237c4762a1bSJed Brown
238d8606c27SBarry Smith      PetscCallA(SNESSetApplicationContext(mysnes,user,ierr))
239c4762a1bSJed Brown
240c4762a1bSJed Brown!  Set function evaluation routine and vector
241d8606c27SBarry Smith      PetscCallA(SNESSetFunction(mysnes,r,FormFunction,user,ierr))
242c4762a1bSJed Brown
243c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
244c4762a1bSJed Brown!  Create matrix data structure; set Jacobian evaluation routine
245c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
246c4762a1bSJed Brown
247c4762a1bSJed Brown!  Set Jacobian matrix data structure and default Jacobian evaluation
248c4762a1bSJed Brown!  routine. User can override with:
249c4762a1bSJed Brown!     -snes_fd : default finite differencing approximation of Jacobian
250c4762a1bSJed Brown!     -snes_mf : matrix-free Newton-Krylov method with no preconditioning
251c4762a1bSJed Brown!                (unless user explicitly sets preconditioner)
252c4762a1bSJed Brown!     -snes_mf_operator : form preconditioning matrix as set by the user,
253c4762a1bSJed Brown!                         but use matrix-free approx for Jacobian-vector
254c4762a1bSJed Brown!                         products within Newton-Krylov method
255c4762a1bSJed Brown!
256c4762a1bSJed Brown!  Note:  For the parallel case, vectors and matrices MUST be partitioned
257c4762a1bSJed Brown!     accordingly.  When using distributed arrays (DMDAs) to create vectors,
258c4762a1bSJed Brown!     the DMDAs determine the problem partitioning.  We must explicitly
259c4762a1bSJed Brown!     specify the local matrix dimensions upon its creation for compatibility
260c4762a1bSJed Brown!     with the vector distribution.  Thus, the generic MatCreate() routine
261c4762a1bSJed Brown!     is NOT sufficient when working with distributed arrays.
262c4762a1bSJed Brown!
263c4762a1bSJed Brown!     Note: Here we only approximately preallocate storage space for the
264c4762a1bSJed Brown!     Jacobian.  See the users manual for a discussion of better techniques
265c4762a1bSJed Brown!     for preallocating matrix memory.
266c4762a1bSJed Brown
267d8606c27SBarry Smith      PetscCallA(PetscOptionsHasName(options,PETSC_NULL_CHARACTER,'-snes_mf',matrix_free,ierr))
268c4762a1bSJed Brown      if (.not. matrix_free) then
269d8606c27SBarry Smith        PetscCallA(DMSetMatType(user%da,MATAIJ,ierr))
270d8606c27SBarry Smith        PetscCallA(DMCreateMatrix(user%da,J,ierr))
271d8606c27SBarry Smith        PetscCallA(SNESSetJacobian(mysnes,J,J,FormJacobian,user,ierr))
272c4762a1bSJed Brown      endif
273c4762a1bSJed Brown
274c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
275c4762a1bSJed Brown!  Customize nonlinear solver; set runtime options
276c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
277c4762a1bSJed Brown!  Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
278d8606c27SBarry Smith      PetscCallA(SNESSetFromOptions(mysnes,ierr))
279c4762a1bSJed Brown
280c4762a1bSJed Brown!     Test Fortran90 wrapper for SNESSet/Get ApplicationContext()
281d8606c27SBarry Smith      PetscCallA(PetscOptionsGetBool(options,PETSC_NULL_CHARACTER,'-test_appctx',flg,set,ierr))
282c4762a1bSJed Brown      if (flg) then
283d8606c27SBarry Smith        PetscCallA(SNESGetApplicationContext(mysnes,puser,ierr))
284c4762a1bSJed Brown      endif
285c4762a1bSJed Brown
286c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
287c4762a1bSJed Brown!  Evaluate initial guess; then solve nonlinear system.
288c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
289c4762a1bSJed Brown!  Note: The user should initialize the vector, x, with the initial guess
290c4762a1bSJed Brown!  for the nonlinear solver prior to calling SNESSolve().  In particular,
291c4762a1bSJed Brown!  to employ an initial guess of zero, the user should explicitly set
292c4762a1bSJed Brown!  this vector to zero by calling VecSet().
293c4762a1bSJed Brown
294d8606c27SBarry Smith      PetscCallA(FormInitialGuess(mysnes,x,ierr))
295d8606c27SBarry Smith      PetscCallA(SNESSolve(mysnes,PETSC_NULL_VEC,x,ierr))
296d8606c27SBarry Smith      PetscCallA(SNESGetIterationNumber(mysnes,its,ierr))
297c4762a1bSJed Brown      if (user%rank .eq. 0) then
298c4762a1bSJed Brown         write(6,100) its
299c4762a1bSJed Brown      endif
300c4762a1bSJed Brown  100 format('Number of SNES iterations = ',i5)
301c4762a1bSJed Brown
302c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
303c4762a1bSJed Brown!  Free work space.  All PETSc objects should be destroyed when they
304c4762a1bSJed Brown!  are no longer needed.
305c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
306d8606c27SBarry Smith      if (.not. matrix_free) PetscCallA(MatDestroy(J,ierr))
307d8606c27SBarry Smith      PetscCallA(VecDestroy(x,ierr))
308d8606c27SBarry Smith      PetscCallA(VecDestroy(r,ierr))
309d8606c27SBarry Smith      PetscCallA(SNESDestroy(mysnes,ierr))
310d8606c27SBarry Smith      PetscCallA(DMDestroy(user%da,ierr))
311c4762a1bSJed Brown
312d8606c27SBarry Smith      PetscCallA(PetscFinalize(ierr))
313c4762a1bSJed Brown      end
314c4762a1bSJed Brown
315c4762a1bSJed Brown! ---------------------------------------------------------------------
316c4762a1bSJed Brown!
317c4762a1bSJed Brown!  FormInitialGuess - Forms initial approximation.
318c4762a1bSJed Brown!
319c4762a1bSJed Brown!  Input Parameters:
320c4762a1bSJed Brown!  X - vector
321c4762a1bSJed Brown!
322c4762a1bSJed Brown!  Output Parameter:
323c4762a1bSJed Brown!  X - vector
324c4762a1bSJed Brown!
325c4762a1bSJed Brown!  Notes:
326c4762a1bSJed Brown!  This routine serves as a wrapper for the lower-level routine
327c4762a1bSJed Brown!  "InitialGuessLocal", where the actual computations are
328c4762a1bSJed Brown!  done using the standard Fortran style of treating the local
329c4762a1bSJed Brown!  vector data as a multidimensional array over the local mesh.
330c4762a1bSJed Brown!  This routine merely handles ghost point scatters and accesses
331c4762a1bSJed Brown!  the local vector data via VecGetArrayF90() and VecRestoreArrayF90().
332c4762a1bSJed Brown!
333c4762a1bSJed Brown      subroutine FormInitialGuess(mysnes,X,ierr)
334c4762a1bSJed Brown#include <petsc/finclude/petscsnes.h>
335c4762a1bSJed Brown      use petscsnes
33677d968b7SBarry Smith      use ex5f90tmodule
337c4762a1bSJed Brown      use f90moduleinterfacest
338c4762a1bSJed Brown!  Input/output variables:
339c4762a1bSJed Brown      type(tSNES)     mysnes
340c4762a1bSJed Brown      type(userctx), pointer:: puser
341c4762a1bSJed Brown      type(tVec)      X
342c4762a1bSJed Brown      PetscErrorCode ierr
343c4762a1bSJed Brown
344c4762a1bSJed Brown!  Declarations for use with local arrays:
345c4762a1bSJed Brown      PetscScalar,pointer :: lx_v(:)
346c4762a1bSJed Brown
347c4762a1bSJed Brown      ierr = 0
348d8606c27SBarry Smith      PetscCallA(SNESGetApplicationContext(mysnes,puser,ierr))
349c4762a1bSJed Brown!  Get a pointer to vector data.
35042ce371bSBarry Smith!    - VecGetArray90() returns a pointer to the data array.
351c4762a1bSJed Brown!    - You MUST call VecRestoreArrayF90() when you no longer need access to
352c4762a1bSJed Brown!      the array.
353c4762a1bSJed Brown
354d8606c27SBarry Smith      PetscCallA(VecGetArrayF90(X,lx_v,ierr))
355c4762a1bSJed Brown
356c4762a1bSJed Brown!  Compute initial guess over the locally owned part of the grid
357d8606c27SBarry Smith      PetscCallA(InitialGuessLocal(puser,lx_v,ierr))
358c4762a1bSJed Brown
359c4762a1bSJed Brown!  Restore vector
360d8606c27SBarry Smith      PetscCallA(VecRestoreArrayF90(X,lx_v,ierr))
361c4762a1bSJed Brown
362c4762a1bSJed Brown!  Insert values into global vector
363c4762a1bSJed Brown
364c4762a1bSJed Brown      return
365c4762a1bSJed Brown      end
366c4762a1bSJed Brown
367c4762a1bSJed Brown! ---------------------------------------------------------------------
368c4762a1bSJed Brown!
369c4762a1bSJed Brown!  InitialGuessLocal - Computes initial approximation, called by
370c4762a1bSJed Brown!  the higher level routine FormInitialGuess().
371c4762a1bSJed Brown!
372c4762a1bSJed Brown!  Input Parameter:
373c4762a1bSJed Brown!  x - local vector data
374c4762a1bSJed Brown!
375c4762a1bSJed Brown!  Output Parameters:
376c4762a1bSJed Brown!  x - local vector data
377c4762a1bSJed Brown!  ierr - error code
378c4762a1bSJed Brown!
379c4762a1bSJed Brown!  Notes:
380c4762a1bSJed Brown!  This routine uses standard Fortran-style computations over a 2-dim array.
381c4762a1bSJed Brown!
382c4762a1bSJed Brown      subroutine InitialGuessLocal(user,x,ierr)
383c4762a1bSJed Brown#include <petsc/finclude/petscsys.h>
384c4762a1bSJed Brown      use petscsys
38577d968b7SBarry Smith      use ex5f90tmodule
386c4762a1bSJed Brown!  Input/output variables:
387c4762a1bSJed Brown      type (userctx)         user
388c4762a1bSJed Brown      PetscScalar  x(user%xs:user%xe,user%ys:user%ye)
389c4762a1bSJed Brown      PetscErrorCode ierr
390c4762a1bSJed Brown
391c4762a1bSJed Brown!  Local variables:
392c4762a1bSJed Brown      PetscInt  i,j
393c4762a1bSJed Brown      PetscScalar   temp1,temp,hx,hy
394c4762a1bSJed Brown      PetscScalar   one
395c4762a1bSJed Brown
396c4762a1bSJed Brown!  Set parameters
397c4762a1bSJed Brown
398c4762a1bSJed Brown      ierr   = 0
399c4762a1bSJed Brown      one    = 1.0
400c4762a1bSJed Brown      hx     = one/(PetscIntToReal(user%mx-1))
401c4762a1bSJed Brown      hy     = one/(PetscIntToReal(user%my-1))
402c4762a1bSJed Brown      temp1  = user%lambda/(user%lambda + one)
403c4762a1bSJed Brown
404c4762a1bSJed Brown      do 20 j=user%ys,user%ye
405c4762a1bSJed Brown         temp = PetscIntToReal(min(j-1,user%my-j))*hy
406c4762a1bSJed Brown         do 10 i=user%xs,user%xe
407c4762a1bSJed Brown            if (i .eq. 1 .or. j .eq. 1 .or. i .eq. user%mx .or. j .eq. user%my) then
408c4762a1bSJed Brown              x(i,j) = 0.0
409c4762a1bSJed Brown            else
410c4762a1bSJed Brown              x(i,j) = temp1 * sqrt(min(PetscIntToReal(min(i-1,user%mx-i)*hx),PetscIntToReal(temp)))
411c4762a1bSJed Brown            endif
412c4762a1bSJed Brown 10      continue
413c4762a1bSJed Brown 20   continue
414c4762a1bSJed Brown
415c4762a1bSJed Brown      return
416c4762a1bSJed Brown      end
417c4762a1bSJed Brown
418c4762a1bSJed Brown! ---------------------------------------------------------------------
419c4762a1bSJed Brown!
420c4762a1bSJed Brown!  FormFunctionLocal - Computes nonlinear function, called by
421c4762a1bSJed Brown!  the higher level routine FormFunction().
422c4762a1bSJed Brown!
423c4762a1bSJed Brown!  Input Parameter:
424c4762a1bSJed Brown!  x - local vector data
425c4762a1bSJed Brown!
426c4762a1bSJed Brown!  Output Parameters:
427c4762a1bSJed Brown!  f - local vector data, f(x)
428c4762a1bSJed Brown!  ierr - error code
429c4762a1bSJed Brown!
430c4762a1bSJed Brown!  Notes:
431c4762a1bSJed Brown!  This routine uses standard Fortran-style computations over a 2-dim array.
432c4762a1bSJed Brown!
433c4762a1bSJed Brown      subroutine FormFunctionLocal(x,f,user,ierr)
434c4762a1bSJed Brown#include <petsc/finclude/petscsys.h>
435c4762a1bSJed Brown      use petscsys
43677d968b7SBarry Smith      use ex5f90tmodule
437c4762a1bSJed Brown!  Input/output variables:
438c4762a1bSJed Brown      type (userctx) user
439c4762a1bSJed Brown      PetscScalar  x(user%gxs:user%gxe,user%gys:user%gye)
440c4762a1bSJed Brown      PetscScalar  f(user%xs:user%xe,user%ys:user%ye)
441c4762a1bSJed Brown      PetscErrorCode ierr
442c4762a1bSJed Brown
443c4762a1bSJed Brown!  Local variables:
444c4762a1bSJed Brown      PetscScalar two,one,hx,hy,hxdhy,hydhx,sc
445c4762a1bSJed Brown      PetscScalar u,uxx,uyy
446c4762a1bSJed Brown      PetscInt  i,j
447c4762a1bSJed Brown
448c4762a1bSJed Brown      one    = 1.0
449c4762a1bSJed Brown      two    = 2.0
450c4762a1bSJed Brown      hx     = one/PetscIntToReal(user%mx-1)
451c4762a1bSJed Brown      hy     = one/PetscIntToReal(user%my-1)
452c4762a1bSJed Brown      sc     = hx*hy*user%lambda
453c4762a1bSJed Brown      hxdhy  = hx/hy
454c4762a1bSJed Brown      hydhx  = hy/hx
455c4762a1bSJed Brown
456c4762a1bSJed Brown!  Compute function over the locally owned part of the grid
457c4762a1bSJed Brown
458c4762a1bSJed Brown      do 20 j=user%ys,user%ye
459c4762a1bSJed Brown         do 10 i=user%xs,user%xe
460c4762a1bSJed Brown            if (i .eq. 1 .or. j .eq. 1 .or. i .eq. user%mx .or. j .eq. user%my) then
461c4762a1bSJed Brown               f(i,j) = x(i,j)
462c4762a1bSJed Brown            else
463c4762a1bSJed Brown               u = x(i,j)
464c4762a1bSJed Brown               uxx = hydhx * (two*u - x(i-1,j) - x(i+1,j))
465c4762a1bSJed Brown               uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1))
466c4762a1bSJed Brown               f(i,j) = uxx + uyy - sc*exp(u)
467c4762a1bSJed Brown            endif
468c4762a1bSJed Brown 10      continue
469c4762a1bSJed Brown 20   continue
470c4762a1bSJed Brown      ierr = 0
471c4762a1bSJed Brown      return
472c4762a1bSJed Brown      end
473c4762a1bSJed Brown
474c4762a1bSJed Brown! ---------------------------------------------------------------------
475c4762a1bSJed Brown!
476c4762a1bSJed Brown!  FormJacobian - Evaluates Jacobian matrix.
477c4762a1bSJed Brown!
478c4762a1bSJed Brown!  Input Parameters:
479c4762a1bSJed Brown!  snes     - the SNES context
480c4762a1bSJed Brown!  x        - input vector
481c4762a1bSJed Brown!  dummy    - optional user-defined context, as set by SNESSetJacobian()
482c4762a1bSJed Brown!             (not used here)
483c4762a1bSJed Brown!
484c4762a1bSJed Brown!  Output Parameters:
485c4762a1bSJed Brown!  jac      - Jacobian matrix
486c4762a1bSJed Brown!  jac_prec - optionally different preconditioning matrix (not used here)
487c4762a1bSJed Brown!  flag     - flag indicating matrix structure
488c4762a1bSJed Brown!
489c4762a1bSJed Brown!  Notes:
490c4762a1bSJed Brown!  This routine serves as a wrapper for the lower-level routine
491c4762a1bSJed Brown!  "FormJacobianLocal", where the actual computations are
492c4762a1bSJed Brown!  done using the standard Fortran style of treating the local
493c4762a1bSJed Brown!  vector data as a multidimensional array over the local mesh.
494c4762a1bSJed Brown!  This routine merely accesses the local vector data via
495c4762a1bSJed Brown!  VecGetArrayF90() and VecRestoreArrayF90().
496c4762a1bSJed Brown!
497c4762a1bSJed Brown!  Notes:
498c4762a1bSJed Brown!  Due to grid point reordering with DMDAs, we must always work
499c4762a1bSJed Brown!  with the local grid points, and then transform them to the new
500c4762a1bSJed Brown!  global numbering with the "ltog" mapping
501c4762a1bSJed Brown!  We cannot work directly with the global numbers for the original
502c4762a1bSJed Brown!  uniprocessor grid!
503c4762a1bSJed Brown!
504c4762a1bSJed Brown!  Two methods are available for imposing this transformation
505c4762a1bSJed Brown!  when setting matrix entries:
506c4762a1bSJed Brown!    (A) MatSetValuesLocal(), using the local ordering (including
507c4762a1bSJed Brown!        ghost points!)
508c4762a1bSJed Brown!        - Set matrix entries using the local ordering
509c4762a1bSJed Brown!          by calling MatSetValuesLocal()
510c4762a1bSJed Brown!    (B) MatSetValues(), using the global ordering
511c4762a1bSJed Brown!        - Use DMGetLocalToGlobalMapping() then
512c4762a1bSJed Brown!          ISLocalToGlobalMappingGetIndicesF90() to extract the local-to-global map
513c4762a1bSJed Brown!        - Then apply this map explicitly yourself
514c4762a1bSJed Brown!        - Set matrix entries using the global ordering by calling
515c4762a1bSJed Brown!          MatSetValues()
516c4762a1bSJed Brown!  Option (A) seems cleaner/easier in many cases, and is the procedure
517c4762a1bSJed Brown!  used in this example.
518c4762a1bSJed Brown!
519c4762a1bSJed Brown      subroutine FormJacobian(mysnes,X,jac,jac_prec,user,ierr)
520c4762a1bSJed Brown#include <petsc/finclude/petscsnes.h>
521c4762a1bSJed Brown      use petscsnes
52277d968b7SBarry Smith      use ex5f90tmodule
523c4762a1bSJed Brown!  Input/output variables:
524c4762a1bSJed Brown      type(tSNES)     mysnes
525c4762a1bSJed Brown      type(tVec)      X
526c4762a1bSJed Brown      type(tMat)      jac,jac_prec
527c4762a1bSJed Brown      type(userctx)  user
528c4762a1bSJed Brown      PetscErrorCode ierr
529c4762a1bSJed Brown
530c4762a1bSJed Brown!  Declarations for use with local arrays:
531c4762a1bSJed Brown      PetscScalar,pointer :: lx_v(:)
532c4762a1bSJed Brown      type(tVec)      localX
533c4762a1bSJed Brown
534c4762a1bSJed Brown!  Scatter ghost points to local vector, using the 2-step process
535c4762a1bSJed Brown!     DMGlobalToLocalBegin(), DMGlobalToLocalEnd()
536c4762a1bSJed Brown!  Computations can be done while messages are in transition,
537c4762a1bSJed Brown!  by placing code between these two statements.
538c4762a1bSJed Brown
539d8606c27SBarry Smith      PetscCallA(DMGetLocalVector(user%da,localX,ierr))
540d8606c27SBarry Smith      PetscCallA(DMGlobalToLocalBegin(user%da,X,INSERT_VALUES,localX,ierr))
541d8606c27SBarry Smith      PetscCallA(DMGlobalToLocalEnd(user%da,X,INSERT_VALUES,localX,ierr))
542c4762a1bSJed Brown
543c4762a1bSJed Brown!  Get a pointer to vector data
544d8606c27SBarry Smith      PetscCallA(VecGetArrayF90(localX,lx_v,ierr))
545c4762a1bSJed Brown
546c4762a1bSJed Brown!  Compute entries for the locally owned part of the Jacobian preconditioner.
547d8606c27SBarry Smith      PetscCallA(FormJacobianLocal(lx_v,jac_prec,user,ierr))
548c4762a1bSJed Brown
549c4762a1bSJed Brown!  Assemble matrix, using the 2-step process:
550c4762a1bSJed Brown!     MatAssemblyBegin(), MatAssemblyEnd()
551c4762a1bSJed Brown!  Computations can be done while messages are in transition,
552c4762a1bSJed Brown!  by placing code between these two statements.
553c4762a1bSJed Brown
554d8606c27SBarry Smith      PetscCallA(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr))
555c4762a1bSJed Brown!      if (jac .ne. jac_prec) then
556d8606c27SBarry Smith         PetscCallA(MatAssemblyBegin(jac_prec,MAT_FINAL_ASSEMBLY,ierr))
557c4762a1bSJed Brown!      endif
558d8606c27SBarry Smith      PetscCallA(VecRestoreArrayF90(localX,lx_v,ierr))
559d8606c27SBarry Smith      PetscCallA(DMRestoreLocalVector(user%da,localX,ierr))
560d8606c27SBarry Smith      PetscCallA(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr))
561c4762a1bSJed Brown!      if (jac .ne. jac_prec) then
562d8606c27SBarry Smith        PetscCallA(MatAssemblyEnd(jac_prec,MAT_FINAL_ASSEMBLY,ierr))
563c4762a1bSJed Brown!      endif
564c4762a1bSJed Brown
565c4762a1bSJed Brown!  Tell the matrix we will never add a new nonzero location to the
566c4762a1bSJed Brown!  matrix. If we do it will generate an error.
567c4762a1bSJed Brown
568d8606c27SBarry Smith      PetscCallA(MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE,ierr))
569c4762a1bSJed Brown
570c4762a1bSJed Brown      return
571c4762a1bSJed Brown      end
572c4762a1bSJed Brown
573c4762a1bSJed Brown! ---------------------------------------------------------------------
574c4762a1bSJed Brown!
575c4762a1bSJed Brown!  FormJacobianLocal - Computes Jacobian preconditioner matrix,
576c4762a1bSJed Brown!  called by the higher level routine FormJacobian().
577c4762a1bSJed Brown!
578c4762a1bSJed Brown!  Input Parameters:
579c4762a1bSJed Brown!  x        - local vector data
580c4762a1bSJed Brown!
581c4762a1bSJed Brown!  Output Parameters:
582c4762a1bSJed Brown!  jac_prec - Jacobian preconditioner matrix
583c4762a1bSJed Brown!  ierr     - error code
584c4762a1bSJed Brown!
585c4762a1bSJed Brown!  Notes:
586c4762a1bSJed Brown!  This routine uses standard Fortran-style computations over a 2-dim array.
587c4762a1bSJed Brown!
588c4762a1bSJed Brown!  Notes:
589c4762a1bSJed Brown!  Due to grid point reordering with DMDAs, we must always work
590c4762a1bSJed Brown!  with the local grid points, and then transform them to the new
591c4762a1bSJed Brown!  global numbering with the "ltog" mapping
592c4762a1bSJed Brown!  We cannot work directly with the global numbers for the original
593c4762a1bSJed Brown!  uniprocessor grid!
594c4762a1bSJed Brown!
595c4762a1bSJed Brown!  Two methods are available for imposing this transformation
596c4762a1bSJed Brown!  when setting matrix entries:
597c4762a1bSJed Brown!    (A) MatSetValuesLocal(), using the local ordering (including
598c4762a1bSJed Brown!        ghost points!)
599c4762a1bSJed Brown!        - Set matrix entries using the local ordering
600c4762a1bSJed Brown!          by calling MatSetValuesLocal()
601c4762a1bSJed Brown!    (B) MatSetValues(), using the global ordering
602c4762a1bSJed Brown!        - Set matrix entries using the global ordering by calling
603c4762a1bSJed Brown!          MatSetValues()
604c4762a1bSJed Brown!  Option (A) seems cleaner/easier in many cases, and is the procedure
605c4762a1bSJed Brown!  used in this example.
606c4762a1bSJed Brown!
607c4762a1bSJed Brown      subroutine FormJacobianLocal(x,jac_prec,user,ierr)
608c4762a1bSJed Brown#include <petsc/finclude/petscmat.h>
609c4762a1bSJed Brown      use petscmat
61077d968b7SBarry Smith      use ex5f90tmodule
611c4762a1bSJed Brown!  Input/output variables:
612c4762a1bSJed Brown      type (userctx) user
613c4762a1bSJed Brown      PetscScalar    x(user%gxs:user%gxe,user%gys:user%gye)
614c4762a1bSJed Brown      type(tMat)      jac_prec
615c4762a1bSJed Brown      PetscErrorCode ierr
616c4762a1bSJed Brown
617c4762a1bSJed Brown!  Local variables:
618c4762a1bSJed Brown      PetscInt    row,col(5),i,j
619c4762a1bSJed Brown      PetscInt    ione,ifive
620c4762a1bSJed Brown      PetscScalar two,one,hx,hy,hxdhy
621c4762a1bSJed Brown      PetscScalar hydhx,sc,v(5)
622c4762a1bSJed Brown
623c4762a1bSJed Brown!  Set parameters
624c4762a1bSJed Brown      ione   = 1
625c4762a1bSJed Brown      ifive  = 5
626c4762a1bSJed Brown      one    = 1.0
627c4762a1bSJed Brown      two    = 2.0
628c4762a1bSJed Brown      hx     = one/PetscIntToReal(user%mx-1)
629c4762a1bSJed Brown      hy     = one/PetscIntToReal(user%my-1)
630c4762a1bSJed Brown      sc     = hx*hy
631c4762a1bSJed Brown      hxdhy  = hx/hy
632c4762a1bSJed Brown      hydhx  = hy/hx
633c4762a1bSJed Brown
634c4762a1bSJed Brown!  Compute entries for the locally owned part of the Jacobian.
635c4762a1bSJed Brown!   - Currently, all PETSc parallel matrix formats are partitioned by
636c4762a1bSJed Brown!     contiguous chunks of rows across the processors.
637c4762a1bSJed Brown!   - Each processor needs to insert only elements that it owns
638c4762a1bSJed Brown!     locally (but any non-local elements will be sent to the
639c4762a1bSJed Brown!     appropriate processor during matrix assembly).
640c4762a1bSJed Brown!   - Here, we set all entries for a particular row at once.
641c4762a1bSJed Brown!   - We can set matrix entries either using either
642c4762a1bSJed Brown!     MatSetValuesLocal() or MatSetValues(), as discussed above.
643c4762a1bSJed Brown!   - Note that MatSetValues() uses 0-based row and column numbers
644c4762a1bSJed Brown!     in Fortran as well as in C.
645c4762a1bSJed Brown
646c4762a1bSJed Brown      do 20 j=user%ys,user%ye
647c4762a1bSJed Brown         row = (j - user%gys)*user%gxm + user%xs - user%gxs - 1
648c4762a1bSJed Brown         do 10 i=user%xs,user%xe
649c4762a1bSJed Brown            row = row + 1
650c4762a1bSJed Brown!           boundary points
651c4762a1bSJed Brown            if (i .eq. 1 .or. j .eq. 1 .or. i .eq. user%mx .or. j .eq. user%my) then
652c4762a1bSJed Brown               col(1) = row
653c4762a1bSJed Brown               v(1)   = one
654d8606c27SBarry Smith               PetscCallA(MatSetValuesLocal(jac_prec,ione,row,ione,col,v,INSERT_VALUES,ierr))
655c4762a1bSJed Brown!           interior grid points
656c4762a1bSJed Brown            else
657c4762a1bSJed Brown               v(1) = -hxdhy
658c4762a1bSJed Brown               v(2) = -hydhx
659c4762a1bSJed Brown               v(3) = two*(hydhx + hxdhy) - sc*user%lambda*exp(x(i,j))
660c4762a1bSJed Brown               v(4) = -hydhx
661c4762a1bSJed Brown               v(5) = -hxdhy
662c4762a1bSJed Brown               col(1) = row - user%gxm
663c4762a1bSJed Brown               col(2) = row - 1
664c4762a1bSJed Brown               col(3) = row
665c4762a1bSJed Brown               col(4) = row + 1
666c4762a1bSJed Brown               col(5) = row + user%gxm
667d8606c27SBarry Smith               PetscCallA(MatSetValuesLocal(jac_prec,ione,row,ifive,col,v,INSERT_VALUES,ierr))
668c4762a1bSJed Brown            endif
669c4762a1bSJed Brown 10      continue
670c4762a1bSJed Brown 20   continue
671c4762a1bSJed Brown      return
672c4762a1bSJed Brown      end
673c4762a1bSJed Brown
674c4762a1bSJed Brown!/*TEST
675c4762a1bSJed Brown!
676c4762a1bSJed Brown!   test:
677c4762a1bSJed Brown!      nsize: 4
678c4762a1bSJed Brown!      args: -snes_mf -pc_type none -da_processors_x 4 -da_processors_y 1 -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
679c4762a1bSJed Brown!
680c4762a1bSJed Brown!TEST*/
681