xref: /petsc/src/snes/tutorials/ex5f90.F90 (revision dfbbaf821b4c49d07b4ce746493b0d955783fdf9)
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
34*dfbbaf82SBarry Smith      module ex5f90module
35c4762a1bSJed Brown      use petscsnes
36*dfbbaf82SBarry Smith      use petscdmda
37c4762a1bSJed Brown#include <petsc/finclude/petscsnes.h>
38c4762a1bSJed Brown      type userctx
39c4762a1bSJed Brown        PetscInt xs,xe,xm,gxs,gxe,gxm
40c4762a1bSJed Brown        PetscInt ys,ye,ym,gys,gye,gym
41c4762a1bSJed Brown        PetscInt mx,my
42c4762a1bSJed Brown        PetscMPIInt rank
43c4762a1bSJed Brown        PetscReal lambda
44c4762a1bSJed Brown      end type userctx
45c4762a1bSJed Brown
46c4762a1bSJed Brown      contains
47c4762a1bSJed Brown! ---------------------------------------------------------------------
48c4762a1bSJed Brown!
49c4762a1bSJed Brown!  FormFunction - Evaluates nonlinear function, F(x).
50c4762a1bSJed Brown!
51c4762a1bSJed Brown!  Input Parameters:
52c4762a1bSJed Brown!  snes - the SNES context
53c4762a1bSJed Brown!  X - input vector
54c4762a1bSJed Brown!  dummy - optional user-defined context, as set by SNESSetFunction()
55c4762a1bSJed Brown!          (not used here)
56c4762a1bSJed Brown!
57c4762a1bSJed Brown!  Output Parameter:
58c4762a1bSJed Brown!  F - function vector
59c4762a1bSJed Brown!
60c4762a1bSJed Brown!  Notes:
61c4762a1bSJed Brown!  This routine serves as a wrapper for the lower-level routine
62c4762a1bSJed Brown!  "FormFunctionLocal", where the actual computations are
63c4762a1bSJed Brown!  done using the standard Fortran style of treating the local
64c4762a1bSJed Brown!  vector data as a multidimensional array over the local mesh.
65c4762a1bSJed Brown!  This routine merely handles ghost point scatters and accesses
66c4762a1bSJed Brown!  the local vector data via VecGetArrayF90() and VecRestoreArrayF90().
67c4762a1bSJed Brown!
68c4762a1bSJed Brown      subroutine FormFunction(snes,X,F,user,ierr)
69c4762a1bSJed Brown      implicit none
70c4762a1bSJed Brown
71c4762a1bSJed Brown!  Input/output variables:
72c4762a1bSJed Brown      SNES           snes
73c4762a1bSJed Brown      Vec            X,F
74c4762a1bSJed Brown      PetscErrorCode ierr
75c4762a1bSJed Brown      type (userctx) user
76c4762a1bSJed Brown      DM             da
77c4762a1bSJed Brown
78c4762a1bSJed Brown!  Declarations for use with local arrays:
79c4762a1bSJed Brown      PetscScalar,pointer :: lx_v(:),lf_v(:)
80c4762a1bSJed Brown      Vec            localX
81c4762a1bSJed Brown
82c4762a1bSJed Brown!  Scatter ghost points to local vector, using the 2-step process
83c4762a1bSJed Brown!     DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
84c4762a1bSJed Brown!  By placing code between these two statements, computations can
85c4762a1bSJed Brown!  be done while messages are in transition.
86d8606c27SBarry Smith      PetscCall(SNESGetDM(snes,da,ierr))
87d8606c27SBarry Smith      PetscCall(DMGetLocalVector(da,localX,ierr))
88d8606c27SBarry Smith      PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX,ierr))
89d8606c27SBarry Smith      PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX,ierr))
90c4762a1bSJed Brown
91c4762a1bSJed Brown!  Get a pointer to vector data.
92c4762a1bSJed Brown!    - For default PETSc vectors, VecGetArray90() returns a pointer to
93c4762a1bSJed Brown!      the data array. Otherwise, the routine is implementation dependent.
94c4762a1bSJed Brown!    - You MUST call VecRestoreArrayF90() when you no longer need access to
95c4762a1bSJed Brown!      the array.
96c4762a1bSJed Brown!    - Note that the interface to VecGetArrayF90() differs from VecGetArray(),
97c4762a1bSJed Brown!      and is useable from Fortran-90 Only.
98c4762a1bSJed Brown
99d8606c27SBarry Smith      PetscCall(VecGetArrayF90(localX,lx_v,ierr))
100d8606c27SBarry Smith      PetscCall(VecGetArrayF90(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
106d8606c27SBarry Smith      PetscCall(VecRestoreArrayF90(localX,lx_v,ierr))
107d8606c27SBarry Smith      PetscCall(VecRestoreArrayF90(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      return
117c4762a1bSJed Brown      end subroutine formfunction
118*dfbbaf82SBarry Smith      end module ex5f90module
119c4762a1bSJed Brown
120*dfbbaf82SBarry Smith      module ex5f90moduleinterfaces
121*dfbbaf82SBarry Smith        use ex5f90module
122c4762a1bSJed Brown
123c4762a1bSJed Brown      Interface SNESSetApplicationContext
124c4762a1bSJed Brown        Subroutine SNESSetApplicationContext(snes,ctx,ierr)
125*dfbbaf82SBarry Smith        use ex5f90module
126c4762a1bSJed Brown          SNES snes
127c4762a1bSJed Brown          type(userctx) ctx
128c4762a1bSJed Brown          PetscErrorCode ierr
129c4762a1bSJed Brown        End Subroutine
130c4762a1bSJed Brown      End Interface SNESSetApplicationContext
131c4762a1bSJed Brown
132c4762a1bSJed Brown      Interface SNESGetApplicationContext
133c4762a1bSJed Brown        Subroutine SNESGetApplicationContext(snes,ctx,ierr)
134*dfbbaf82SBarry Smith        use ex5f90module
135c4762a1bSJed Brown          SNES snes
136c4762a1bSJed Brown          type(userctx), pointer :: ctx
137c4762a1bSJed Brown          PetscErrorCode ierr
138c4762a1bSJed Brown        End Subroutine
139c4762a1bSJed Brown      End Interface SNESGetApplicationContext
140*dfbbaf82SBarry Smith      end module ex5f90moduleinterfaces
141c4762a1bSJed Brown
142c4762a1bSJed Brown      program main
143*dfbbaf82SBarry Smith      use ex5f90module
144*dfbbaf82SBarry Smith      use ex5f90moduleinterfaces
145c4762a1bSJed Brown      implicit none
146c4762a1bSJed Brown!
147c4762a1bSJed Brown
148c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149c4762a1bSJed Brown!                   Variable declarations
150c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151c4762a1bSJed Brown!
152c4762a1bSJed Brown!  Variables:
153c4762a1bSJed Brown!     snes        - nonlinear solver
154c4762a1bSJed Brown!     x, r        - solution, residual vectors
155c4762a1bSJed Brown!     J           - Jacobian matrix
156c4762a1bSJed Brown!     its         - iterations for convergence
157c4762a1bSJed Brown!     Nx, Ny      - number of preocessors in x- and y- directions
158c4762a1bSJed Brown!     matrix_free - flag - 1 indicates matrix-free version
159c4762a1bSJed Brown!
160c4762a1bSJed Brown      SNES             snes
161c4762a1bSJed Brown      Vec              x,r
162c4762a1bSJed Brown      Mat              J
163c4762a1bSJed Brown      PetscErrorCode   ierr
164c4762a1bSJed Brown      PetscInt         its
165c4762a1bSJed Brown      PetscBool        flg,matrix_free
166c4762a1bSJed Brown      PetscInt         ione,nfour
167c4762a1bSJed Brown      PetscReal lambda_max,lambda_min
168c4762a1bSJed Brown      type (userctx)   user
169c4762a1bSJed Brown      DM               da
170c4762a1bSJed Brown
171c4762a1bSJed Brown!  Note: Any user-defined Fortran routines (such as FormJacobian)
172c4762a1bSJed Brown!  MUST be declared as external.
173c4762a1bSJed Brown      external FormInitialGuess,FormJacobian
174c4762a1bSJed Brown
175c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
176c4762a1bSJed Brown!  Initialize program
177c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178d8606c27SBarry Smith      PetscCallA(PetscInitialize(ierr))
179d8606c27SBarry Smith      PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,user%rank,ierr))
180c4762a1bSJed Brown
181c4762a1bSJed Brown!  Initialize problem parameters
182c4762a1bSJed Brown      lambda_max  = 6.81
183c4762a1bSJed Brown      lambda_min  = 0.0
184c4762a1bSJed Brown      user%lambda = 6.0
185c4762a1bSJed Brown      ione = 1
186c4762a1bSJed Brown      nfour = 4
187d8606c27SBarry Smith      PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-par',user%lambda,flg,ierr))
188d8606c27SBarry Smith      if (user%lambda .ge. lambda_max .or. user%lambda .le. lambda_min) then
189d8606c27SBarry Smith         SETERRA(PETSC_COMM_SELF,PETSC_ERR_USER,'Lambda provided with -par is out of range ')
190d8606c27SBarry Smith      endif
191c4762a1bSJed Brown
192c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
193c4762a1bSJed Brown!  Create nonlinear solver context
194c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
195d8606c27SBarry Smith      PetscCallA(SNESCreate(PETSC_COMM_WORLD,snes,ierr))
196c4762a1bSJed Brown
197c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198c4762a1bSJed Brown!  Create vector data structures; set function evaluation routine
199c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
200c4762a1bSJed Brown
201c4762a1bSJed Brown!  Create distributed array (DMDA) to manage parallel grid and vectors
202c4762a1bSJed Brown
203c4762a1bSJed Brown! This really needs only the star-type stencil, but we use the box
204c4762a1bSJed Brown! stencil temporarily.
205d8606c27SBarry 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,da,ierr))
206d8606c27SBarry Smith      PetscCallA(DMSetFromOptions(da,ierr))
207d8606c27SBarry Smith      PetscCallA(DMSetUp(da,ierr))
208c4762a1bSJed Brown
209d8606c27SBarry 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_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ierr))
210c4762a1bSJed Brown
211c4762a1bSJed Brown!
212c4762a1bSJed Brown!   Visualize the distribution of the array across the processors
213c4762a1bSJed Brown!
214d8606c27SBarry Smith!     PetscCallA(DMView(da,PETSC_VIEWER_DRAW_WORLD,ierr))
215c4762a1bSJed Brown
216c4762a1bSJed Brown!  Extract global and local vectors from DMDA; then duplicate for remaining
217c4762a1bSJed Brown!  vectors that are the same types
218d8606c27SBarry Smith      PetscCallA(DMCreateGlobalVector(da,x,ierr))
219d8606c27SBarry Smith      PetscCallA(VecDuplicate(x,r,ierr))
220c4762a1bSJed Brown
221c4762a1bSJed Brown!  Get local grid boundaries (for 2-dimensional DMDA)
222d8606c27SBarry Smith      PetscCallA(DMDAGetCorners(da,user%xs,user%ys,PETSC_NULL_INTEGER,user%xm,user%ym,PETSC_NULL_INTEGER,ierr))
223d8606c27SBarry Smith      PetscCallA(DMDAGetGhostCorners(da,user%gxs,user%gys,PETSC_NULL_INTEGER,user%gxm,user%gym,PETSC_NULL_INTEGER,ierr))
224c4762a1bSJed Brown
225c4762a1bSJed Brown!  Here we shift the starting indices up by one so that we can easily
226c4762a1bSJed Brown!  use the Fortran convention of 1-based indices (rather 0-based indices).
227c4762a1bSJed Brown      user%xs  = user%xs+1
228c4762a1bSJed Brown      user%ys  = user%ys+1
229c4762a1bSJed Brown      user%gxs = user%gxs+1
230c4762a1bSJed Brown      user%gys = user%gys+1
231c4762a1bSJed Brown
232c4762a1bSJed Brown      user%ye  = user%ys+user%ym-1
233c4762a1bSJed Brown      user%xe  = user%xs+user%xm-1
234c4762a1bSJed Brown      user%gye = user%gys+user%gym-1
235c4762a1bSJed Brown      user%gxe = user%gxs+user%gxm-1
236c4762a1bSJed Brown
237d8606c27SBarry Smith      PetscCallA(SNESSetApplicationContext(snes,user,ierr))
238c4762a1bSJed Brown
239c4762a1bSJed Brown!  Set function evaluation routine and vector
240d8606c27SBarry Smith      PetscCallA(SNESSetFunction(snes,r,FormFunction,user,ierr))
241c4762a1bSJed Brown
242c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
243c4762a1bSJed Brown!  Create matrix data structure; set Jacobian evaluation routine
244c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
245c4762a1bSJed Brown
246c4762a1bSJed Brown!  Set Jacobian matrix data structure and default Jacobian evaluation
247c4762a1bSJed Brown!  routine. User can override with:
248c4762a1bSJed Brown!     -snes_fd : default finite differencing approximation of Jacobian
249c4762a1bSJed Brown!     -snes_mf : matrix-free Newton-Krylov method with no preconditioning
250c4762a1bSJed Brown!                (unless user explicitly sets preconditioner)
251c4762a1bSJed Brown!     -snes_mf_operator : form preconditioning matrix as set by the user,
252c4762a1bSJed Brown!                         but use matrix-free approx for Jacobian-vector
253c4762a1bSJed Brown!                         products within Newton-Krylov method
254c4762a1bSJed Brown!
255c4762a1bSJed Brown!  Note:  For the parallel case, vectors and matrices MUST be partitioned
256c4762a1bSJed Brown!     accordingly.  When using distributed arrays (DMDAs) to create vectors,
257c4762a1bSJed Brown!     the DMDAs determine the problem partitioning.  We must explicitly
258c4762a1bSJed Brown!     specify the local matrix dimensions upon its creation for compatibility
259c4762a1bSJed Brown!     with the vector distribution.  Thus, the generic MatCreate() routine
260c4762a1bSJed Brown!     is NOT sufficient when working with distributed arrays.
261c4762a1bSJed Brown!
262c4762a1bSJed Brown!     Note: Here we only approximately preallocate storage space for the
263c4762a1bSJed Brown!     Jacobian.  See the users manual for a discussion of better techniques
264c4762a1bSJed Brown!     for preallocating matrix memory.
265c4762a1bSJed Brown
266d8606c27SBarry Smith      PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-snes_mf',matrix_free,ierr))
267c4762a1bSJed Brown      if (.not. matrix_free) then
268d8606c27SBarry Smith        PetscCallA(DMSetMatType(da,MATAIJ,ierr))
269d8606c27SBarry Smith        PetscCallA(DMCreateMatrix(da,J,ierr))
270d8606c27SBarry Smith        PetscCallA(SNESSetJacobian(snes,J,J,FormJacobian,user,ierr))
271c4762a1bSJed Brown      endif
272c4762a1bSJed Brown
273c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
274c4762a1bSJed Brown!  Customize nonlinear solver; set runtime options
275c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
276c4762a1bSJed Brown!  Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
277d8606c27SBarry Smith      PetscCallA(SNESSetDM(snes,da,ierr))
278d8606c27SBarry Smith      PetscCallA(SNESSetFromOptions(snes,ierr))
279c4762a1bSJed Brown
280c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
281c4762a1bSJed Brown!  Evaluate initial guess; then solve nonlinear system.
282c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
283c4762a1bSJed Brown!  Note: The user should initialize the vector, x, with the initial guess
284c4762a1bSJed Brown!  for the nonlinear solver prior to calling SNESSolve().  In particular,
285c4762a1bSJed Brown!  to employ an initial guess of zero, the user should explicitly set
286c4762a1bSJed Brown!  this vector to zero by calling VecSet().
287c4762a1bSJed Brown
288d8606c27SBarry Smith      PetscCallA(FormInitialGuess(snes,x,ierr))
289d8606c27SBarry Smith      PetscCallA(SNESSolve(snes,PETSC_NULL_VEC,x,ierr))
290d8606c27SBarry Smith      PetscCallA(SNESGetIterationNumber(snes,its,ierr))
291c4762a1bSJed Brown      if (user%rank .eq. 0) then
292c4762a1bSJed Brown         write(6,100) its
293c4762a1bSJed Brown      endif
294c4762a1bSJed Brown  100 format('Number of SNES iterations = ',i5)
295c4762a1bSJed Brown
296c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
297c4762a1bSJed Brown!  Free work space.  All PETSc objects should be destroyed when they
298c4762a1bSJed Brown!  are no longer needed.
299c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
300d8606c27SBarry Smith      if (.not. matrix_free) PetscCallA(MatDestroy(J,ierr))
301d8606c27SBarry Smith      PetscCallA(VecDestroy(x,ierr))
302d8606c27SBarry Smith      PetscCallA(VecDestroy(r,ierr))
303d8606c27SBarry Smith      PetscCallA(SNESDestroy(snes,ierr))
304d8606c27SBarry Smith      PetscCallA(DMDestroy(da,ierr))
305c4762a1bSJed Brown
306d8606c27SBarry Smith      PetscCallA(PetscFinalize(ierr))
307c4762a1bSJed Brown      end
308c4762a1bSJed Brown
309c4762a1bSJed Brown! ---------------------------------------------------------------------
310c4762a1bSJed Brown!
311c4762a1bSJed Brown!  FormInitialGuess - Forms initial approximation.
312c4762a1bSJed Brown!
313c4762a1bSJed Brown!  Input Parameters:
314c4762a1bSJed Brown!  X - vector
315c4762a1bSJed Brown!
316c4762a1bSJed Brown!  Output Parameter:
317c4762a1bSJed Brown!  X - vector
318c4762a1bSJed Brown!
319c4762a1bSJed Brown!  Notes:
320c4762a1bSJed Brown!  This routine serves as a wrapper for the lower-level routine
321c4762a1bSJed Brown!  "InitialGuessLocal", where the actual computations are
322c4762a1bSJed Brown!  done using the standard Fortran style of treating the local
323c4762a1bSJed Brown!  vector data as a multidimensional array over the local mesh.
324c4762a1bSJed Brown!  This routine merely handles ghost point scatters and accesses
325c4762a1bSJed Brown!  the local vector data via VecGetArrayF90() and VecRestoreArrayF90().
326c4762a1bSJed Brown!
327c4762a1bSJed Brown      subroutine FormInitialGuess(snes,X,ierr)
328*dfbbaf82SBarry Smith      use ex5f90module
329*dfbbaf82SBarry Smith      use ex5f90moduleinterfaces
330c4762a1bSJed Brown      implicit none
331c4762a1bSJed Brown
332c4762a1bSJed Brown!  Input/output variables:
333c4762a1bSJed Brown      SNES           snes
334c4762a1bSJed Brown      type(userctx), pointer:: puser
335c4762a1bSJed Brown      Vec            X
336c4762a1bSJed Brown      PetscErrorCode ierr
337c4762a1bSJed Brown      DM             da
338c4762a1bSJed Brown
339c4762a1bSJed Brown!  Declarations for use with local arrays:
340c4762a1bSJed Brown      PetscScalar,pointer :: lx_v(:)
341c4762a1bSJed Brown
342c4762a1bSJed Brown      ierr = 0
343d8606c27SBarry Smith      PetscCallA(SNESGetDM(snes,da,ierr))
344d8606c27SBarry Smith      PetscCallA(SNESGetApplicationContext(snes,puser,ierr))
345c4762a1bSJed Brown!  Get a pointer to vector data.
346c4762a1bSJed Brown!    - For default PETSc vectors, VecGetArray90() returns a pointer to
347c4762a1bSJed Brown!      the data array. Otherwise, the routine is implementation dependent.
348c4762a1bSJed Brown!    - You MUST call VecRestoreArrayF90() when you no longer need access to
349c4762a1bSJed Brown!      the array.
350c4762a1bSJed Brown!    - Note that the interface to VecGetArrayF90() differs from VecGetArray(),
351c4762a1bSJed Brown!      and is useable from Fortran-90 Only.
352c4762a1bSJed Brown
353d8606c27SBarry Smith      PetscCallA(VecGetArrayF90(X,lx_v,ierr))
354c4762a1bSJed Brown
355c4762a1bSJed Brown!  Compute initial guess over the locally owned part of the grid
356d8606c27SBarry Smith      PetscCallA(InitialGuessLocal(puser,lx_v,ierr))
357c4762a1bSJed Brown
358c4762a1bSJed Brown!  Restore vector
359d8606c27SBarry Smith      PetscCallA(VecRestoreArrayF90(X,lx_v,ierr))
360c4762a1bSJed Brown
361c4762a1bSJed Brown!  Insert values into global vector
362c4762a1bSJed Brown
363c4762a1bSJed Brown      return
364c4762a1bSJed Brown      end
365c4762a1bSJed Brown
366c4762a1bSJed Brown! ---------------------------------------------------------------------
367c4762a1bSJed Brown!
368c4762a1bSJed Brown!  InitialGuessLocal - Computes initial approximation, called by
369c4762a1bSJed Brown!  the higher level routine FormInitialGuess().
370c4762a1bSJed Brown!
371c4762a1bSJed Brown!  Input Parameter:
372c4762a1bSJed Brown!  x - local vector data
373c4762a1bSJed Brown!
374c4762a1bSJed Brown!  Output Parameters:
375c4762a1bSJed Brown!  x - local vector data
376c4762a1bSJed Brown!  ierr - error code
377c4762a1bSJed Brown!
378c4762a1bSJed Brown!  Notes:
379c4762a1bSJed Brown!  This routine uses standard Fortran-style computations over a 2-dim array.
380c4762a1bSJed Brown!
381c4762a1bSJed Brown      subroutine InitialGuessLocal(user,x,ierr)
382*dfbbaf82SBarry Smith      use ex5f90module
383c4762a1bSJed Brown      implicit none
384c4762a1bSJed Brown
385c4762a1bSJed Brown!  Input/output variables:
386c4762a1bSJed Brown      type (userctx)         user
387c4762a1bSJed Brown      PetscScalar  x(user%xs:user%xe,user%ys:user%ye)
388c4762a1bSJed Brown      PetscErrorCode ierr
389c4762a1bSJed Brown
390c4762a1bSJed Brown!  Local variables:
391c4762a1bSJed Brown      PetscInt  i,j
392c4762a1bSJed Brown      PetscReal   temp1,temp,hx,hy
393c4762a1bSJed Brown      PetscReal   one
394c4762a1bSJed Brown
395c4762a1bSJed Brown!  Set parameters
396c4762a1bSJed Brown
397c4762a1bSJed Brown      ierr   = 0
398c4762a1bSJed Brown      one    = 1.0
399c4762a1bSJed Brown      hx     = one/(user%mx-1)
400c4762a1bSJed Brown      hy     = one/(user%my-1)
401c4762a1bSJed Brown      temp1  = user%lambda/(user%lambda + one)
402c4762a1bSJed Brown
403c4762a1bSJed Brown      do 20 j=user%ys,user%ye
404c4762a1bSJed Brown         temp = min(j-1,user%my-j)*hy
405c4762a1bSJed Brown         do 10 i=user%xs,user%xe
406c4762a1bSJed Brown            if (i .eq. 1 .or. j .eq. 1 .or. i .eq. user%mx .or. j .eq. user%my) then
407c4762a1bSJed Brown              x(i,j) = 0.0
408c4762a1bSJed Brown            else
409c4762a1bSJed Brown              x(i,j) = temp1 * sqrt(min(hx*min(i-1,user%mx-i),temp))
410c4762a1bSJed Brown            endif
411c4762a1bSJed Brown 10      continue
412c4762a1bSJed Brown 20   continue
413c4762a1bSJed Brown
414c4762a1bSJed Brown      return
415c4762a1bSJed Brown      end
416c4762a1bSJed Brown
417c4762a1bSJed Brown! ---------------------------------------------------------------------
418c4762a1bSJed Brown!
419c4762a1bSJed Brown!  FormFunctionLocal - Computes nonlinear function, called by
420c4762a1bSJed Brown!  the higher level routine FormFunction().
421c4762a1bSJed Brown!
422c4762a1bSJed Brown!  Input Parameter:
423c4762a1bSJed Brown!  x - local vector data
424c4762a1bSJed Brown!
425c4762a1bSJed Brown!  Output Parameters:
426c4762a1bSJed Brown!  f - local vector data, f(x)
427c4762a1bSJed Brown!  ierr - error code
428c4762a1bSJed Brown!
429c4762a1bSJed Brown!  Notes:
430c4762a1bSJed Brown!  This routine uses standard Fortran-style computations over a 2-dim array.
431c4762a1bSJed Brown!
432c4762a1bSJed Brown      subroutine FormFunctionLocal(x,f,user,ierr)
433*dfbbaf82SBarry Smith      use ex5f90module
434c4762a1bSJed Brown
435c4762a1bSJed Brown      implicit none
436c4762a1bSJed Brown
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/(user%mx-1)
451c4762a1bSJed Brown      hy     = one/(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
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
512c4762a1bSJed Brown!        - Set matrix entries using the global ordering by calling
513c4762a1bSJed Brown!          MatSetValues()
514c4762a1bSJed Brown!  Option (A) seems cleaner/easier in many cases, and is the procedure
515c4762a1bSJed Brown!  used in this example.
516c4762a1bSJed Brown!
517c4762a1bSJed Brown      subroutine FormJacobian(snes,X,jac,jac_prec,user,ierr)
518*dfbbaf82SBarry Smith      use ex5f90module
519c4762a1bSJed Brown      implicit none
520c4762a1bSJed Brown
521c4762a1bSJed Brown!  Input/output variables:
522c4762a1bSJed Brown      SNES         snes
523c4762a1bSJed Brown      Vec          X
524c4762a1bSJed Brown      Mat          jac,jac_prec
525c4762a1bSJed Brown      type(userctx)  user
526c4762a1bSJed Brown      PetscErrorCode ierr
527c4762a1bSJed Brown      DM             da
528c4762a1bSJed Brown
529c4762a1bSJed Brown!  Declarations for use with local arrays:
530c4762a1bSJed Brown      PetscScalar,pointer :: lx_v(:)
531c4762a1bSJed Brown      Vec            localX
532c4762a1bSJed Brown
533c4762a1bSJed Brown!  Scatter ghost points to local vector, using the 2-step process
534c4762a1bSJed Brown!     DMGlobalToLocalBegin(), DMGlobalToLocalEnd()
535c4762a1bSJed Brown!  Computations can be done while messages are in transition,
536c4762a1bSJed Brown!  by placing code between these two statements.
537c4762a1bSJed Brown
538d8606c27SBarry Smith      PetscCallA(SNESGetDM(snes,da,ierr))
539d8606c27SBarry Smith      PetscCallA(DMGetLocalVector(da,localX,ierr))
540d8606c27SBarry Smith      PetscCallA(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX,ierr))
541d8606c27SBarry Smith      PetscCallA(DMGlobalToLocalEnd(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(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!        - Then apply this map explicitly yourself
603c4762a1bSJed Brown!        - Set matrix entries using the global ordering by calling
604c4762a1bSJed Brown!          MatSetValues()
605c4762a1bSJed Brown!  Option (A) seems cleaner/easier in many cases, and is the procedure
606c4762a1bSJed Brown!  used in this example.
607c4762a1bSJed Brown!
608c4762a1bSJed Brown      subroutine FormJacobianLocal(x,jac_prec,user,ierr)
609*dfbbaf82SBarry Smith      use ex5f90module
610c4762a1bSJed Brown      implicit none
611c4762a1bSJed Brown
612c4762a1bSJed Brown!  Input/output variables:
613c4762a1bSJed Brown      type (userctx) user
614c4762a1bSJed Brown      PetscScalar    x(user%gxs:user%gxe,user%gys:user%gye)
615c4762a1bSJed Brown      Mat            jac_prec
616c4762a1bSJed Brown      PetscErrorCode ierr
617c4762a1bSJed Brown
618c4762a1bSJed Brown!  Local variables:
619c4762a1bSJed Brown      PetscInt    row,col(5),i,j
620c4762a1bSJed Brown      PetscInt    ione,ifive
621c4762a1bSJed Brown      PetscScalar two,one,hx,hy,hxdhy
622c4762a1bSJed Brown      PetscScalar hydhx,sc,v(5)
623c4762a1bSJed Brown
624c4762a1bSJed Brown!  Set parameters
625c4762a1bSJed Brown      ione   = 1
626c4762a1bSJed Brown      ifive  = 5
627c4762a1bSJed Brown      one    = 1.0
628c4762a1bSJed Brown      two    = 2.0
629c4762a1bSJed Brown      hx     = one/(user%mx-1)
630c4762a1bSJed Brown      hy     = one/(user%my-1)
631c4762a1bSJed Brown      sc     = hx*hy
632c4762a1bSJed Brown      hxdhy  = hx/hy
633c4762a1bSJed Brown      hydhx  = hy/hx
634c4762a1bSJed Brown
635c4762a1bSJed Brown!  Compute entries for the locally owned part of the Jacobian.
636c4762a1bSJed Brown!   - Currently, all PETSc parallel matrix formats are partitioned by
637c4762a1bSJed Brown!     contiguous chunks of rows across the processors.
638c4762a1bSJed Brown!   - Each processor needs to insert only elements that it owns
639c4762a1bSJed Brown!     locally (but any non-local elements will be sent to the
640c4762a1bSJed Brown!     appropriate processor during matrix assembly).
641c4762a1bSJed Brown!   - Here, we set all entries for a particular row at once.
642c4762a1bSJed Brown!   - We can set matrix entries either using either
643c4762a1bSJed Brown!     MatSetValuesLocal() or MatSetValues(), as discussed above.
644c4762a1bSJed Brown!   - Note that MatSetValues() uses 0-based row and column numbers
645c4762a1bSJed Brown!     in Fortran as well as in C.
646c4762a1bSJed Brown
647c4762a1bSJed Brown      do 20 j=user%ys,user%ye
648c4762a1bSJed Brown         row = (j - user%gys)*user%gxm + user%xs - user%gxs - 1
649c4762a1bSJed Brown         do 10 i=user%xs,user%xe
650c4762a1bSJed Brown            row = row + 1
651c4762a1bSJed Brown!           boundary points
652c4762a1bSJed Brown            if (i .eq. 1 .or. j .eq. 1 .or. i .eq. user%mx .or. j .eq. user%my) then
653c4762a1bSJed Brown               col(1) = row
654c4762a1bSJed Brown               v(1)   = one
655d8606c27SBarry Smith               PetscCallA(MatSetValuesLocal(jac_prec,ione,row,ione,col,v,INSERT_VALUES,ierr))
656c4762a1bSJed Brown!           interior grid points
657c4762a1bSJed Brown            else
658c4762a1bSJed Brown               v(1) = -hxdhy
659c4762a1bSJed Brown               v(2) = -hydhx
660c4762a1bSJed Brown               v(3) = two*(hydhx + hxdhy) - sc*user%lambda*exp(x(i,j))
661c4762a1bSJed Brown               v(4) = -hydhx
662c4762a1bSJed Brown               v(5) = -hxdhy
663c4762a1bSJed Brown               col(1) = row - user%gxm
664c4762a1bSJed Brown               col(2) = row - 1
665c4762a1bSJed Brown               col(3) = row
666c4762a1bSJed Brown               col(4) = row + 1
667c4762a1bSJed Brown               col(5) = row + user%gxm
668d8606c27SBarry Smith               PetscCallA(MatSetValuesLocal(jac_prec,ione,row,ifive,col,v,INSERT_VALUES,ierr))
669c4762a1bSJed Brown            endif
670c4762a1bSJed Brown 10      continue
671c4762a1bSJed Brown 20   continue
672c4762a1bSJed Brown
673c4762a1bSJed Brown      return
674c4762a1bSJed Brown      end
675c4762a1bSJed Brown
676c4762a1bSJed Brown!
677c4762a1bSJed Brown!/*TEST
678c4762a1bSJed Brown!
679c4762a1bSJed Brown!   test:
680c4762a1bSJed Brown!      nsize: 4
681c4762a1bSJed 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
682c4762a1bSJed Brown!      requires: !single
683c4762a1bSJed Brown!
684c4762a1bSJed Brown!   test:
685c4762a1bSJed Brown!      suffix: 2
686c4762a1bSJed Brown!      nsize: 4
687c4762a1bSJed Brown!      args: -da_processors_x 2 -da_processors_y 2 -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
688c4762a1bSJed Brown!      requires: !single
689c4762a1bSJed Brown!
690c4762a1bSJed Brown!   test:
691c4762a1bSJed Brown!      suffix: 3
692c4762a1bSJed Brown!      nsize: 3
693c4762a1bSJed Brown!      args: -snes_fd -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
694c4762a1bSJed Brown!      requires: !single
695c4762a1bSJed Brown!
696c4762a1bSJed Brown!   test:
697c4762a1bSJed Brown!      suffix: 4
698c4762a1bSJed Brown!      nsize: 3
699c4762a1bSJed Brown!      args: -snes_mf_operator -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
700c4762a1bSJed Brown!      requires: !single
701c4762a1bSJed Brown!
702c4762a1bSJed Brown!   test:
703c4762a1bSJed Brown!      suffix: 5
704c4762a1bSJed Brown!      requires: !single
705c4762a1bSJed Brown!
706c4762a1bSJed Brown!TEST*/
707