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