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