xref: /petsc/src/snes/tutorials/ex5f.F90 (revision 42ce371b2bd7d45eb85bb2bb31075ac1967f9fc8)
1c4762a1bSJed Brown!
2*42ce371bSBarry Smith!  This example shows how to avoid Fortran line lengths larger than 132 characters.
3*42ce371bSBarry Smith!  We recommend starting from src/snes/tutorials/ex5f90.F90 instead of this example
4*42ce371bSBarry Smith!
5c4762a1bSJed Brown!  Description: This example solves a nonlinear system in parallel with SNES.
6c4762a1bSJed Brown!  We solve the  Bratu (SFI - solid fuel ignition) problem in a 2D rectangular
7c4762a1bSJed Brown!  domain, using distributed arrays (DMDAs) to partition the parallel grid.
8c4762a1bSJed Brown!  The command line options include:
9c4762a1bSJed Brown!    -par <param>, where <param> indicates the nonlinearity of the problem
10c4762a1bSJed Brown!       problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)
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
7565afa91aSSatish Balay      call PetscInitialize(ierr)
7665afa91aSSatish Balay      CHKERRA(ierr)
7765afa91aSSatish Balay      call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
7865afa91aSSatish Balay      CHKERRMPIA(ierr)
7965afa91aSSatish Balay      call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
8065afa91aSSatish 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
8865afa91aSSatish Balay      call PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-par',lambda,PETSC_NULL_BOOL,ierr)
8965afa91aSSatish Balay      CHKERRA(ierr)
9065afa91aSSatish 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
10065afa91aSSatish Balay      call SNESCreate(PETSC_COMM_WORLD,snes,ierr)
10165afa91aSSatish Balay      CHKERRA(ierr)
102c4762a1bSJed Brown
103c4762a1bSJed Brown!  Set convergence test routine if desired
104c4762a1bSJed Brown
10565afa91aSSatish Balay      call PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-my_snes_convergence',flg,ierr)
10665afa91aSSatish Balay      CHKERRA(ierr)
107c4762a1bSJed Brown      if (flg) then
10865afa91aSSatish Balay        call SNESSetConvergenceTest(snes,MySNESConverged,0,PETSC_NULL_FUNCTION,ierr)
10965afa91aSSatish 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
118*42ce371bSBarry Smith!     This really needs only the star-type stencil, but we use the box stencil
11960cf0239SBarry Smith
12065afa91aSSatish Balay      call DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,i4,i4,PETSC_DECIDE,PETSC_DECIDE, &
12165afa91aSSatish Balay                        i1,i1, PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da,ierr)
12265afa91aSSatish Balay      CHKERRA(ierr)
12365afa91aSSatish Balay      call DMSetFromOptions(da,ierr)
12465afa91aSSatish Balay      CHKERRA(ierr)
12565afa91aSSatish Balay      call DMSetUp(da,ierr)
12665afa91aSSatish Balay      CHKERRA(ierr)
127c4762a1bSJed Brown
128c4762a1bSJed Brown!  Extract global and local vectors from DMDA; then duplicate for remaining
129c4762a1bSJed Brown!  vectors that are the same types
130c4762a1bSJed Brown
13165afa91aSSatish Balay      call DMCreateGlobalVector(da,x,ierr)
13265afa91aSSatish Balay      CHKERRA(ierr)
13365afa91aSSatish Balay      call VecDuplicate(x,r,ierr)
13465afa91aSSatish Balay      CHKERRA(ierr)
135c4762a1bSJed Brown
136c4762a1bSJed Brown!  Get local grid boundaries (for 2-dimensional DMDA)
137c4762a1bSJed Brown
13860cf0239SBarry Smith      call DMDAGetInfo(da,PETSC_NULL_INTEGER,mx,my,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
13960cf0239SBarry Smith                       PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
14060cf0239SBarry Smith                       PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ierr)
14165afa91aSSatish Balay      CHKERRA(ierr)
14265afa91aSSatish Balay      call DMDAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym,PETSC_NULL_INTEGER,ierr)
14365afa91aSSatish Balay      CHKERRA(ierr)
14465afa91aSSatish Balay      call DMDAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,gxm,gym,PETSC_NULL_INTEGER,ierr)
14565afa91aSSatish Balay      CHKERRA(ierr)
146c4762a1bSJed Brown
147c4762a1bSJed Brown!  Here we shift the starting indices up by one so that we can easily
148c4762a1bSJed Brown!  use the Fortran convention of 1-based indices (rather 0-based indices).
149c4762a1bSJed Brown
150c4762a1bSJed Brown      xs  = xs+1
151c4762a1bSJed Brown      ys  = ys+1
152c4762a1bSJed Brown      gxs = gxs+1
153c4762a1bSJed Brown      gys = gys+1
154c4762a1bSJed Brown
155c4762a1bSJed Brown      ye  = ys+ym-1
156c4762a1bSJed Brown      xe  = xs+xm-1
157c4762a1bSJed Brown      gye = gys+gym-1
158c4762a1bSJed Brown      gxe = gxs+gxm-1
159c4762a1bSJed Brown
160c4762a1bSJed Brown!  Set function evaluation routine and vector
161c4762a1bSJed Brown
16265afa91aSSatish Balay      call DMDASNESSetFunctionLocal(da,INSERT_VALUES,FormFunctionLocal,da,ierr)
16365afa91aSSatish Balay      CHKERRA(ierr)
16465afa91aSSatish Balay      call DMDASNESSetJacobianLocal(da,FormJacobianLocal,da,ierr)
16565afa91aSSatish Balay      CHKERRA(ierr)
16665afa91aSSatish Balay      call SNESSetDM(snes,da,ierr)
16765afa91aSSatish Balay      CHKERRA(ierr)
168c4762a1bSJed Brown
169c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
170c4762a1bSJed Brown!  Customize nonlinear solver; set runtime options
171c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
172c4762a1bSJed Brown
173c4762a1bSJed Brown!  Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
174c4762a1bSJed Brown
17565afa91aSSatish Balay      call SNESSetFromOptions(snes,ierr)
17665afa91aSSatish Balay      CHKERRA(ierr)
177c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178c4762a1bSJed Brown!  Evaluate initial guess; then solve nonlinear system.
179c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
180c4762a1bSJed Brown
181c4762a1bSJed Brown!  Note: The user should initialize the vector, x, with the initial guess
182c4762a1bSJed Brown!  for the nonlinear solver prior to calling SNESSolve().  In particular,
183c4762a1bSJed Brown!  to employ an initial guess of zero, the user should explicitly set
184c4762a1bSJed Brown!  this vector to zero by calling VecSet().
185c4762a1bSJed Brown
18665afa91aSSatish Balay      call FormInitialGuess(x,ierr)
18765afa91aSSatish Balay      CHKERRA(ierr)
18865afa91aSSatish Balay      call SNESSolve(snes,PETSC_NULL_VEC,x,ierr)
18965afa91aSSatish Balay      CHKERRA(ierr)
19065afa91aSSatish Balay      call SNESGetIterationNumber(snes,its,ierr)
19165afa91aSSatish Balay      CHKERRA(ierr)
192c4762a1bSJed Brown      if (rank .eq. 0) then
193c4762a1bSJed Brown         write(6,100) its
194c4762a1bSJed Brown      endif
195c4762a1bSJed Brown  100 format('Number of SNES iterations = ',i5)
196c4762a1bSJed Brown
197c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198c4762a1bSJed Brown!  Free work space.  All PETSc objects should be destroyed when they
199c4762a1bSJed Brown!  are no longer needed.
200c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
201c4762a1bSJed Brown
20265afa91aSSatish Balay      call VecDestroy(x,ierr)
20365afa91aSSatish Balay      CHKERRA(ierr)
20465afa91aSSatish Balay      call VecDestroy(r,ierr)
20565afa91aSSatish Balay      CHKERRA(ierr)
20665afa91aSSatish Balay      call SNESDestroy(snes,ierr)
20765afa91aSSatish Balay      CHKERRA(ierr)
20865afa91aSSatish Balay      call DMDestroy(da,ierr)
20965afa91aSSatish Balay      CHKERRA(ierr)
21065afa91aSSatish Balay      call PetscFinalize(ierr)
21165afa91aSSatish Balay      CHKERRA(ierr)
212c4762a1bSJed Brown      end
213c4762a1bSJed Brown
214c4762a1bSJed Brown! ---------------------------------------------------------------------
215c4762a1bSJed Brown!
216c4762a1bSJed Brown!  FormInitialGuess - Forms initial approximation.
217c4762a1bSJed Brown!
218c4762a1bSJed Brown!  Input Parameters:
219c4762a1bSJed Brown!  X - vector
220c4762a1bSJed Brown!
221c4762a1bSJed Brown!  Output Parameter:
222c4762a1bSJed Brown!  X - vector
223c4762a1bSJed Brown!
224c4762a1bSJed Brown!  Notes:
225c4762a1bSJed Brown!  This routine serves as a wrapper for the lower-level routine
226c4762a1bSJed Brown!  "ApplicationInitialGuess", where the actual computations are
227c4762a1bSJed Brown!  done using the standard Fortran style of treating the local
228c4762a1bSJed Brown!  vector data as a multidimensional array over the local mesh.
229c4762a1bSJed Brown!  This routine merely handles ghost point scatters and accesses
230*42ce371bSBarry Smith!  the local vector data via VecGetArrayF90() and VecRestoreArrayF90().
231c4762a1bSJed Brown!
232c4762a1bSJed Brown      subroutine FormInitialGuess(X,ierr)
233dfbbaf82SBarry Smith      use ex5fmodule
234c4762a1bSJed Brown      implicit none
235c4762a1bSJed Brown
236c4762a1bSJed Brown!  Input/output variables:
237c4762a1bSJed Brown      Vec      X
238c4762a1bSJed Brown      PetscErrorCode  ierr
239c4762a1bSJed Brown!  Declarations for use with local arrays:
240*42ce371bSBarry Smith      PetscScalar, pointer :: lx_v(:)
241c4762a1bSJed Brown
242c4762a1bSJed Brown      ierr = 0
243c4762a1bSJed Brown
244c4762a1bSJed Brown!  Get a pointer to vector data.
245c4762a1bSJed Brown!    - For default PETSc vectors, VecGetArray() returns a pointer to
246c4762a1bSJed Brown!      the data array.  Otherwise, the routine is implementation dependent.
247c4762a1bSJed Brown!    - You MUST call VecRestoreArray() when you no longer need access to
248c4762a1bSJed Brown!      the array.
249c4762a1bSJed Brown!    - Note that the Fortran interface to VecGetArray() differs from the
250c4762a1bSJed Brown!      C version.  See the users manual for details.
251c4762a1bSJed Brown
252*42ce371bSBarry Smith      call VecGetArrayF90(X,lx_v,ierr)
25365afa91aSSatish Balay      CHKERRQ(ierr)
254c4762a1bSJed Brown
255c4762a1bSJed Brown!  Compute initial guess over the locally owned part of the grid
256c4762a1bSJed Brown
257*42ce371bSBarry Smith      call InitialGuessLocal(lx_v,ierr)
25865afa91aSSatish Balay      CHKERRQ(ierr)
259c4762a1bSJed Brown
260c4762a1bSJed Brown!  Restore vector
261c4762a1bSJed Brown
262*42ce371bSBarry Smith      call VecRestoreArrayF90(X,lx_v,ierr)
26365afa91aSSatish Balay      CHKERRQ(ierr)
264c4762a1bSJed Brown
265c4762a1bSJed Brown      return
266c4762a1bSJed Brown      end
267c4762a1bSJed Brown
268c4762a1bSJed Brown! ---------------------------------------------------------------------
269c4762a1bSJed Brown!
270c4762a1bSJed Brown!  InitialGuessLocal - Computes initial approximation, called by
271c4762a1bSJed Brown!  the higher level routine FormInitialGuess().
272c4762a1bSJed Brown!
273c4762a1bSJed Brown!  Input Parameter:
274c4762a1bSJed Brown!  x - local vector data
275c4762a1bSJed Brown!
276c4762a1bSJed Brown!  Output Parameters:
277c4762a1bSJed Brown!  x - local vector data
278c4762a1bSJed Brown!  ierr - error code
279c4762a1bSJed Brown!
280c4762a1bSJed Brown!  Notes:
281c4762a1bSJed Brown!  This routine uses standard Fortran-style computations over a 2-dim array.
282c4762a1bSJed Brown!
283c4762a1bSJed Brown      subroutine InitialGuessLocal(x,ierr)
284dfbbaf82SBarry Smith      use ex5fmodule
285c4762a1bSJed Brown      implicit none
286c4762a1bSJed Brown
287c4762a1bSJed Brown!  Input/output variables:
288c4762a1bSJed Brown      PetscScalar    x(xs:xe,ys:ye)
289c4762a1bSJed Brown      PetscErrorCode ierr
290c4762a1bSJed Brown
291c4762a1bSJed Brown!  Local variables:
292c4762a1bSJed Brown      PetscInt  i,j
293c4762a1bSJed Brown      PetscReal temp1,temp,one,hx,hy
294c4762a1bSJed Brown
295c4762a1bSJed Brown!  Set parameters
296c4762a1bSJed Brown
297c4762a1bSJed Brown      ierr   = 0
298c4762a1bSJed Brown      one    = 1.0
299c4762a1bSJed Brown      hx     = one/((mx-1))
300c4762a1bSJed Brown      hy     = one/((my-1))
301c4762a1bSJed Brown      temp1  = lambda/(lambda + one)
302c4762a1bSJed Brown
303c4762a1bSJed Brown      do 20 j=ys,ye
304c4762a1bSJed Brown         temp = (min(j-1,my-j))*hy
305c4762a1bSJed Brown         do 10 i=xs,xe
306c4762a1bSJed Brown            if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then
307c4762a1bSJed Brown              x(i,j) = 0.0
308c4762a1bSJed Brown            else
309c4762a1bSJed Brown              x(i,j) = temp1 * sqrt(min(min(i-1,mx-i)*hx,(temp)))
310c4762a1bSJed Brown            endif
311c4762a1bSJed Brown 10      continue
312c4762a1bSJed Brown 20   continue
313c4762a1bSJed Brown
314c4762a1bSJed Brown      return
315c4762a1bSJed Brown      end
316c4762a1bSJed Brown
317c4762a1bSJed Brown! ---------------------------------------------------------------------
318c4762a1bSJed Brown!
319c4762a1bSJed Brown!  FormFunctionLocal - Computes nonlinear function, called by
320c4762a1bSJed Brown!  the higher level routine FormFunction().
321c4762a1bSJed Brown!
322c4762a1bSJed Brown!  Input Parameter:
323c4762a1bSJed Brown!  x - local vector data
324c4762a1bSJed Brown!
325c4762a1bSJed Brown!  Output Parameters:
326c4762a1bSJed Brown!  f - local vector data, f(x)
327c4762a1bSJed Brown!  ierr - error code
328c4762a1bSJed Brown!
329c4762a1bSJed Brown!  Notes:
330c4762a1bSJed Brown!  This routine uses standard Fortran-style computations over a 2-dim array.
331c4762a1bSJed Brown!
332c4762a1bSJed Brown!
333c4762a1bSJed Brown      subroutine FormFunctionLocal(info,x,f,da,ierr)
334dfbbaf82SBarry Smith      use ex5fmodule
335c4762a1bSJed Brown      implicit none
336c4762a1bSJed Brown
337c4762a1bSJed Brown      DM da
338c4762a1bSJed Brown
339c4762a1bSJed Brown!  Input/output variables:
340c4762a1bSJed Brown      DMDALocalInfo info(DMDA_LOCAL_INFO_SIZE)
341c4762a1bSJed Brown      PetscScalar x(gxs:gxe,gys:gye)
342c4762a1bSJed Brown      PetscScalar f(xs:xe,ys:ye)
343c4762a1bSJed Brown      PetscErrorCode     ierr
344c4762a1bSJed Brown
345c4762a1bSJed Brown!  Local variables:
346c4762a1bSJed Brown      PetscScalar two,one,hx,hy
347c4762a1bSJed Brown      PetscScalar hxdhy,hydhx,sc
348c4762a1bSJed Brown      PetscScalar u,uxx,uyy
349c4762a1bSJed Brown      PetscInt  i,j
350c4762a1bSJed Brown
351c4762a1bSJed Brown      xs     = info(DMDA_LOCAL_INFO_XS)+1
352c4762a1bSJed Brown      xe     = xs+info(DMDA_LOCAL_INFO_XM)-1
353c4762a1bSJed Brown      ys     = info(DMDA_LOCAL_INFO_YS)+1
354c4762a1bSJed Brown      ye     = ys+info(DMDA_LOCAL_INFO_YM)-1
355c4762a1bSJed Brown      mx     = info(DMDA_LOCAL_INFO_MX)
356c4762a1bSJed Brown      my     = info(DMDA_LOCAL_INFO_MY)
357c4762a1bSJed Brown
358c4762a1bSJed Brown      one    = 1.0
359c4762a1bSJed Brown      two    = 2.0
360c4762a1bSJed Brown      hx     = one/(mx-1)
361c4762a1bSJed Brown      hy     = one/(my-1)
362c4762a1bSJed Brown      sc     = hx*hy*lambda
363c4762a1bSJed Brown      hxdhy  = hx/hy
364c4762a1bSJed Brown      hydhx  = hy/hx
365c4762a1bSJed Brown
366c4762a1bSJed Brown!  Compute function over the locally owned part of the grid
367c4762a1bSJed Brown
368c4762a1bSJed Brown      do 20 j=ys,ye
369c4762a1bSJed Brown         do 10 i=xs,xe
370c4762a1bSJed Brown            if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then
371c4762a1bSJed Brown               f(i,j) = x(i,j)
372c4762a1bSJed Brown            else
373c4762a1bSJed Brown               u = x(i,j)
374c4762a1bSJed Brown               uxx = hydhx * (two*u - x(i-1,j) - x(i+1,j))
375c4762a1bSJed Brown               uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1))
376c4762a1bSJed Brown               f(i,j) = uxx + uyy - sc*exp(u)
377c4762a1bSJed Brown            endif
378c4762a1bSJed Brown 10      continue
379c4762a1bSJed Brown 20   continue
380c4762a1bSJed Brown
38165afa91aSSatish Balay      call PetscLogFlops(11.0d0*ym*xm,ierr)
38265afa91aSSatish Balay      CHKERRQ(ierr)
383c4762a1bSJed Brown
384c4762a1bSJed Brown      return
385c4762a1bSJed Brown      end
386c4762a1bSJed Brown
387c4762a1bSJed Brown! ---------------------------------------------------------------------
388c4762a1bSJed Brown!
389c4762a1bSJed Brown!  FormJacobianLocal - Computes Jacobian matrix, called by
390c4762a1bSJed Brown!  the higher level routine FormJacobian().
391c4762a1bSJed Brown!
392c4762a1bSJed Brown!  Input Parameters:
393c4762a1bSJed Brown!  x        - local vector data
394c4762a1bSJed Brown!
395c4762a1bSJed Brown!  Output Parameters:
396c4762a1bSJed Brown!  jac      - Jacobian matrix
397c4762a1bSJed Brown!  jac_prec - optionally different preconditioning matrix (not used here)
398c4762a1bSJed Brown!  ierr     - error code
399c4762a1bSJed Brown!
400c4762a1bSJed Brown!  Notes:
401c4762a1bSJed Brown!  This routine uses standard Fortran-style computations over a 2-dim array.
402c4762a1bSJed Brown!
403c4762a1bSJed Brown!  Notes:
404c4762a1bSJed Brown!  Due to grid point reordering with DMDAs, we must always work
405c4762a1bSJed Brown!  with the local grid points, and then transform them to the new
406c4762a1bSJed Brown!  global numbering with the "ltog" mapping
407c4762a1bSJed Brown!  We cannot work directly with the global numbers for the original
408c4762a1bSJed Brown!  uniprocessor grid!
409c4762a1bSJed Brown!
410c4762a1bSJed Brown!  Two methods are available for imposing this transformation
411c4762a1bSJed Brown!  when setting matrix entries:
412c4762a1bSJed Brown!    (A) MatSetValuesLocal(), using the local ordering (including
413c4762a1bSJed Brown!        ghost points!)
414c4762a1bSJed Brown!          by calling MatSetValuesLocal()
415c4762a1bSJed Brown!    (B) MatSetValues(), using the global ordering
416c4762a1bSJed Brown!        - Use DMDAGetGlobalIndices() to extract the local-to-global map
417c4762a1bSJed Brown!        - Then apply this map explicitly yourself
418c4762a1bSJed Brown!        - Set matrix entries using the global ordering by calling
419c4762a1bSJed Brown!          MatSetValues()
420c4762a1bSJed Brown!  Option (A) seems cleaner/easier in many cases, and is the procedure
421c4762a1bSJed Brown!  used in this example.
422c4762a1bSJed Brown!
423c4762a1bSJed Brown      subroutine FormJacobianLocal(info,x,A,jac,da,ierr)
424dfbbaf82SBarry Smith      use ex5fmodule
425c4762a1bSJed Brown      implicit none
426c4762a1bSJed Brown
427c4762a1bSJed Brown      DM da
428c4762a1bSJed Brown
429c4762a1bSJed Brown!  Input/output variables:
430c4762a1bSJed Brown      PetscScalar x(gxs:gxe,gys:gye)
431c4762a1bSJed Brown      Mat         A,jac
432c4762a1bSJed Brown      PetscErrorCode  ierr
433c4762a1bSJed Brown      DMDALocalInfo info(DMDA_LOCAL_INFO_SIZE)
434c4762a1bSJed Brown
435c4762a1bSJed Brown!  Local variables:
436c4762a1bSJed Brown      PetscInt  row,col(5),i,j,i1,i5
437c4762a1bSJed Brown      PetscScalar two,one,hx,hy,v(5)
438c4762a1bSJed Brown      PetscScalar hxdhy,hydhx,sc
439c4762a1bSJed Brown
440c4762a1bSJed Brown!  Set parameters
441c4762a1bSJed Brown
442c4762a1bSJed Brown      i1     = 1
443c4762a1bSJed Brown      i5     = 5
444c4762a1bSJed Brown      one    = 1.0
445c4762a1bSJed Brown      two    = 2.0
446c4762a1bSJed Brown      hx     = one/(mx-1)
447c4762a1bSJed Brown      hy     = one/(my-1)
448c4762a1bSJed Brown      sc     = hx*hy
449c4762a1bSJed Brown      hxdhy  = hx/hy
450c4762a1bSJed Brown      hydhx  = hy/hx
451c4762a1bSJed Brown
452c4762a1bSJed Brown!  Compute entries for the locally owned part of the Jacobian.
453c4762a1bSJed Brown!   - Currently, all PETSc parallel matrix formats are partitioned by
454c4762a1bSJed Brown!     contiguous chunks of rows across the processors.
455c4762a1bSJed Brown!   - Each processor needs to insert only elements that it owns
456c4762a1bSJed Brown!     locally (but any non-local elements will be sent to the
457c4762a1bSJed Brown!     appropriate processor during matrix assembly).
458c4762a1bSJed Brown!   - Here, we set all entries for a particular row at once.
459c4762a1bSJed Brown!   - We can set matrix entries either using either
460c4762a1bSJed Brown!     MatSetValuesLocal() or MatSetValues(), as discussed above.
461c4762a1bSJed Brown!   - Note that MatSetValues() uses 0-based row and column numbers
462c4762a1bSJed Brown!     in Fortran as well as in C.
463c4762a1bSJed Brown
464c4762a1bSJed Brown      do 20 j=ys,ye
465c4762a1bSJed Brown         row = (j - gys)*gxm + xs - gxs - 1
466c4762a1bSJed Brown         do 10 i=xs,xe
467c4762a1bSJed Brown            row = row + 1
468c4762a1bSJed Brown!           boundary points
469c4762a1bSJed Brown            if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then
470c4762a1bSJed Brown!       Some f90 compilers need 4th arg to be of same type in both calls
471c4762a1bSJed Brown               col(1) = row
472c4762a1bSJed Brown               v(1)   = one
47365afa91aSSatish Balay               call MatSetValuesLocal(jac,i1,row,i1,col,v,INSERT_VALUES,ierr)
47465afa91aSSatish Balay               CHKERRQ(ierr)
475c4762a1bSJed Brown!           interior grid points
476c4762a1bSJed Brown            else
477c4762a1bSJed Brown               v(1) = -hxdhy
478c4762a1bSJed Brown               v(2) = -hydhx
479c4762a1bSJed Brown               v(3) = two*(hydhx + hxdhy) - sc*lambda*exp(x(i,j))
480c4762a1bSJed Brown               v(4) = -hydhx
481c4762a1bSJed Brown               v(5) = -hxdhy
482c4762a1bSJed Brown               col(1) = row - gxm
483c4762a1bSJed Brown               col(2) = row - 1
484c4762a1bSJed Brown               col(3) = row
485c4762a1bSJed Brown               col(4) = row + 1
486c4762a1bSJed Brown               col(5) = row + gxm
48765afa91aSSatish Balay               call MatSetValuesLocal(jac,i1,row,i5,col,v, INSERT_VALUES,ierr)
48865afa91aSSatish Balay               CHKERRQ(ierr)
489c4762a1bSJed Brown            endif
490c4762a1bSJed Brown 10      continue
491c4762a1bSJed Brown 20   continue
49265afa91aSSatish Balay      call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
49365afa91aSSatish Balay      CHKERRQ(ierr)
49465afa91aSSatish Balay      call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
49565afa91aSSatish Balay      CHKERRQ(ierr)
496c4762a1bSJed Brown      if (A .ne. jac) then
49765afa91aSSatish Balay         call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
49865afa91aSSatish Balay         CHKERRQ(ierr)
49965afa91aSSatish Balay         call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)
50065afa91aSSatish Balay         CHKERRQ(ierr)
501c4762a1bSJed Brown      endif
502c4762a1bSJed Brown      return
503c4762a1bSJed Brown      end
504c4762a1bSJed Brown
505c4762a1bSJed Brown!
506c4762a1bSJed Brown!     Simple convergence test based on the infinity norm of the residual being small
507c4762a1bSJed Brown!
508c4762a1bSJed Brown      subroutine MySNESConverged(snes,it,xnorm,snorm,fnorm,reason,dummy,ierr)
509dfbbaf82SBarry Smith      use ex5fmodule
510c4762a1bSJed Brown      implicit none
511c4762a1bSJed Brown
512c4762a1bSJed Brown      SNES snes
513c4762a1bSJed Brown      PetscInt it,dummy
514c4762a1bSJed Brown      PetscReal xnorm,snorm,fnorm,nrm
515c4762a1bSJed Brown      SNESConvergedReason reason
516c4762a1bSJed Brown      Vec f
517c4762a1bSJed Brown      PetscErrorCode ierr
518c4762a1bSJed Brown
51965afa91aSSatish Balay      call SNESGetFunction(snes,f,PETSC_NULL_FUNCTION,dummy,ierr)
52065afa91aSSatish Balay      CHKERRQ(ierr)
52165afa91aSSatish Balay      call VecNorm(f,NORM_INFINITY,nrm,ierr)
52265afa91aSSatish Balay      CHKERRQ(ierr)
523c4762a1bSJed Brown      if (nrm .le. 1.e-5) reason = SNES_CONVERGED_FNORM_ABS
524c4762a1bSJed Brown
525c4762a1bSJed Brown      end
526c4762a1bSJed Brown
527c4762a1bSJed Brown!/*TEST
528c4762a1bSJed Brown!
529c4762a1bSJed Brown!   build:
530c4762a1bSJed Brown!      requires: !complex !single
531c4762a1bSJed Brown!
532c4762a1bSJed Brown!   test:
533c4762a1bSJed Brown!      nsize: 4
5348f8b3c79SStefano Zampini!      args: -snes_mf -pc_type none -da_processors_x 4 -da_processors_y 1 -snes_monitor_short \
5358f8b3c79SStefano Zampini!            -ksp_gmres_cgs_refinement_type refine_always
536c4762a1bSJed Brown!
537c4762a1bSJed Brown!   test:
538c4762a1bSJed Brown!      suffix: 2
539c4762a1bSJed Brown!      nsize: 4
540c4762a1bSJed Brown!      args: -da_processors_x 2 -da_processors_y 2 -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
541c4762a1bSJed Brown!
542c4762a1bSJed Brown!   test:
543c4762a1bSJed Brown!      suffix: 3
544c4762a1bSJed Brown!      nsize: 3
545c4762a1bSJed Brown!      args: -snes_fd -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
546c4762a1bSJed Brown!
547c4762a1bSJed Brown!   test:
548c4762a1bSJed Brown!      suffix: 6
549c4762a1bSJed Brown!      nsize: 1
550c4762a1bSJed Brown!      args: -snes_monitor_short -my_snes_convergence
551c4762a1bSJed Brown!
552c4762a1bSJed Brown!TEST*/
553