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