xref: /petsc/src/snes/tutorials/ex5f.F90 (revision 65afa91a19bc4c6a7a312154dcc5d3d6fd9186ca)
1c4762a1bSJed Brown!
2c4762a1bSJed Brown!  Description: This example 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 <param>, where <param> 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!
14c4762a1bSJed Brown!  Solid Fuel Ignition (SFI) problem.  This problem is modeled by
15c4762a1bSJed Brown!  the partial differential equation
16c4762a1bSJed Brown!
17c4762a1bSJed Brown!          -Laplacian u - lambda*exp(u) = 0,  0 < x,y < 1,
18c4762a1bSJed Brown!
19c4762a1bSJed Brown!  with boundary conditions
20c4762a1bSJed Brown!
21c4762a1bSJed Brown!           u = 0  for  x = 0, x = 1, y = 0, y = 1.
22c4762a1bSJed Brown!
23c4762a1bSJed Brown!  A finite difference approximation with the usual 5-point stencil
24c4762a1bSJed Brown!  is used to discretize the boundary value problem to obtain a nonlinear
25c4762a1bSJed Brown!  system of equations.
26c4762a1bSJed Brown!
27c4762a1bSJed Brown!  --------------------------------------------------------------------------
28dfbbaf82SBarry Smith      module ex5fmodule
29dfbbaf82SBarry Smith      use petscsnes
30dfbbaf82SBarry Smith      use petscdmda
31dfbbaf82SBarry Smith#include <petsc/finclude/petscsnes.h>
32dfbbaf82SBarry Smith#include <petsc/finclude/petscdm.h>
33dfbbaf82SBarry Smith#include <petsc/finclude/petscdmda.h>
34dfbbaf82SBarry Smith      PetscInt xs,xe,xm,gxs,gxe,gxm
35dfbbaf82SBarry Smith      PetscInt ys,ye,ym,gys,gye,gym
36dfbbaf82SBarry Smith      PetscInt mx,my
37dfbbaf82SBarry Smith      PetscMPIInt rank,size
38dfbbaf82SBarry Smith      PetscReal lambda
39dfbbaf82SBarry Smith      end module ex5fmodule
40c4762a1bSJed Brown
41c4762a1bSJed Brown      program main
42dfbbaf82SBarry Smith      use ex5fmodule
43c4762a1bSJed Brown      implicit none
44c4762a1bSJed Brown
45c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
46c4762a1bSJed Brown!                   Variable declarations
47c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48c4762a1bSJed Brown!
49c4762a1bSJed Brown!  Variables:
50c4762a1bSJed Brown!     snes        - nonlinear solver
51c4762a1bSJed Brown!     x, r        - solution, residual vectors
52c4762a1bSJed Brown!     its         - iterations for convergence
53c4762a1bSJed Brown!
54c4762a1bSJed Brown!  See additional variable declarations in the file ex5f.h
55c4762a1bSJed Brown!
56c4762a1bSJed Brown      SNES           snes
57c4762a1bSJed Brown      Vec            x,r
58c4762a1bSJed Brown      PetscInt       its,i1,i4
59c4762a1bSJed Brown      PetscErrorCode ierr
60c4762a1bSJed Brown      PetscReal      lambda_max,lambda_min
61c4762a1bSJed Brown      PetscBool      flg
62c4762a1bSJed Brown      DM             da
63c4762a1bSJed Brown
64c4762a1bSJed Brown!  Note: Any user-defined Fortran routines (such as FormJacobianLocal)
65c4762a1bSJed Brown!  MUST be declared as external.
66c4762a1bSJed Brown
67c4762a1bSJed Brown      external FormInitialGuess
68c4762a1bSJed Brown      external FormFunctionLocal,FormJacobianLocal
69c4762a1bSJed Brown      external MySNESConverged
70c4762a1bSJed Brown
71c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
72c4762a1bSJed Brown!  Initialize program
73c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
74c4762a1bSJed Brown
75*65afa91aSSatish Balay      call PetscInitialize(ierr)
76*65afa91aSSatish Balay      CHKERRA(ierr)
77*65afa91aSSatish Balay      call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
78*65afa91aSSatish Balay      CHKERRMPIA(ierr)
79*65afa91aSSatish Balay      call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
80*65afa91aSSatish Balay      CHKERRMPIA(ierr)
81c4762a1bSJed Brown!  Initialize problem parameters
82c4762a1bSJed Brown
83c4762a1bSJed Brown      i1 = 1
84c4762a1bSJed Brown      i4 = 4
85c4762a1bSJed Brown      lambda_max = 6.81
86c4762a1bSJed Brown      lambda_min = 0.0
87c4762a1bSJed Brown      lambda     = 6.0
88*65afa91aSSatish Balay      call PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-par',lambda,PETSC_NULL_BOOL,ierr)
89*65afa91aSSatish Balay      CHKERRA(ierr)
90*65afa91aSSatish Balay
91c4762a1bSJed Brown! this statement is split into multiple-lines to keep lines under 132 char limit - required by 'make check'
92c4762a1bSJed Brown      if (lambda .ge. lambda_max .or. lambda .le. lambda_min) then
93c4762a1bSJed Brown        ierr = PETSC_ERR_ARG_OUTOFRANGE; SETERRA(PETSC_COMM_WORLD,ierr,'Lambda')
94c4762a1bSJed Brown      endif
95c4762a1bSJed Brown
96c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97c4762a1bSJed Brown!  Create nonlinear solver context
98c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
99c4762a1bSJed Brown
100*65afa91aSSatish Balay      call SNESCreate(PETSC_COMM_WORLD,snes,ierr)
101*65afa91aSSatish Balay      CHKERRA(ierr)
102c4762a1bSJed Brown
103c4762a1bSJed Brown!  Set convergence test routine if desired
104c4762a1bSJed Brown
105*65afa91aSSatish Balay      call PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-my_snes_convergence',flg,ierr)
106*65afa91aSSatish Balay      CHKERRA(ierr)
107c4762a1bSJed Brown      if (flg) then
108*65afa91aSSatish Balay        call SNESSetConvergenceTest(snes,MySNESConverged,0,PETSC_NULL_FUNCTION,ierr)
109*65afa91aSSatish Balay        CHKERRA(ierr)
110c4762a1bSJed Brown      endif
111c4762a1bSJed Brown
112c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113c4762a1bSJed Brown!  Create vector data structures; set function evaluation routine
114c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115c4762a1bSJed Brown
116c4762a1bSJed Brown!  Create distributed array (DMDA) to manage parallel grid and vectors
117c4762a1bSJed Brown
11860cf0239SBarry Smith!     This really needs only the star-type stencil, but we use the box stencil temporarily.
11960cf0239SBarry Smith
12060cf0239SBarry Smith#if defined(PETSC_HAVE_FORTRAN_FREE_LINE_LENGTH_NONE)
121d8606c27SBarry Smith      PetscCallA(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,i4,i4,PETSC_DECIDE,PETSC_DECIDE,i1,i1,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da,ierr))
12260cf0239SBarry Smith#else
123*65afa91aSSatish Balay      call DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,i4,i4,PETSC_DECIDE,PETSC_DECIDE, &
124*65afa91aSSatish Balay                        i1,i1, PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da,ierr)
125*65afa91aSSatish Balay      CHKERRA(ierr)
12660cf0239SBarry Smith#endif
127*65afa91aSSatish Balay      call DMSetFromOptions(da,ierr)
128*65afa91aSSatish Balay      CHKERRA(ierr)
129*65afa91aSSatish Balay      call DMSetUp(da,ierr)
130*65afa91aSSatish Balay      CHKERRA(ierr)
131c4762a1bSJed Brown
132c4762a1bSJed Brown!  Extract global and local vectors from DMDA; then duplicate for remaining
133c4762a1bSJed Brown!  vectors that are the same types
134c4762a1bSJed Brown
135*65afa91aSSatish Balay      call DMCreateGlobalVector(da,x,ierr)
136*65afa91aSSatish Balay      CHKERRA(ierr)
137*65afa91aSSatish Balay      call VecDuplicate(x,r,ierr)
138*65afa91aSSatish Balay      CHKERRA(ierr)
139c4762a1bSJed Brown
140c4762a1bSJed Brown!  Get local grid boundaries (for 2-dimensional DMDA)
141c4762a1bSJed Brown
14260cf0239SBarry Smith#if defined(PETSC_HAVE_FORTRAN_FREE_LINE_LENGTH_NONE)
143d8606c27SBarry Smith      PetscCallA(DMDAGetInfo(da,PETSC_NULL_INTEGER,mx,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))
14460cf0239SBarry Smith#else
14560cf0239SBarry Smith      call DMDAGetInfo(da,PETSC_NULL_INTEGER,mx,my,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
14660cf0239SBarry Smith                       PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
14760cf0239SBarry Smith                       PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ierr)
148*65afa91aSSatish Balay      CHKERRA(ierr)
14960cf0239SBarry Smith#endif
150*65afa91aSSatish Balay      call DMDAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym,PETSC_NULL_INTEGER,ierr)
151*65afa91aSSatish Balay      CHKERRA(ierr)
152*65afa91aSSatish Balay      call DMDAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,gxm,gym,PETSC_NULL_INTEGER,ierr)
153*65afa91aSSatish Balay      CHKERRA(ierr)
154c4762a1bSJed Brown
155c4762a1bSJed Brown!  Here we shift the starting indices up by one so that we can easily
156c4762a1bSJed Brown!  use the Fortran convention of 1-based indices (rather 0-based indices).
157c4762a1bSJed Brown
158c4762a1bSJed Brown      xs  = xs+1
159c4762a1bSJed Brown      ys  = ys+1
160c4762a1bSJed Brown      gxs = gxs+1
161c4762a1bSJed Brown      gys = gys+1
162c4762a1bSJed Brown
163c4762a1bSJed Brown      ye  = ys+ym-1
164c4762a1bSJed Brown      xe  = xs+xm-1
165c4762a1bSJed Brown      gye = gys+gym-1
166c4762a1bSJed Brown      gxe = gxs+gxm-1
167c4762a1bSJed Brown
168c4762a1bSJed Brown!  Set function evaluation routine and vector
169c4762a1bSJed Brown
170*65afa91aSSatish Balay      call DMDASNESSetFunctionLocal(da,INSERT_VALUES,FormFunctionLocal,da,ierr)
171*65afa91aSSatish Balay      CHKERRA(ierr)
172*65afa91aSSatish Balay      call DMDASNESSetJacobianLocal(da,FormJacobianLocal,da,ierr)
173*65afa91aSSatish Balay      CHKERRA(ierr)
174*65afa91aSSatish Balay      call SNESSetDM(snes,da,ierr)
175*65afa91aSSatish Balay      CHKERRA(ierr)
176c4762a1bSJed Brown
177c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178c4762a1bSJed Brown!  Customize nonlinear solver; set runtime options
179c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
180c4762a1bSJed Brown
181c4762a1bSJed Brown!  Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
182c4762a1bSJed Brown
183*65afa91aSSatish Balay      call SNESSetFromOptions(snes,ierr)
184*65afa91aSSatish Balay      CHKERRA(ierr)
185c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
186c4762a1bSJed Brown!  Evaluate initial guess; then solve nonlinear system.
187c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
188c4762a1bSJed Brown
189c4762a1bSJed Brown!  Note: The user should initialize the vector, x, with the initial guess
190c4762a1bSJed Brown!  for the nonlinear solver prior to calling SNESSolve().  In particular,
191c4762a1bSJed Brown!  to employ an initial guess of zero, the user should explicitly set
192c4762a1bSJed Brown!  this vector to zero by calling VecSet().
193c4762a1bSJed Brown
194*65afa91aSSatish Balay      call FormInitialGuess(x,ierr)
195*65afa91aSSatish Balay      CHKERRA(ierr)
196*65afa91aSSatish Balay      call SNESSolve(snes,PETSC_NULL_VEC,x,ierr)
197*65afa91aSSatish Balay      CHKERRA(ierr)
198*65afa91aSSatish Balay      call SNESGetIterationNumber(snes,its,ierr)
199*65afa91aSSatish Balay      CHKERRA(ierr)
200c4762a1bSJed Brown      if (rank .eq. 0) then
201c4762a1bSJed Brown         write(6,100) its
202c4762a1bSJed Brown      endif
203c4762a1bSJed Brown  100 format('Number of SNES iterations = ',i5)
204c4762a1bSJed Brown
205c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
206c4762a1bSJed Brown!  Free work space.  All PETSc objects should be destroyed when they
207c4762a1bSJed Brown!  are no longer needed.
208c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
209c4762a1bSJed Brown
210*65afa91aSSatish Balay      call VecDestroy(x,ierr)
211*65afa91aSSatish Balay      CHKERRA(ierr)
212*65afa91aSSatish Balay      call VecDestroy(r,ierr)
213*65afa91aSSatish Balay      CHKERRA(ierr)
214*65afa91aSSatish Balay      call SNESDestroy(snes,ierr)
215*65afa91aSSatish Balay      CHKERRA(ierr)
216*65afa91aSSatish Balay      call DMDestroy(da,ierr)
217*65afa91aSSatish Balay      CHKERRA(ierr)
218*65afa91aSSatish Balay      call PetscFinalize(ierr)
219*65afa91aSSatish Balay      CHKERRA(ierr)
220c4762a1bSJed Brown      end
221c4762a1bSJed Brown
222c4762a1bSJed Brown! ---------------------------------------------------------------------
223c4762a1bSJed Brown!
224c4762a1bSJed Brown!  FormInitialGuess - Forms initial approximation.
225c4762a1bSJed Brown!
226c4762a1bSJed Brown!  Input Parameters:
227c4762a1bSJed Brown!  X - vector
228c4762a1bSJed Brown!
229c4762a1bSJed Brown!  Output Parameter:
230c4762a1bSJed Brown!  X - vector
231c4762a1bSJed Brown!
232c4762a1bSJed Brown!  Notes:
233c4762a1bSJed Brown!  This routine serves as a wrapper for the lower-level routine
234c4762a1bSJed Brown!  "ApplicationInitialGuess", where the actual computations are
235c4762a1bSJed Brown!  done using the standard Fortran style of treating the local
236c4762a1bSJed Brown!  vector data as a multidimensional array over the local mesh.
237c4762a1bSJed Brown!  This routine merely handles ghost point scatters and accesses
238c4762a1bSJed Brown!  the local vector data via VecGetArray() and VecRestoreArray().
239c4762a1bSJed Brown!
240c4762a1bSJed Brown      subroutine FormInitialGuess(X,ierr)
241dfbbaf82SBarry Smith      use ex5fmodule
242c4762a1bSJed Brown      implicit none
243c4762a1bSJed Brown
244c4762a1bSJed Brown!  Input/output variables:
245c4762a1bSJed Brown      Vec      X
246c4762a1bSJed Brown      PetscErrorCode  ierr
247c4762a1bSJed Brown!  Declarations for use with local arrays:
248c4762a1bSJed Brown      PetscScalar lx_v(0:1)
249c4762a1bSJed Brown      PetscOffset lx_i
250c4762a1bSJed Brown
251c4762a1bSJed Brown      ierr = 0
252c4762a1bSJed Brown
253c4762a1bSJed Brown!  Get a pointer to vector data.
254c4762a1bSJed Brown!    - For default PETSc vectors, VecGetArray() returns a pointer to
255c4762a1bSJed Brown!      the data array.  Otherwise, the routine is implementation dependent.
256c4762a1bSJed Brown!    - You MUST call VecRestoreArray() when you no longer need access to
257c4762a1bSJed Brown!      the array.
258c4762a1bSJed Brown!    - Note that the Fortran interface to VecGetArray() differs from the
259c4762a1bSJed Brown!      C version.  See the users manual for details.
260c4762a1bSJed Brown
261*65afa91aSSatish Balay      call VecGetArray(X,lx_v,lx_i,ierr)
262*65afa91aSSatish Balay      CHKERRQ(ierr)
263c4762a1bSJed Brown
264c4762a1bSJed Brown!  Compute initial guess over the locally owned part of the grid
265c4762a1bSJed Brown
266*65afa91aSSatish Balay      call InitialGuessLocal(lx_v(lx_i),ierr)
267*65afa91aSSatish Balay      CHKERRQ(ierr)
268c4762a1bSJed Brown
269c4762a1bSJed Brown!  Restore vector
270c4762a1bSJed Brown
271*65afa91aSSatish Balay      call VecRestoreArray(X,lx_v,lx_i,ierr)
272*65afa91aSSatish Balay      CHKERRQ(ierr)
273c4762a1bSJed Brown
274c4762a1bSJed Brown      return
275c4762a1bSJed Brown      end
276c4762a1bSJed Brown
277c4762a1bSJed Brown! ---------------------------------------------------------------------
278c4762a1bSJed Brown!
279c4762a1bSJed Brown!  InitialGuessLocal - Computes initial approximation, called by
280c4762a1bSJed Brown!  the higher level routine FormInitialGuess().
281c4762a1bSJed Brown!
282c4762a1bSJed Brown!  Input Parameter:
283c4762a1bSJed Brown!  x - local vector data
284c4762a1bSJed Brown!
285c4762a1bSJed Brown!  Output Parameters:
286c4762a1bSJed Brown!  x - local vector data
287c4762a1bSJed Brown!  ierr - error code
288c4762a1bSJed Brown!
289c4762a1bSJed Brown!  Notes:
290c4762a1bSJed Brown!  This routine uses standard Fortran-style computations over a 2-dim array.
291c4762a1bSJed Brown!
292c4762a1bSJed Brown      subroutine InitialGuessLocal(x,ierr)
293dfbbaf82SBarry Smith      use ex5fmodule
294c4762a1bSJed Brown      implicit none
295c4762a1bSJed Brown
296c4762a1bSJed Brown!  Input/output variables:
297c4762a1bSJed Brown      PetscScalar    x(xs:xe,ys:ye)
298c4762a1bSJed Brown      PetscErrorCode ierr
299c4762a1bSJed Brown
300c4762a1bSJed Brown!  Local variables:
301c4762a1bSJed Brown      PetscInt  i,j
302c4762a1bSJed Brown      PetscReal temp1,temp,one,hx,hy
303c4762a1bSJed Brown
304c4762a1bSJed Brown!  Set parameters
305c4762a1bSJed Brown
306c4762a1bSJed Brown      ierr   = 0
307c4762a1bSJed Brown      one    = 1.0
308c4762a1bSJed Brown      hx     = one/((mx-1))
309c4762a1bSJed Brown      hy     = one/((my-1))
310c4762a1bSJed Brown      temp1  = lambda/(lambda + one)
311c4762a1bSJed Brown
312c4762a1bSJed Brown      do 20 j=ys,ye
313c4762a1bSJed Brown         temp = (min(j-1,my-j))*hy
314c4762a1bSJed Brown         do 10 i=xs,xe
315c4762a1bSJed Brown            if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then
316c4762a1bSJed Brown              x(i,j) = 0.0
317c4762a1bSJed Brown            else
318c4762a1bSJed Brown              x(i,j) = temp1 * sqrt(min(min(i-1,mx-i)*hx,(temp)))
319c4762a1bSJed Brown            endif
320c4762a1bSJed Brown 10      continue
321c4762a1bSJed Brown 20   continue
322c4762a1bSJed Brown
323c4762a1bSJed Brown      return
324c4762a1bSJed Brown      end
325c4762a1bSJed Brown
326c4762a1bSJed Brown! ---------------------------------------------------------------------
327c4762a1bSJed Brown!
328c4762a1bSJed Brown!  FormFunctionLocal - Computes nonlinear function, called by
329c4762a1bSJed Brown!  the higher level routine FormFunction().
330c4762a1bSJed Brown!
331c4762a1bSJed Brown!  Input Parameter:
332c4762a1bSJed Brown!  x - local vector data
333c4762a1bSJed Brown!
334c4762a1bSJed Brown!  Output Parameters:
335c4762a1bSJed Brown!  f - local vector data, f(x)
336c4762a1bSJed Brown!  ierr - error code
337c4762a1bSJed Brown!
338c4762a1bSJed Brown!  Notes:
339c4762a1bSJed Brown!  This routine uses standard Fortran-style computations over a 2-dim array.
340c4762a1bSJed Brown!
341c4762a1bSJed Brown!
342c4762a1bSJed Brown      subroutine FormFunctionLocal(info,x,f,da,ierr)
343dfbbaf82SBarry Smith      use ex5fmodule
344c4762a1bSJed Brown      implicit none
345c4762a1bSJed Brown
346c4762a1bSJed Brown      DM da
347c4762a1bSJed Brown
348c4762a1bSJed Brown!  Input/output variables:
349c4762a1bSJed Brown      DMDALocalInfo info(DMDA_LOCAL_INFO_SIZE)
350c4762a1bSJed Brown      PetscScalar x(gxs:gxe,gys:gye)
351c4762a1bSJed Brown      PetscScalar f(xs:xe,ys:ye)
352c4762a1bSJed Brown      PetscErrorCode     ierr
353c4762a1bSJed Brown
354c4762a1bSJed Brown!  Local variables:
355c4762a1bSJed Brown      PetscScalar two,one,hx,hy
356c4762a1bSJed Brown      PetscScalar hxdhy,hydhx,sc
357c4762a1bSJed Brown      PetscScalar u,uxx,uyy
358c4762a1bSJed Brown      PetscInt  i,j
359c4762a1bSJed Brown
360c4762a1bSJed Brown      xs     = info(DMDA_LOCAL_INFO_XS)+1
361c4762a1bSJed Brown      xe     = xs+info(DMDA_LOCAL_INFO_XM)-1
362c4762a1bSJed Brown      ys     = info(DMDA_LOCAL_INFO_YS)+1
363c4762a1bSJed Brown      ye     = ys+info(DMDA_LOCAL_INFO_YM)-1
364c4762a1bSJed Brown      mx     = info(DMDA_LOCAL_INFO_MX)
365c4762a1bSJed Brown      my     = info(DMDA_LOCAL_INFO_MY)
366c4762a1bSJed Brown
367c4762a1bSJed Brown      one    = 1.0
368c4762a1bSJed Brown      two    = 2.0
369c4762a1bSJed Brown      hx     = one/(mx-1)
370c4762a1bSJed Brown      hy     = one/(my-1)
371c4762a1bSJed Brown      sc     = hx*hy*lambda
372c4762a1bSJed Brown      hxdhy  = hx/hy
373c4762a1bSJed Brown      hydhx  = hy/hx
374c4762a1bSJed Brown
375c4762a1bSJed Brown!  Compute function over the locally owned part of the grid
376c4762a1bSJed Brown
377c4762a1bSJed Brown      do 20 j=ys,ye
378c4762a1bSJed Brown         do 10 i=xs,xe
379c4762a1bSJed Brown            if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then
380c4762a1bSJed Brown               f(i,j) = x(i,j)
381c4762a1bSJed Brown            else
382c4762a1bSJed Brown               u = x(i,j)
383c4762a1bSJed Brown               uxx = hydhx * (two*u - x(i-1,j) - x(i+1,j))
384c4762a1bSJed Brown               uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1))
385c4762a1bSJed Brown               f(i,j) = uxx + uyy - sc*exp(u)
386c4762a1bSJed Brown            endif
387c4762a1bSJed Brown 10      continue
388c4762a1bSJed Brown 20   continue
389c4762a1bSJed Brown
390*65afa91aSSatish Balay      call PetscLogFlops(11.0d0*ym*xm,ierr)
391*65afa91aSSatish Balay      CHKERRQ(ierr)
392c4762a1bSJed Brown
393c4762a1bSJed Brown      return
394c4762a1bSJed Brown      end
395c4762a1bSJed Brown
396c4762a1bSJed Brown! ---------------------------------------------------------------------
397c4762a1bSJed Brown!
398c4762a1bSJed Brown!  FormJacobianLocal - Computes Jacobian matrix, called by
399c4762a1bSJed Brown!  the higher level routine FormJacobian().
400c4762a1bSJed Brown!
401c4762a1bSJed Brown!  Input Parameters:
402c4762a1bSJed Brown!  x        - local vector data
403c4762a1bSJed Brown!
404c4762a1bSJed Brown!  Output Parameters:
405c4762a1bSJed Brown!  jac      - Jacobian matrix
406c4762a1bSJed Brown!  jac_prec - optionally different preconditioning matrix (not used here)
407c4762a1bSJed Brown!  ierr     - error code
408c4762a1bSJed Brown!
409c4762a1bSJed Brown!  Notes:
410c4762a1bSJed Brown!  This routine uses standard Fortran-style computations over a 2-dim array.
411c4762a1bSJed Brown!
412c4762a1bSJed Brown!  Notes:
413c4762a1bSJed Brown!  Due to grid point reordering with DMDAs, we must always work
414c4762a1bSJed Brown!  with the local grid points, and then transform them to the new
415c4762a1bSJed Brown!  global numbering with the "ltog" mapping
416c4762a1bSJed Brown!  We cannot work directly with the global numbers for the original
417c4762a1bSJed Brown!  uniprocessor grid!
418c4762a1bSJed Brown!
419c4762a1bSJed Brown!  Two methods are available for imposing this transformation
420c4762a1bSJed Brown!  when setting matrix entries:
421c4762a1bSJed Brown!    (A) MatSetValuesLocal(), using the local ordering (including
422c4762a1bSJed Brown!        ghost points!)
423c4762a1bSJed Brown!          by calling MatSetValuesLocal()
424c4762a1bSJed Brown!    (B) MatSetValues(), using the global ordering
425c4762a1bSJed Brown!        - Use DMDAGetGlobalIndices() to extract the local-to-global map
426c4762a1bSJed Brown!        - Then apply this map explicitly yourself
427c4762a1bSJed Brown!        - Set matrix entries using the global ordering by calling
428c4762a1bSJed Brown!          MatSetValues()
429c4762a1bSJed Brown!  Option (A) seems cleaner/easier in many cases, and is the procedure
430c4762a1bSJed Brown!  used in this example.
431c4762a1bSJed Brown!
432c4762a1bSJed Brown      subroutine FormJacobianLocal(info,x,A,jac,da,ierr)
433dfbbaf82SBarry Smith      use ex5fmodule
434c4762a1bSJed Brown      implicit none
435c4762a1bSJed Brown
436c4762a1bSJed Brown      DM da
437c4762a1bSJed Brown
438c4762a1bSJed Brown!  Input/output variables:
439c4762a1bSJed Brown      PetscScalar x(gxs:gxe,gys:gye)
440c4762a1bSJed Brown      Mat         A,jac
441c4762a1bSJed Brown      PetscErrorCode  ierr
442c4762a1bSJed Brown      DMDALocalInfo info(DMDA_LOCAL_INFO_SIZE)
443c4762a1bSJed Brown
444c4762a1bSJed Brown!  Local variables:
445c4762a1bSJed Brown      PetscInt  row,col(5),i,j,i1,i5
446c4762a1bSJed Brown      PetscScalar two,one,hx,hy,v(5)
447c4762a1bSJed Brown      PetscScalar hxdhy,hydhx,sc
448c4762a1bSJed Brown
449c4762a1bSJed Brown!  Set parameters
450c4762a1bSJed Brown
451c4762a1bSJed Brown      i1     = 1
452c4762a1bSJed Brown      i5     = 5
453c4762a1bSJed Brown      one    = 1.0
454c4762a1bSJed Brown      two    = 2.0
455c4762a1bSJed Brown      hx     = one/(mx-1)
456c4762a1bSJed Brown      hy     = one/(my-1)
457c4762a1bSJed Brown      sc     = hx*hy
458c4762a1bSJed Brown      hxdhy  = hx/hy
459c4762a1bSJed Brown      hydhx  = hy/hx
460c4762a1bSJed Brown
461c4762a1bSJed Brown!  Compute entries for the locally owned part of the Jacobian.
462c4762a1bSJed Brown!   - Currently, all PETSc parallel matrix formats are partitioned by
463c4762a1bSJed Brown!     contiguous chunks of rows across the processors.
464c4762a1bSJed Brown!   - Each processor needs to insert only elements that it owns
465c4762a1bSJed Brown!     locally (but any non-local elements will be sent to the
466c4762a1bSJed Brown!     appropriate processor during matrix assembly).
467c4762a1bSJed Brown!   - Here, we set all entries for a particular row at once.
468c4762a1bSJed Brown!   - We can set matrix entries either using either
469c4762a1bSJed Brown!     MatSetValuesLocal() or MatSetValues(), as discussed above.
470c4762a1bSJed Brown!   - Note that MatSetValues() uses 0-based row and column numbers
471c4762a1bSJed Brown!     in Fortran as well as in C.
472c4762a1bSJed Brown
473c4762a1bSJed Brown      do 20 j=ys,ye
474c4762a1bSJed Brown         row = (j - gys)*gxm + xs - gxs - 1
475c4762a1bSJed Brown         do 10 i=xs,xe
476c4762a1bSJed Brown            row = row + 1
477c4762a1bSJed Brown!           boundary points
478c4762a1bSJed Brown            if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then
479c4762a1bSJed Brown!       Some f90 compilers need 4th arg to be of same type in both calls
480c4762a1bSJed Brown               col(1) = row
481c4762a1bSJed Brown               v(1)   = one
482*65afa91aSSatish Balay               call MatSetValuesLocal(jac,i1,row,i1,col,v,INSERT_VALUES,ierr)
483*65afa91aSSatish Balay               CHKERRQ(ierr)
484c4762a1bSJed Brown!           interior grid points
485c4762a1bSJed Brown            else
486c4762a1bSJed Brown               v(1) = -hxdhy
487c4762a1bSJed Brown               v(2) = -hydhx
488c4762a1bSJed Brown               v(3) = two*(hydhx + hxdhy) - sc*lambda*exp(x(i,j))
489c4762a1bSJed Brown               v(4) = -hydhx
490c4762a1bSJed Brown               v(5) = -hxdhy
491c4762a1bSJed Brown               col(1) = row - gxm
492c4762a1bSJed Brown               col(2) = row - 1
493c4762a1bSJed Brown               col(3) = row
494c4762a1bSJed Brown               col(4) = row + 1
495c4762a1bSJed Brown               col(5) = row + gxm
496*65afa91aSSatish Balay               call MatSetValuesLocal(jac,i1,row,i5,col,v, INSERT_VALUES,ierr)
497*65afa91aSSatish Balay               CHKERRQ(ierr)
498c4762a1bSJed Brown            endif
499c4762a1bSJed Brown 10      continue
500c4762a1bSJed Brown 20   continue
501*65afa91aSSatish Balay      call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
502*65afa91aSSatish Balay      CHKERRQ(ierr)
503*65afa91aSSatish Balay      call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
504*65afa91aSSatish Balay      CHKERRQ(ierr)
505c4762a1bSJed Brown      if (A .ne. jac) then
506*65afa91aSSatish Balay         call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
507*65afa91aSSatish Balay         CHKERRQ(ierr)
508*65afa91aSSatish Balay         call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)
509*65afa91aSSatish Balay         CHKERRQ(ierr)
510c4762a1bSJed Brown      endif
511c4762a1bSJed Brown      return
512c4762a1bSJed Brown      end
513c4762a1bSJed Brown
514c4762a1bSJed Brown!
515c4762a1bSJed Brown!     Simple convergence test based on the infinity norm of the residual being small
516c4762a1bSJed Brown!
517c4762a1bSJed Brown      subroutine MySNESConverged(snes,it,xnorm,snorm,fnorm,reason,dummy,ierr)
518dfbbaf82SBarry Smith      use ex5fmodule
519c4762a1bSJed Brown      implicit none
520c4762a1bSJed Brown
521c4762a1bSJed Brown      SNES snes
522c4762a1bSJed Brown      PetscInt it,dummy
523c4762a1bSJed Brown      PetscReal xnorm,snorm,fnorm,nrm
524c4762a1bSJed Brown      SNESConvergedReason reason
525c4762a1bSJed Brown      Vec f
526c4762a1bSJed Brown      PetscErrorCode ierr
527c4762a1bSJed Brown
528*65afa91aSSatish Balay      call SNESGetFunction(snes,f,PETSC_NULL_FUNCTION,dummy,ierr)
529*65afa91aSSatish Balay      CHKERRQ(ierr)
530*65afa91aSSatish Balay      call VecNorm(f,NORM_INFINITY,nrm,ierr)
531*65afa91aSSatish Balay      CHKERRQ(ierr)
532c4762a1bSJed Brown      if (nrm .le. 1.e-5) reason = SNES_CONVERGED_FNORM_ABS
533c4762a1bSJed Brown
534c4762a1bSJed Brown      end
535c4762a1bSJed Brown
536c4762a1bSJed Brown!/*TEST
537c4762a1bSJed Brown!
538c4762a1bSJed Brown!   build:
539c4762a1bSJed Brown!      requires: !complex !single
540c4762a1bSJed Brown!
541c4762a1bSJed Brown!   test:
542c4762a1bSJed Brown!      nsize: 4
5438f8b3c79SStefano Zampini!      args: -snes_mf -pc_type none -da_processors_x 4 -da_processors_y 1 -snes_monitor_short \
5448f8b3c79SStefano Zampini!            -ksp_gmres_cgs_refinement_type refine_always
545c4762a1bSJed Brown!
546c4762a1bSJed Brown!   test:
547c4762a1bSJed Brown!      suffix: 2
548c4762a1bSJed Brown!      nsize: 4
549c4762a1bSJed Brown!      args: -da_processors_x 2 -da_processors_y 2 -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
550c4762a1bSJed Brown!
551c4762a1bSJed Brown!   test:
552c4762a1bSJed Brown!      suffix: 3
553c4762a1bSJed Brown!      nsize: 3
554c4762a1bSJed Brown!      args: -snes_fd -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
555c4762a1bSJed Brown!
556c4762a1bSJed Brown!   test:
557c4762a1bSJed Brown!      suffix: 6
558c4762a1bSJed Brown!      nsize: 1
559c4762a1bSJed Brown!      args: -snes_monitor_short -my_snes_convergence
560c4762a1bSJed Brown!
561c4762a1bSJed Brown!TEST*/
562