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