xref: /petsc/src/snes/tests/ex1f.F90 (revision d8606c274c09e255c003062beb17b1be973467bc)
1c4762a1bSJed Brown!
2c4762a1bSJed Brown!  Description: This example solves a nonlinear system on 1 processor with SNES.
3c4762a1bSJed Brown!  We solve the  Bratu (SFI - solid fuel ignition) problem in a 2D rectangular
4c4762a1bSJed Brown!  domain.  The command line options include:
5c4762a1bSJed Brown!    -par <parameter>, where <parameter> indicates the nonlinearity of the problem
6c4762a1bSJed Brown!       problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)
7c4762a1bSJed Brown!    -mx <xg>, where <xg> = number of grid points in the x-direction
8c4762a1bSJed Brown!    -my <yg>, where <yg> = number of grid points in the y-direction
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!  The parallel version of this code is snes/tutorials/ex5f.F
28c4762a1bSJed Brown!
29c4762a1bSJed Brown!  --------------------------------------------------------------------------
30c4762a1bSJed Brown      subroutine postcheck(snes,x,y,w,changed_y,changed_w,ctx,ierr)
31c4762a1bSJed Brown#include <petsc/finclude/petscsnes.h>
32c4762a1bSJed Brown      use petscsnes
33c4762a1bSJed Brown      implicit none
34c4762a1bSJed Brown      SNES           snes
35c4762a1bSJed Brown      PetscReal      norm
36c4762a1bSJed Brown      Vec            tmp,x,y,w
37c4762a1bSJed Brown      PetscBool      changed_w,changed_y
38c4762a1bSJed Brown      PetscErrorCode ierr
39c4762a1bSJed Brown      PetscInt       ctx
40c4762a1bSJed Brown      PetscScalar    mone
41c4762a1bSJed Brown
42*d8606c27SBarry Smith      PetscCallA(VecDuplicate(x,tmp,ierr))
43c4762a1bSJed Brown      mone = -1.0
44*d8606c27SBarry Smith      PetscCallA(VecWAXPY(tmp,mone,x,w,ierr))
45*d8606c27SBarry Smith      PetscCallA(VecNorm(tmp,NORM_2,norm,ierr))
46*d8606c27SBarry Smith      PetscCallA(VecDestroy(tmp,ierr))
47c4762a1bSJed Brown      print*, 'Norm of search step ',norm
48c4762a1bSJed Brown      changed_y = PETSC_FALSE
49c4762a1bSJed Brown      changed_w = PETSC_FALSE
50c4762a1bSJed Brown      return
51c4762a1bSJed Brown      end
52c4762a1bSJed Brown
53c4762a1bSJed Brown      program main
54c4762a1bSJed Brown#include <petsc/finclude/petscdraw.h>
55c4762a1bSJed Brown      use petscsnes
56c4762a1bSJed Brown      implicit none
5717a42bb7SSatish Balay      interface SNESSetJacobian
5817a42bb7SSatish Balay      subroutine SNESSetJacobian1(a,b,c,d,e,z)
5917a42bb7SSatish Balay       use petscsnes
6017a42bb7SSatish Balay       SNES a
6117a42bb7SSatish Balay       Mat b
6217a42bb7SSatish Balay       Mat c
6317a42bb7SSatish Balay       external d
6417a42bb7SSatish Balay       MatFDColoring e
6517a42bb7SSatish Balay       PetscErrorCode z
6617a42bb7SSatish Balay      end subroutine
6717a42bb7SSatish Balay      subroutine SNESSetJacobian2(a,b,c,d,e,z)
6817a42bb7SSatish Balay       use petscsnes
6917a42bb7SSatish Balay       SNES a
7017a42bb7SSatish Balay       Mat b
7117a42bb7SSatish Balay       Mat c
7217a42bb7SSatish Balay       external d
7317a42bb7SSatish Balay       integer e
7417a42bb7SSatish Balay       PetscErrorCode z
7517a42bb7SSatish Balay      end subroutine
7617a42bb7SSatish Balay      end interface
77c4762a1bSJed Brown!
78c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79c4762a1bSJed Brown!                   Variable declarations
80c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
81c4762a1bSJed Brown!
82c4762a1bSJed Brown!  Variables:
83c4762a1bSJed Brown!     snes        - nonlinear solver
84c4762a1bSJed Brown!     x, r        - solution, residual vectors
85c4762a1bSJed Brown!     J           - Jacobian matrix
86c4762a1bSJed Brown!     its         - iterations for convergence
87c4762a1bSJed Brown!     matrix_free - flag - 1 indicates matrix-free version
88c4762a1bSJed Brown!     lambda      - nonlinearity parameter
89c4762a1bSJed Brown!     draw        - drawing context
90c4762a1bSJed Brown!
91c4762a1bSJed Brown      SNES               snes
92c4762a1bSJed Brown      MatColoring        mc
93c4762a1bSJed Brown      Vec                x,r
94c4762a1bSJed Brown      PetscDraw               draw
95c4762a1bSJed Brown      Mat                J
96c4762a1bSJed Brown      PetscBool  matrix_free,flg,fd_coloring
97c4762a1bSJed Brown      PetscErrorCode ierr
98c4762a1bSJed Brown      PetscInt   its,N, mx,my,i5
99c4762a1bSJed Brown      PetscMPIInt size,rank
100c4762a1bSJed Brown      PetscReal   lambda_max,lambda_min,lambda
101c4762a1bSJed Brown      MatFDColoring      fdcoloring
102c4762a1bSJed Brown      ISColoring         iscoloring
103c4762a1bSJed Brown      PetscBool          pc
104c4762a1bSJed Brown      external           postcheck
105c4762a1bSJed Brown
106c4762a1bSJed Brown      PetscScalar        lx_v(0:1)
107c4762a1bSJed Brown      PetscOffset        lx_i
108c4762a1bSJed Brown
109c4762a1bSJed Brown!  Store parameters in common block
110c4762a1bSJed Brown
111c4762a1bSJed Brown      common /params/ lambda,mx,my,fd_coloring
112c4762a1bSJed Brown
113c4762a1bSJed Brown!  Note: Any user-defined Fortran routines (such as FormJacobian)
114c4762a1bSJed Brown!  MUST be declared as external.
115c4762a1bSJed Brown
116c4762a1bSJed Brown      external FormFunction,FormInitialGuess,FormJacobian
117c4762a1bSJed Brown
118c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
119c4762a1bSJed Brown!  Initialize program
120c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121c4762a1bSJed Brown
122*d8606c27SBarry Smith      PetscCallA(PetscInitialize(ierr))
123*d8606c27SBarry Smith      PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD,size,ierr))
124*d8606c27SBarry Smith      PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr))
125c4762a1bSJed Brown
126c4762a1bSJed Brown      if (size .ne. 1) then; SETERRA(PETSC_COMM_SELF,PETSC_ERR_WRONG_MPI_SIZE,'This is a uniprocessor example only'); endif
127c4762a1bSJed Brown
128c4762a1bSJed Brown!  Initialize problem parameters
129c4762a1bSJed Brown      i5 = 5
130c4762a1bSJed Brown      lambda_max = 6.81
131c4762a1bSJed Brown      lambda_min = 0.0
132c4762a1bSJed Brown      lambda     = 6.0
133c4762a1bSJed Brown      mx         = 4
134c4762a1bSJed Brown      my         = 4
135*d8606c27SBarry Smith      PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-mx',mx,flg,ierr))
136*d8606c27SBarry Smith      PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-my',my,flg,ierr))
137*d8606c27SBarry Smith      PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-par',lambda,flg,ierr))
138c4762a1bSJed Brown      if (lambda .ge. lambda_max .or. lambda .le. lambda_min) then; SETERRA(PETSC_COMM_SELF,PETSC_ERR_USER,'Lambda out of range '); endif
139c4762a1bSJed Brown      N  = mx*my
140c4762a1bSJed Brown      pc = PETSC_FALSE
141*d8606c27SBarry Smith      PetscCallA(PetscOptionsGetBool(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-pc',pc,PETSC_NULL_BOOL,ierr))
142c4762a1bSJed Brown
143c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144c4762a1bSJed Brown!  Create nonlinear solver context
145c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146c4762a1bSJed Brown
147*d8606c27SBarry Smith      PetscCallA(SNESCreate(PETSC_COMM_WORLD,snes,ierr))
148c4762a1bSJed Brown
149c4762a1bSJed Brown      if (pc .eqv. PETSC_TRUE) then
150*d8606c27SBarry Smith        PetscCallA(SNESSetType(snes,SNESNEWTONTR,ierr))
151*d8606c27SBarry Smith        PetscCallA(SNESNewtonTRSetPostCheck(snes, postcheck,snes,ierr))
152c4762a1bSJed Brown      endif
153c4762a1bSJed Brown
154c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
155c4762a1bSJed Brown!  Create vector data structures; set function evaluation routine
156c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157c4762a1bSJed Brown
158*d8606c27SBarry Smith      PetscCallA(VecCreate(PETSC_COMM_WORLD,x,ierr))
159*d8606c27SBarry Smith      PetscCallA(VecSetSizes(x,PETSC_DECIDE,N,ierr))
160*d8606c27SBarry Smith      PetscCallA(VecSetFromOptions(x,ierr))
161*d8606c27SBarry Smith      PetscCallA(VecDuplicate(x,r,ierr))
162c4762a1bSJed Brown
163c4762a1bSJed Brown!  Set function evaluation routine and vector.  Whenever the nonlinear
164c4762a1bSJed Brown!  solver needs to evaluate the nonlinear function, it will call this
165c4762a1bSJed Brown!  routine.
166c4762a1bSJed Brown!   - Note that the final routine argument is the user-defined
167c4762a1bSJed Brown!     context that provides application-specific data for the
168c4762a1bSJed Brown!     function evaluation routine.
169c4762a1bSJed Brown
170*d8606c27SBarry Smith      PetscCallA(SNESSetFunction(snes,r,FormFunction,fdcoloring,ierr))
171c4762a1bSJed Brown
172c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
173c4762a1bSJed Brown!  Create matrix data structure; set Jacobian evaluation routine
174c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
175c4762a1bSJed Brown
176c4762a1bSJed Brown!  Create matrix. Here we only approximately preallocate storage space
177c4762a1bSJed Brown!  for the Jacobian.  See the users manual for a discussion of better
178c4762a1bSJed Brown!  techniques for preallocating matrix memory.
179c4762a1bSJed Brown
180*d8606c27SBarry Smith      PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-snes_mf',matrix_free,ierr))
181c4762a1bSJed Brown      if (.not. matrix_free) then
182*d8606c27SBarry Smith        PetscCallA(MatCreateSeqAIJ(PETSC_COMM_WORLD,N,N,i5,PETSC_NULL_INTEGER,J,ierr))
183c4762a1bSJed Brown      endif
184c4762a1bSJed Brown
185c4762a1bSJed Brown!
186c4762a1bSJed Brown!     This option will cause the Jacobian to be computed via finite differences
187c4762a1bSJed Brown!    efficiently using a coloring of the columns of the matrix.
188c4762a1bSJed Brown!
189c4762a1bSJed Brown      fd_coloring = .false.
190*d8606c27SBarry Smith      PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-snes_fd_coloring',fd_coloring,ierr))
191c4762a1bSJed Brown      if (fd_coloring) then
192c4762a1bSJed Brown
193c4762a1bSJed Brown!
194c4762a1bSJed Brown!      This initializes the nonzero structure of the Jacobian. This is artificial
195c4762a1bSJed Brown!      because clearly if we had a routine to compute the Jacobian we won't need
196c4762a1bSJed Brown!      to use finite differences.
197c4762a1bSJed Brown!
198*d8606c27SBarry Smith        PetscCallA(FormJacobian(snes,x,J,J,0,ierr))
199c4762a1bSJed Brown!
200c4762a1bSJed Brown!       Color the matrix, i.e. determine groups of columns that share no common
201a5b23f4aSJose E. Roman!      rows. These columns in the Jacobian can all be computed simultaneously.
202c4762a1bSJed Brown!
203*d8606c27SBarry Smith        PetscCallA(MatColoringCreate(J,mc,ierr))
204*d8606c27SBarry Smith        PetscCallA(MatColoringSetType(mc,MATCOLORINGNATURAL,ierr))
205*d8606c27SBarry Smith        PetscCallA(MatColoringSetFromOptions(mc,ierr))
206*d8606c27SBarry Smith        PetscCallA(MatColoringApply(mc,iscoloring,ierr))
207*d8606c27SBarry Smith        PetscCallA(MatColoringDestroy(mc,ierr))
208c4762a1bSJed Brown!
209c4762a1bSJed Brown!       Create the data structure that SNESComputeJacobianDefaultColor() uses
210c4762a1bSJed Brown!       to compute the actual Jacobians via finite differences.
211c4762a1bSJed Brown!
212*d8606c27SBarry Smith        PetscCallA(MatFDColoringCreate(J,iscoloring,fdcoloring,ierr))
213*d8606c27SBarry Smith        PetscCallA(MatFDColoringSetFunction(fdcoloring,FormFunction,fdcoloring,ierr))
214*d8606c27SBarry Smith        PetscCallA(MatFDColoringSetFromOptions(fdcoloring,ierr))
215*d8606c27SBarry Smith        PetscCallA(MatFDColoringSetUp(J,iscoloring,fdcoloring,ierr))
216c4762a1bSJed Brown!
217c4762a1bSJed Brown!        Tell SNES to use the routine SNESComputeJacobianDefaultColor()
218c4762a1bSJed Brown!      to compute Jacobians.
219c4762a1bSJed Brown!
220*d8606c27SBarry Smith        PetscCallA(SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,fdcoloring,ierr))
221*d8606c27SBarry Smith        PetscCallA(ISColoringDestroy(iscoloring,ierr))
222c4762a1bSJed Brown
223c4762a1bSJed Brown      else if (.not. matrix_free) then
224c4762a1bSJed Brown
225c4762a1bSJed Brown!  Set Jacobian matrix data structure and default Jacobian evaluation
226c4762a1bSJed Brown!  routine.  Whenever the nonlinear solver needs to compute the
227c4762a1bSJed Brown!  Jacobian matrix, it will call this routine.
228c4762a1bSJed Brown!   - Note that the final routine argument is the user-defined
229c4762a1bSJed Brown!     context that provides application-specific data for the
230c4762a1bSJed Brown!     Jacobian evaluation routine.
231c4762a1bSJed Brown!   - The user can override with:
232c4762a1bSJed Brown!      -snes_fd : default finite differencing approximation of Jacobian
233c4762a1bSJed Brown!      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
234c4762a1bSJed Brown!                 (unless user explicitly sets preconditioner)
235c4762a1bSJed Brown!      -snes_mf_operator : form preconditioning matrix as set by the user,
236c4762a1bSJed Brown!                          but use matrix-free approx for Jacobian-vector
237c4762a1bSJed Brown!                          products within Newton-Krylov method
238c4762a1bSJed Brown!
239*d8606c27SBarry Smith        PetscCallA(SNESSetJacobian(snes,J,J,FormJacobian,0,ierr))
240c4762a1bSJed Brown      endif
241c4762a1bSJed Brown
242c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
243c4762a1bSJed Brown!  Customize nonlinear solver; set runtime options
244c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
245c4762a1bSJed Brown
246c4762a1bSJed Brown!  Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
247c4762a1bSJed Brown
248*d8606c27SBarry Smith      PetscCallA(SNESSetFromOptions(snes,ierr))
249c4762a1bSJed Brown
250c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
251c4762a1bSJed Brown!  Evaluate initial guess; then solve nonlinear system.
252c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
253c4762a1bSJed Brown
254c4762a1bSJed Brown!  Note: The user should initialize the vector, x, with the initial guess
255c4762a1bSJed Brown!  for the nonlinear solver prior to calling SNESSolve().  In particular,
256c4762a1bSJed Brown!  to employ an initial guess of zero, the user should explicitly set
257c4762a1bSJed Brown!  this vector to zero by calling VecSet().
258c4762a1bSJed Brown
259*d8606c27SBarry Smith      PetscCallA(FormInitialGuess(x,ierr))
260*d8606c27SBarry Smith      PetscCallA(SNESSolve(snes,PETSC_NULL_VEC,x,ierr))
261*d8606c27SBarry Smith      PetscCallA(SNESGetIterationNumber(snes,its,ierr))
262c4762a1bSJed Brown      if (rank .eq. 0) then
263c4762a1bSJed Brown         write(6,100) its
264c4762a1bSJed Brown      endif
265c4762a1bSJed Brown  100 format('Number of SNES iterations = ',i1)
266c4762a1bSJed Brown
267c4762a1bSJed Brown!  PetscDraw contour plot of solution
268c4762a1bSJed Brown
269*d8606c27SBarry Smith      PetscCallA(PetscDrawCreate(PETSC_COMM_WORLD,PETSC_NULL_CHARACTER,'Solution',300,0,300,300,draw,ierr))
270*d8606c27SBarry Smith      PetscCallA(PetscDrawSetFromOptions(draw,ierr))
271c4762a1bSJed Brown
272*d8606c27SBarry Smith      PetscCallA(VecGetArrayRead(x,lx_v,lx_i,ierr))
273*d8606c27SBarry Smith      PetscCallA(PetscDrawTensorContour(draw,mx,my,PETSC_NULL_REAL,PETSC_NULL_REAL,lx_v(lx_i+1),ierr))
274*d8606c27SBarry Smith      PetscCallA(VecRestoreArrayRead(x,lx_v,lx_i,ierr))
275c4762a1bSJed Brown
276c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
277c4762a1bSJed Brown!  Free work space.  All PETSc objects should be destroyed when they
278c4762a1bSJed Brown!  are no longer needed.
279c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
280c4762a1bSJed Brown
281*d8606c27SBarry Smith      if (.not. matrix_free) PetscCallA(MatDestroy(J,ierr))
282*d8606c27SBarry Smith      if (fd_coloring) PetscCallA(MatFDColoringDestroy(fdcoloring,ierr))
283c4762a1bSJed Brown
284*d8606c27SBarry Smith      PetscCallA(VecDestroy(x,ierr))
285*d8606c27SBarry Smith      PetscCallA(VecDestroy(r,ierr))
286*d8606c27SBarry Smith      PetscCallA(SNESDestroy(snes,ierr))
287*d8606c27SBarry Smith      PetscCallA(PetscDrawDestroy(draw,ierr))
288*d8606c27SBarry Smith      PetscCallA(PetscFinalize(ierr))
289c4762a1bSJed Brown      end
290c4762a1bSJed Brown
291c4762a1bSJed Brown! ---------------------------------------------------------------------
292c4762a1bSJed Brown!
293c4762a1bSJed Brown!  FormInitialGuess - Forms initial approximation.
294c4762a1bSJed Brown!
295c4762a1bSJed Brown!  Input Parameter:
296c4762a1bSJed Brown!  X - vector
297c4762a1bSJed Brown!
298c4762a1bSJed Brown!  Output Parameters:
299c4762a1bSJed Brown!  X - vector
300c4762a1bSJed Brown!  ierr - error code
301c4762a1bSJed Brown!
302c4762a1bSJed Brown!  Notes:
303c4762a1bSJed Brown!  This routine serves as a wrapper for the lower-level routine
304c4762a1bSJed Brown!  "ApplicationInitialGuess", where the actual computations are
305c4762a1bSJed Brown!  done using the standard Fortran style of treating the local
306c4762a1bSJed Brown!  vector data as a multidimensional array over the local mesh.
307c4762a1bSJed Brown!  This routine merely accesses the local vector data via
308c4762a1bSJed Brown!  VecGetArray() and VecRestoreArray().
309c4762a1bSJed Brown!
310c4762a1bSJed Brown      subroutine FormInitialGuess(X,ierr)
311c4762a1bSJed Brown      use petscsnes
312c4762a1bSJed Brown      implicit none
313c4762a1bSJed Brown
314c4762a1bSJed Brown!  Input/output variables:
315c4762a1bSJed Brown      Vec           X
316c4762a1bSJed Brown      PetscErrorCode    ierr
317c4762a1bSJed Brown
318c4762a1bSJed Brown!  Declarations for use with local arrays:
319c4762a1bSJed Brown      PetscScalar   lx_v(0:1)
320c4762a1bSJed Brown      PetscOffset   lx_i
321c4762a1bSJed Brown
322c4762a1bSJed Brown      ierr   = 0
323c4762a1bSJed Brown
324c4762a1bSJed Brown!  Get a pointer to vector data.
325c4762a1bSJed Brown!    - For default PETSc vectors, VecGetArray() returns a pointer to
326c4762a1bSJed Brown!      the data array.  Otherwise, the routine is implementation dependent.
327c4762a1bSJed Brown!    - You MUST call VecRestoreArray() when you no longer need access to
328c4762a1bSJed Brown!      the array.
329c4762a1bSJed Brown!    - Note that the Fortran interface to VecGetArray() differs from the
330c4762a1bSJed Brown!      C version.  See the users manual for details.
331c4762a1bSJed Brown
332*d8606c27SBarry Smith      PetscCallA(VecGetArray(X,lx_v,lx_i,ierr))
333c4762a1bSJed Brown
334c4762a1bSJed Brown!  Compute initial guess
335c4762a1bSJed Brown
336*d8606c27SBarry Smith      PetscCallA(ApplicationInitialGuess(lx_v(lx_i),ierr))
337c4762a1bSJed Brown
338c4762a1bSJed Brown!  Restore vector
339c4762a1bSJed Brown
340*d8606c27SBarry Smith      PetscCallA(VecRestoreArray(X,lx_v,lx_i,ierr))
341c4762a1bSJed Brown
342c4762a1bSJed Brown      return
343c4762a1bSJed Brown      end
344c4762a1bSJed Brown
345c4762a1bSJed Brown! ---------------------------------------------------------------------
346c4762a1bSJed Brown!
347c4762a1bSJed Brown!  ApplicationInitialGuess - Computes initial approximation, called by
348c4762a1bSJed Brown!  the higher level routine FormInitialGuess().
349c4762a1bSJed Brown!
350c4762a1bSJed Brown!  Input Parameter:
351c4762a1bSJed Brown!  x - local vector data
352c4762a1bSJed Brown!
353c4762a1bSJed Brown!  Output Parameters:
354c4762a1bSJed Brown!  f - local vector data, f(x)
355c4762a1bSJed Brown!  ierr - error code
356c4762a1bSJed Brown!
357c4762a1bSJed Brown!  Notes:
358c4762a1bSJed Brown!  This routine uses standard Fortran-style computations over a 2-dim array.
359c4762a1bSJed Brown!
360c4762a1bSJed Brown      subroutine ApplicationInitialGuess(x,ierr)
361c4762a1bSJed Brown      use petscksp
362c4762a1bSJed Brown      implicit none
363c4762a1bSJed Brown
364c4762a1bSJed Brown!  Common blocks:
365c4762a1bSJed Brown      PetscReal   lambda
366c4762a1bSJed Brown      PetscInt     mx,my
367c4762a1bSJed Brown      PetscBool         fd_coloring
368c4762a1bSJed Brown      common      /params/ lambda,mx,my,fd_coloring
369c4762a1bSJed Brown
370c4762a1bSJed Brown!  Input/output variables:
371c4762a1bSJed Brown      PetscScalar x(mx,my)
372c4762a1bSJed Brown      PetscErrorCode     ierr
373c4762a1bSJed Brown
374c4762a1bSJed Brown!  Local variables:
375c4762a1bSJed Brown      PetscInt     i,j
376c4762a1bSJed Brown      PetscReal temp1,temp,hx,hy,one
377c4762a1bSJed Brown
378c4762a1bSJed Brown!  Set parameters
379c4762a1bSJed Brown
380c4762a1bSJed Brown      ierr   = 0
381c4762a1bSJed Brown      one    = 1.0
382c4762a1bSJed Brown      hx     = one/(mx-1)
383c4762a1bSJed Brown      hy     = one/(my-1)
384c4762a1bSJed Brown      temp1  = lambda/(lambda + one)
385c4762a1bSJed Brown
386c4762a1bSJed Brown      do 20 j=1,my
387c4762a1bSJed Brown         temp = min(j-1,my-j)*hy
388c4762a1bSJed Brown         do 10 i=1,mx
389c4762a1bSJed Brown            if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then
390c4762a1bSJed Brown              x(i,j) = 0.0
391c4762a1bSJed Brown            else
392c4762a1bSJed Brown              x(i,j) = temp1 * sqrt(min(min(i-1,mx-i)*hx,temp))
393c4762a1bSJed Brown            endif
394c4762a1bSJed Brown 10      continue
395c4762a1bSJed Brown 20   continue
396c4762a1bSJed Brown
397c4762a1bSJed Brown      return
398c4762a1bSJed Brown      end
399c4762a1bSJed Brown
400c4762a1bSJed Brown! ---------------------------------------------------------------------
401c4762a1bSJed Brown!
402c4762a1bSJed Brown!  FormFunction - Evaluates nonlinear function, F(x).
403c4762a1bSJed Brown!
404c4762a1bSJed Brown!  Input Parameters:
405c4762a1bSJed Brown!  snes  - the SNES context
406c4762a1bSJed Brown!  X     - input vector
407c4762a1bSJed Brown!  dummy - optional user-defined context, as set by SNESSetFunction()
408c4762a1bSJed Brown!          (not used here)
409c4762a1bSJed Brown!
410c4762a1bSJed Brown!  Output Parameter:
411c4762a1bSJed Brown!  F     - vector with newly computed function
412c4762a1bSJed Brown!
413c4762a1bSJed Brown!  Notes:
414c4762a1bSJed Brown!  This routine serves as a wrapper for the lower-level routine
415c4762a1bSJed Brown!  "ApplicationFunction", where the actual computations are
416c4762a1bSJed Brown!  done using the standard Fortran style of treating the local
417c4762a1bSJed Brown!  vector data as a multidimensional array over the local mesh.
418c4762a1bSJed Brown!  This routine merely accesses the local vector data via
419c4762a1bSJed Brown!  VecGetArray() and VecRestoreArray().
420c4762a1bSJed Brown!
421c4762a1bSJed Brown      subroutine FormFunction(snes,X,F,fdcoloring,ierr)
422c4762a1bSJed Brown      use petscsnes
423c4762a1bSJed Brown      implicit none
424c4762a1bSJed Brown
425c4762a1bSJed Brown!  Input/output variables:
426c4762a1bSJed Brown      SNES              snes
427c4762a1bSJed Brown      Vec               X,F
428c4762a1bSJed Brown      PetscErrorCode          ierr
429c4762a1bSJed Brown      MatFDColoring fdcoloring
430c4762a1bSJed Brown
431c4762a1bSJed Brown!  Common blocks:
432c4762a1bSJed Brown      PetscReal         lambda
433c4762a1bSJed Brown      PetscInt          mx,my
434c4762a1bSJed Brown      PetscBool         fd_coloring
435c4762a1bSJed Brown      common            /params/ lambda,mx,my,fd_coloring
436c4762a1bSJed Brown
437c4762a1bSJed Brown!  Declarations for use with local arrays:
438c4762a1bSJed Brown      PetscScalar       lx_v(0:1),lf_v(0:1)
439c4762a1bSJed Brown      PetscOffset       lx_i,lf_i
440c4762a1bSJed Brown
441c4762a1bSJed Brown      PetscInt, pointer :: indices(:)
442c4762a1bSJed Brown
443c4762a1bSJed Brown!  Get pointers to vector data.
444c4762a1bSJed Brown!    - For default PETSc vectors, VecGetArray() returns a pointer to
445c4762a1bSJed Brown!      the data array.  Otherwise, the routine is implementation dependent.
446c4762a1bSJed Brown!    - You MUST call VecRestoreArray() when you no longer need access to
447c4762a1bSJed Brown!      the array.
448c4762a1bSJed Brown!    - Note that the Fortran interface to VecGetArray() differs from the
449c4762a1bSJed Brown!      C version.  See the Fortran chapter of the users manual for details.
450c4762a1bSJed Brown
451*d8606c27SBarry Smith      PetscCallA(VecGetArrayRead(X,lx_v,lx_i,ierr))
452*d8606c27SBarry Smith      PetscCallA(VecGetArray(F,lf_v,lf_i,ierr))
453c4762a1bSJed Brown
454c4762a1bSJed Brown!  Compute function
455c4762a1bSJed Brown
456*d8606c27SBarry Smith      PetscCallA(ApplicationFunction(lx_v(lx_i),lf_v(lf_i),ierr))
457c4762a1bSJed Brown
458c4762a1bSJed Brown!  Restore vectors
459c4762a1bSJed Brown
460*d8606c27SBarry Smith      PetscCallA(VecRestoreArrayRead(X,lx_v,lx_i,ierr))
461*d8606c27SBarry Smith      PetscCallA(VecRestoreArray(F,lf_v,lf_i,ierr))
462c4762a1bSJed Brown
463*d8606c27SBarry Smith      PetscCallA(PetscLogFlops(11.0d0*mx*my,ierr))
464c4762a1bSJed Brown!
465c4762a1bSJed Brown!     fdcoloring is in the common block and used here ONLY to test the
466c4762a1bSJed Brown!     calls to MatFDColoringGetPerturbedColumnsF90() and  MatFDColoringRestorePerturbedColumnsF90()
467c4762a1bSJed Brown!
468c4762a1bSJed Brown      if (fd_coloring) then
469*d8606c27SBarry Smith         PetscCallA(MatFDColoringGetPerturbedColumnsF90(fdcoloring,indices,ierr))
470c4762a1bSJed Brown         print*,'Indices from GetPerturbedColumnsF90'
471c4762a1bSJed Brown         write(*,1000) indices
472c4762a1bSJed Brown 1000    format(50i4)
473*d8606c27SBarry Smith         PetscCallA(MatFDColoringRestorePerturbedColumnsF90(fdcoloring,indices,ierr))
474c4762a1bSJed Brown      endif
475c4762a1bSJed Brown      return
476c4762a1bSJed Brown      end
477c4762a1bSJed Brown
478c4762a1bSJed Brown! ---------------------------------------------------------------------
479c4762a1bSJed Brown!
480c4762a1bSJed Brown!  ApplicationFunction - Computes nonlinear function, called by
481c4762a1bSJed Brown!  the higher level routine FormFunction().
482c4762a1bSJed Brown!
483c4762a1bSJed Brown!  Input Parameter:
484c4762a1bSJed Brown!  x    - local vector data
485c4762a1bSJed Brown!
486c4762a1bSJed Brown!  Output Parameters:
487c4762a1bSJed Brown!  f    - local vector data, f(x)
488c4762a1bSJed Brown!  ierr - error code
489c4762a1bSJed Brown!
490c4762a1bSJed Brown!  Notes:
491c4762a1bSJed Brown!  This routine uses standard Fortran-style computations over a 2-dim array.
492c4762a1bSJed Brown!
493c4762a1bSJed Brown      subroutine ApplicationFunction(x,f,ierr)
494c4762a1bSJed Brown      use petscsnes
495c4762a1bSJed Brown      implicit none
496c4762a1bSJed Brown
497c4762a1bSJed Brown!  Common blocks:
498c4762a1bSJed Brown      PetscReal      lambda
499c4762a1bSJed Brown      PetscInt        mx,my
500c4762a1bSJed Brown      PetscBool         fd_coloring
501c4762a1bSJed Brown      common         /params/ lambda,mx,my,fd_coloring
502c4762a1bSJed Brown
503c4762a1bSJed Brown!  Input/output variables:
504c4762a1bSJed Brown      PetscScalar    x(mx,my),f(mx,my)
505c4762a1bSJed Brown      PetscErrorCode       ierr
506c4762a1bSJed Brown
507c4762a1bSJed Brown!  Local variables:
508c4762a1bSJed Brown      PetscScalar    two,one,hx,hy
509c4762a1bSJed Brown      PetscScalar    hxdhy,hydhx,sc
510c4762a1bSJed Brown      PetscScalar    u,uxx,uyy
511c4762a1bSJed Brown      PetscInt        i,j
512c4762a1bSJed Brown
513c4762a1bSJed Brown      ierr   = 0
514c4762a1bSJed Brown      one    = 1.0
515c4762a1bSJed Brown      two    = 2.0
516c4762a1bSJed Brown      hx     = one/(mx-1)
517c4762a1bSJed Brown      hy     = one/(my-1)
518c4762a1bSJed Brown      sc     = hx*hy*lambda
519c4762a1bSJed Brown      hxdhy  = hx/hy
520c4762a1bSJed Brown      hydhx  = hy/hx
521c4762a1bSJed Brown
522c4762a1bSJed Brown!  Compute function
523c4762a1bSJed Brown
524c4762a1bSJed Brown      do 20 j=1,my
525c4762a1bSJed Brown         do 10 i=1,mx
526c4762a1bSJed Brown            if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then
527c4762a1bSJed Brown               f(i,j) = x(i,j)
528c4762a1bSJed Brown            else
529c4762a1bSJed Brown               u = x(i,j)
530c4762a1bSJed Brown               uxx = hydhx * (two*u - x(i-1,j) - x(i+1,j))
531c4762a1bSJed Brown               uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1))
532c4762a1bSJed Brown               f(i,j) = uxx + uyy - sc*exp(u)
533c4762a1bSJed Brown            endif
534c4762a1bSJed Brown 10      continue
535c4762a1bSJed Brown 20   continue
536c4762a1bSJed Brown
537c4762a1bSJed Brown      return
538c4762a1bSJed Brown      end
539c4762a1bSJed Brown
540c4762a1bSJed Brown! ---------------------------------------------------------------------
541c4762a1bSJed Brown!
542c4762a1bSJed Brown!  FormJacobian - Evaluates Jacobian matrix.
543c4762a1bSJed Brown!
544c4762a1bSJed Brown!  Input Parameters:
545c4762a1bSJed Brown!  snes    - the SNES context
546c4762a1bSJed Brown!  x       - input vector
547c4762a1bSJed Brown!  dummy   - optional user-defined context, as set by SNESSetJacobian()
548c4762a1bSJed Brown!            (not used here)
549c4762a1bSJed Brown!
550c4762a1bSJed Brown!  Output Parameters:
551c4762a1bSJed Brown!  jac      - Jacobian matrix
552c4762a1bSJed Brown!  jac_prec - optionally different preconditioning matrix (not used here)
553c4762a1bSJed Brown!  flag     - flag indicating matrix structure
554c4762a1bSJed Brown!
555c4762a1bSJed Brown!  Notes:
556c4762a1bSJed Brown!  This routine serves as a wrapper for the lower-level routine
557c4762a1bSJed Brown!  "ApplicationJacobian", where the actual computations are
558c4762a1bSJed Brown!  done using the standard Fortran style of treating the local
559c4762a1bSJed Brown!  vector data as a multidimensional array over the local mesh.
560c4762a1bSJed Brown!  This routine merely accesses the local vector data via
561c4762a1bSJed Brown!  VecGetArray() and VecRestoreArray().
562c4762a1bSJed Brown!
563c4762a1bSJed Brown      subroutine FormJacobian(snes,X,jac,jac_prec,dummy,ierr)
564c4762a1bSJed Brown      use petscsnes
565c4762a1bSJed Brown      implicit none
566c4762a1bSJed Brown
567c4762a1bSJed Brown!  Input/output variables:
568c4762a1bSJed Brown      SNES          snes
569c4762a1bSJed Brown      Vec           X
570c4762a1bSJed Brown      Mat           jac,jac_prec
571c4762a1bSJed Brown      PetscErrorCode      ierr
572c4762a1bSJed Brown      integer dummy
573c4762a1bSJed Brown
574c4762a1bSJed Brown!  Common blocks:
575c4762a1bSJed Brown      PetscReal     lambda
576c4762a1bSJed Brown      PetscInt       mx,my
577c4762a1bSJed Brown      PetscBool         fd_coloring
578c4762a1bSJed Brown      common        /params/ lambda,mx,my,fd_coloring
579c4762a1bSJed Brown
580c4762a1bSJed Brown!  Declarations for use with local array:
581c4762a1bSJed Brown      PetscScalar   lx_v(0:1)
582c4762a1bSJed Brown      PetscOffset   lx_i
583c4762a1bSJed Brown
584c4762a1bSJed Brown!  Get a pointer to vector data
585c4762a1bSJed Brown
586*d8606c27SBarry Smith      PetscCallA(VecGetArrayRead(X,lx_v,lx_i,ierr))
587c4762a1bSJed Brown
588c4762a1bSJed Brown!  Compute Jacobian entries
589c4762a1bSJed Brown
590*d8606c27SBarry Smith      PetscCallA(ApplicationJacobian(lx_v(lx_i),jac,jac_prec,ierr))
591c4762a1bSJed Brown
592c4762a1bSJed Brown!  Restore vector
593c4762a1bSJed Brown
594*d8606c27SBarry Smith      PetscCallA(VecRestoreArrayRead(X,lx_v,lx_i,ierr))
595c4762a1bSJed Brown
596c4762a1bSJed Brown!  Assemble matrix
597c4762a1bSJed Brown
598*d8606c27SBarry Smith      PetscCallA(MatAssemblyBegin(jac_prec,MAT_FINAL_ASSEMBLY,ierr))
599*d8606c27SBarry Smith      PetscCallA(MatAssemblyEnd(jac_prec,MAT_FINAL_ASSEMBLY,ierr))
600c4762a1bSJed Brown
601c4762a1bSJed Brown      return
602c4762a1bSJed Brown      end
603c4762a1bSJed Brown
604c4762a1bSJed Brown! ---------------------------------------------------------------------
605c4762a1bSJed Brown!
606c4762a1bSJed Brown!  ApplicationJacobian - Computes Jacobian matrix, called by
607c4762a1bSJed Brown!  the higher level routine FormJacobian().
608c4762a1bSJed Brown!
609c4762a1bSJed Brown!  Input Parameters:
610c4762a1bSJed Brown!  x        - local vector data
611c4762a1bSJed Brown!
612c4762a1bSJed Brown!  Output Parameters:
613c4762a1bSJed Brown!  jac      - Jacobian matrix
614c4762a1bSJed Brown!  jac_prec - optionally different preconditioning matrix (not used here)
615c4762a1bSJed Brown!  ierr     - error code
616c4762a1bSJed Brown!
617c4762a1bSJed Brown!  Notes:
618c4762a1bSJed Brown!  This routine uses standard Fortran-style computations over a 2-dim array.
619c4762a1bSJed Brown!
620c4762a1bSJed Brown      subroutine ApplicationJacobian(x,jac,jac_prec,ierr)
621c4762a1bSJed Brown      use petscsnes
622c4762a1bSJed Brown      implicit none
623c4762a1bSJed Brown
624c4762a1bSJed Brown!  Common blocks:
625c4762a1bSJed Brown      PetscReal    lambda
626c4762a1bSJed Brown      PetscInt      mx,my
627c4762a1bSJed Brown      PetscBool         fd_coloring
628c4762a1bSJed Brown      common       /params/ lambda,mx,my,fd_coloring
629c4762a1bSJed Brown
630c4762a1bSJed Brown!  Input/output variables:
631c4762a1bSJed Brown      PetscScalar  x(mx,my)
632c4762a1bSJed Brown      Mat          jac,jac_prec
633c4762a1bSJed Brown      PetscErrorCode      ierr
634c4762a1bSJed Brown
635c4762a1bSJed Brown!  Local variables:
636c4762a1bSJed Brown      PetscInt      i,j,row(1),col(5),i1,i5
637c4762a1bSJed Brown      PetscScalar  two,one, hx,hy
638c4762a1bSJed Brown      PetscScalar  hxdhy,hydhx,sc,v(5)
639c4762a1bSJed Brown
640c4762a1bSJed Brown!  Set parameters
641c4762a1bSJed Brown
642c4762a1bSJed Brown      i1     = 1
643c4762a1bSJed Brown      i5     = 5
644c4762a1bSJed Brown      one    = 1.0
645c4762a1bSJed Brown      two    = 2.0
646c4762a1bSJed Brown      hx     = one/(mx-1)
647c4762a1bSJed Brown      hy     = one/(my-1)
648c4762a1bSJed Brown      sc     = hx*hy
649c4762a1bSJed Brown      hxdhy  = hx/hy
650c4762a1bSJed Brown      hydhx  = hy/hx
651c4762a1bSJed Brown
652c4762a1bSJed Brown!  Compute entries of the Jacobian matrix
653c4762a1bSJed Brown!   - Here, we set all entries for a particular row at once.
654c4762a1bSJed Brown!   - Note that MatSetValues() uses 0-based row and column numbers
655c4762a1bSJed Brown!     in Fortran as well as in C.
656c4762a1bSJed Brown
657c4762a1bSJed Brown      do 20 j=1,my
658c4762a1bSJed Brown         row(1) = (j-1)*mx - 1
659c4762a1bSJed Brown         do 10 i=1,mx
660c4762a1bSJed Brown            row(1) = row(1) + 1
661c4762a1bSJed Brown!           boundary points
662c4762a1bSJed Brown            if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then
663*d8606c27SBarry Smith               PetscCallA(MatSetValues(jac_prec,i1,row,i1,row,one,INSERT_VALUES,ierr))
664c4762a1bSJed Brown!           interior grid points
665c4762a1bSJed Brown            else
666c4762a1bSJed Brown               v(1) = -hxdhy
667c4762a1bSJed Brown               v(2) = -hydhx
668c4762a1bSJed Brown               v(3) = two*(hydhx + hxdhy) - sc*lambda*exp(x(i,j))
669c4762a1bSJed Brown               v(4) = -hydhx
670c4762a1bSJed Brown               v(5) = -hxdhy
671c4762a1bSJed Brown               col(1) = row(1) - mx
672c4762a1bSJed Brown               col(2) = row(1) - 1
673c4762a1bSJed Brown               col(3) = row(1)
674c4762a1bSJed Brown               col(4) = row(1) + 1
675c4762a1bSJed Brown               col(5) = row(1) + mx
676*d8606c27SBarry Smith               PetscCallA(MatSetValues(jac_prec,i1,row,i5,col,v,INSERT_VALUES,ierr))
677c4762a1bSJed Brown            endif
678c4762a1bSJed Brown 10      continue
679c4762a1bSJed Brown 20   continue
680c4762a1bSJed Brown
681c4762a1bSJed Brown      return
682c4762a1bSJed Brown      end
683c4762a1bSJed Brown
684c4762a1bSJed Brown!
685c4762a1bSJed Brown!/*TEST
686c4762a1bSJed Brown!
687c4762a1bSJed Brown!   build:
688c4762a1bSJed Brown!      requires: !single
689c4762a1bSJed Brown!
690c4762a1bSJed Brown!   test:
691c4762a1bSJed Brown!      args: -snes_monitor_short -nox -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always
692c4762a1bSJed Brown!
693c4762a1bSJed Brown!   test:
694c4762a1bSJed Brown!      suffix: 2
695c4762a1bSJed Brown!      args: -snes_monitor_short -nox -snes_fd -ksp_gmres_cgs_refinement_type refine_always
696c4762a1bSJed Brown!
697c4762a1bSJed Brown!   test:
698c4762a1bSJed Brown!      suffix: 3
699c4762a1bSJed Brown!      args: -snes_monitor_short -nox -snes_fd_coloring -mat_coloring_type sl -ksp_gmres_cgs_refinement_type refine_always
700c4762a1bSJed Brown!      filter: sort -b
701c4762a1bSJed Brown!      filter_output: sort -b
702c4762a1bSJed Brown!
703c4762a1bSJed Brown!   test:
704c4762a1bSJed Brown!     suffix: 4
705c4762a1bSJed Brown!     args: -pc -par 6.807 -nox
706c4762a1bSJed Brown!
707c4762a1bSJed Brown!TEST*/
708