xref: /petsc/src/snes/tests/ex1.c (revision ffc4695bcb29f4b022f59a5fd6bc99fc280ff6d8)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Solves the nonlinear system, the Bratu (SFI - solid fuel ignition) problem in a 2D rectangular domain.\n\
3c4762a1bSJed Brown This example also illustrates the use of matrix coloring.  Runtime options include:\n\
4c4762a1bSJed Brown   -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\
5c4762a1bSJed Brown      problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)\n\
6c4762a1bSJed Brown   -mx <xg>, where <xg> = number of grid points in the x-direction\n\
7c4762a1bSJed Brown   -my <yg>, where <yg> = number of grid points in the y-direction\n\n";
8c4762a1bSJed Brown 
9c4762a1bSJed Brown /*T
10c4762a1bSJed Brown    Concepts: SNES^sequential Bratu example
11c4762a1bSJed Brown    Processors: 1
12c4762a1bSJed Brown T*/
13c4762a1bSJed Brown 
14c4762a1bSJed Brown 
15c4762a1bSJed Brown 
16c4762a1bSJed Brown /* ------------------------------------------------------------------------
17c4762a1bSJed Brown 
18c4762a1bSJed Brown     Solid Fuel Ignition (SFI) problem.  This problem is modeled by
19c4762a1bSJed Brown     the partial differential equation
20c4762a1bSJed Brown 
21c4762a1bSJed Brown             -Laplacian u - lambda*exp(u) = 0,  0 < x,y < 1,
22c4762a1bSJed Brown 
23c4762a1bSJed Brown     with boundary conditions
24c4762a1bSJed Brown 
25c4762a1bSJed Brown              u = 0  for  x = 0, x = 1, y = 0, y = 1.
26c4762a1bSJed Brown 
27c4762a1bSJed Brown     A finite difference approximation with the usual 5-point stencil
28c4762a1bSJed Brown     is used to discretize the boundary value problem to obtain a nonlinear
29c4762a1bSJed Brown     system of equations.
30c4762a1bSJed Brown 
31c4762a1bSJed Brown     The parallel version of this code is snes/tutorials/ex5.c
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   ------------------------------------------------------------------------- */
34c4762a1bSJed Brown 
35c4762a1bSJed Brown /*
36c4762a1bSJed Brown    Include "petscsnes.h" so that we can use SNES solvers.  Note that
37c4762a1bSJed Brown    this file automatically includes:
38c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
39c4762a1bSJed Brown      petscmat.h - matrices
40c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
41c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
42c4762a1bSJed Brown      petscksp.h   - linear solvers
43c4762a1bSJed Brown */
44c4762a1bSJed Brown 
45c4762a1bSJed Brown #include <petscsnes.h>
46c4762a1bSJed Brown 
47c4762a1bSJed Brown /*
48c4762a1bSJed Brown    User-defined application context - contains data needed by the
49c4762a1bSJed Brown    application-provided call-back routines, FormJacobian() and
50c4762a1bSJed Brown    FormFunction().
51c4762a1bSJed Brown */
52c4762a1bSJed Brown typedef struct {
53c4762a1bSJed Brown   PetscReal param;              /* test problem parameter */
54c4762a1bSJed Brown   PetscInt  mx;                 /* Discretization in x-direction */
55c4762a1bSJed Brown   PetscInt  my;                 /* Discretization in y-direction */
56c4762a1bSJed Brown } AppCtx;
57c4762a1bSJed Brown 
58c4762a1bSJed Brown /*
59c4762a1bSJed Brown    User-defined routines
60c4762a1bSJed Brown */
61c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);
62c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
63c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(AppCtx*,Vec);
64c4762a1bSJed Brown extern PetscErrorCode ConvergenceTest(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
65c4762a1bSJed Brown extern PetscErrorCode ConvergenceDestroy(void*);
66c4762a1bSJed Brown extern PetscErrorCode postcheck(SNES,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
67c4762a1bSJed Brown 
68c4762a1bSJed Brown int main(int argc,char **argv)
69c4762a1bSJed Brown {
70c4762a1bSJed Brown   SNES           snes;                 /* nonlinear solver context */
71c4762a1bSJed Brown   Vec            x,r;                 /* solution, residual vectors */
72c4762a1bSJed Brown   Mat            J;                    /* Jacobian matrix */
73c4762a1bSJed Brown   AppCtx         user;                 /* user-defined application context */
74c4762a1bSJed Brown   PetscErrorCode ierr;
75c4762a1bSJed Brown   PetscInt       i,its,N,hist_its[50];
76c4762a1bSJed Brown   PetscMPIInt    size;
77c4762a1bSJed Brown   PetscReal      bratu_lambda_max = 6.81,bratu_lambda_min = 0.,history[50];
78c4762a1bSJed Brown   MatFDColoring  fdcoloring;
79c4762a1bSJed Brown   PetscBool      matrix_free = PETSC_FALSE,flg,fd_coloring = PETSC_FALSE, use_convergence_test = PETSC_FALSE,pc = PETSC_FALSE;
80c4762a1bSJed Brown   KSP            ksp;
81c4762a1bSJed Brown   PetscInt       *testarray;
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
84*ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
85c4762a1bSJed Brown   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   /*
88c4762a1bSJed Brown      Initialize problem parameters
89c4762a1bSJed Brown   */
90c4762a1bSJed Brown   user.mx = 4; user.my = 4; user.param = 6.0;
91c4762a1bSJed Brown   ierr    = PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,NULL);CHKERRQ(ierr);
92c4762a1bSJed Brown   ierr    = PetscOptionsGetInt(NULL,NULL,"-my",&user.my,NULL);CHKERRQ(ierr);
93c4762a1bSJed Brown   ierr    = PetscOptionsGetReal(NULL,NULL,"-par",&user.param,NULL);CHKERRQ(ierr);
94c4762a1bSJed Brown   ierr    = PetscOptionsGetBool(NULL,NULL,"-pc",&pc,NULL);CHKERRQ(ierr);
95c4762a1bSJed Brown   if (user.param >= bratu_lambda_max || user.param <= bratu_lambda_min) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lambda is out of range");
96c4762a1bSJed Brown   N = user.mx*user.my;
97c4762a1bSJed Brown   ierr    = PetscOptionsGetBool(NULL,NULL,"-use_convergence_test",&use_convergence_test,NULL);CHKERRQ(ierr);
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100c4762a1bSJed Brown      Create nonlinear solver context
101c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr);
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   if (pc) {
106c4762a1bSJed Brown     ierr = SNESSetType(snes,SNESNEWTONTR);CHKERRQ(ierr);
107c4762a1bSJed Brown     ierr = SNESNewtonTRSetPostCheck(snes, postcheck,NULL);CHKERRQ(ierr);
108c4762a1bSJed Brown   }
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111c4762a1bSJed Brown      Create vector data structures; set function evaluation routine
112c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr);
115c4762a1bSJed Brown   ierr = VecSetSizes(x,PETSC_DECIDE,N);CHKERRQ(ierr);
116c4762a1bSJed Brown   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
117c4762a1bSJed Brown   ierr = VecDuplicate(x,&r);CHKERRQ(ierr);
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   /*
120c4762a1bSJed Brown      Set function evaluation routine and vector.  Whenever the nonlinear
121c4762a1bSJed Brown      solver needs to evaluate the nonlinear function, it will call this
122c4762a1bSJed Brown      routine.
123c4762a1bSJed Brown       - Note that the final routine argument is the user-defined
124c4762a1bSJed Brown         context that provides application-specific data for the
125c4762a1bSJed Brown         function evaluation routine.
126c4762a1bSJed Brown   */
127c4762a1bSJed Brown   ierr = SNESSetFunction(snes,r,FormFunction,(void*)&user);CHKERRQ(ierr);
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine
131c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   /*
134c4762a1bSJed Brown      Create matrix. Here we only approximately preallocate storage space
135c4762a1bSJed Brown      for the Jacobian.  See the users manual for a discussion of better
136c4762a1bSJed Brown      techniques for preallocating matrix memory.
137c4762a1bSJed Brown   */
138c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-snes_mf",&matrix_free,NULL);CHKERRQ(ierr);
139c4762a1bSJed Brown   if (!matrix_free) {
140c4762a1bSJed Brown     PetscBool matrix_free_operator = PETSC_FALSE;
141c4762a1bSJed Brown     ierr = PetscOptionsGetBool(NULL,NULL,"-snes_mf_operator",&matrix_free_operator,NULL);CHKERRQ(ierr);
142c4762a1bSJed Brown     if (matrix_free_operator) matrix_free = PETSC_FALSE;
143c4762a1bSJed Brown   }
144c4762a1bSJed Brown   if (!matrix_free) {
145c4762a1bSJed Brown     ierr = MatCreateSeqAIJ(PETSC_COMM_WORLD,N,N,5,NULL,&J);CHKERRQ(ierr);
146c4762a1bSJed Brown   }
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   /*
149c4762a1bSJed Brown      This option will cause the Jacobian to be computed via finite differences
150c4762a1bSJed Brown     efficiently using a coloring of the columns of the matrix.
151c4762a1bSJed Brown   */
152c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-snes_fd_coloring",&fd_coloring,NULL);CHKERRQ(ierr);
153c4762a1bSJed Brown   if (matrix_free && fd_coloring) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"Use only one of -snes_mf, -snes_fd_coloring options!\nYou can do -snes_mf_operator -snes_fd_coloring");
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   if (fd_coloring) {
156c4762a1bSJed Brown     ISColoring   iscoloring;
157c4762a1bSJed Brown     MatColoring  mc;
158c4762a1bSJed Brown 
159c4762a1bSJed Brown     /*
160c4762a1bSJed Brown       This initializes the nonzero structure of the Jacobian. This is artificial
161c4762a1bSJed Brown       because clearly if we had a routine to compute the Jacobian we won't need
162c4762a1bSJed Brown       to use finite differences.
163c4762a1bSJed Brown     */
164c4762a1bSJed Brown     ierr = FormJacobian(snes,x,J,J,&user);CHKERRQ(ierr);
165c4762a1bSJed Brown 
166c4762a1bSJed Brown     /*
167c4762a1bSJed Brown        Color the matrix, i.e. determine groups of columns that share no common
168c4762a1bSJed Brown       rows. These columns in the Jacobian can all be computed simulataneously.
169c4762a1bSJed Brown     */
170c4762a1bSJed Brown     ierr = MatColoringCreate(J,&mc);CHKERRQ(ierr);
171c4762a1bSJed Brown     ierr = MatColoringSetType(mc,MATCOLORINGSL);CHKERRQ(ierr);
172c4762a1bSJed Brown     ierr = MatColoringSetFromOptions(mc);CHKERRQ(ierr);
173c4762a1bSJed Brown     ierr = MatColoringApply(mc,&iscoloring);CHKERRQ(ierr);
174c4762a1bSJed Brown     ierr = MatColoringDestroy(&mc);CHKERRQ(ierr);
175c4762a1bSJed Brown     /*
176c4762a1bSJed Brown        Create the data structure that SNESComputeJacobianDefaultColor() uses
177c4762a1bSJed Brown        to compute the actual Jacobians via finite differences.
178c4762a1bSJed Brown     */
179c4762a1bSJed Brown     ierr = MatFDColoringCreate(J,iscoloring,&fdcoloring);CHKERRQ(ierr);
180c4762a1bSJed Brown     ierr = MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))FormFunction,&user);CHKERRQ(ierr);
181c4762a1bSJed Brown     ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr);
182c4762a1bSJed Brown     ierr = MatFDColoringSetUp(J,iscoloring,fdcoloring);CHKERRQ(ierr);
183c4762a1bSJed Brown     /*
184c4762a1bSJed Brown         Tell SNES to use the routine SNESComputeJacobianDefaultColor()
185c4762a1bSJed Brown       to compute Jacobians.
186c4762a1bSJed Brown     */
187c4762a1bSJed Brown     ierr = SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,fdcoloring);CHKERRQ(ierr);
188c4762a1bSJed Brown     ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
189c4762a1bSJed Brown   }
190c4762a1bSJed Brown   /*
191c4762a1bSJed Brown      Set Jacobian matrix data structure and default Jacobian evaluation
192c4762a1bSJed Brown      routine.  Whenever the nonlinear solver needs to compute the
193c4762a1bSJed Brown      Jacobian matrix, it will call this routine.
194c4762a1bSJed Brown       - Note that the final routine argument is the user-defined
195c4762a1bSJed Brown         context that provides application-specific data for the
196c4762a1bSJed Brown         Jacobian evaluation routine.
197c4762a1bSJed Brown       - The user can override with:
198c4762a1bSJed Brown          -snes_fd : default finite differencing approximation of Jacobian
199c4762a1bSJed Brown          -snes_mf : matrix-free Newton-Krylov method with no preconditioning
200c4762a1bSJed Brown                     (unless user explicitly sets preconditioner)
201c4762a1bSJed Brown          -snes_mf_operator : form preconditioning matrix as set by the user,
202c4762a1bSJed Brown                              but use matrix-free approx for Jacobian-vector
203c4762a1bSJed Brown                              products within Newton-Krylov method
204c4762a1bSJed Brown   */
205c4762a1bSJed Brown   else if (!matrix_free) {
206c4762a1bSJed Brown     ierr = SNESSetJacobian(snes,J,J,FormJacobian,(void*)&user);CHKERRQ(ierr);
207c4762a1bSJed Brown   }
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
210c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
211c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
212c4762a1bSJed Brown 
213c4762a1bSJed Brown   /*
214c4762a1bSJed Brown      Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
215c4762a1bSJed Brown   */
216c4762a1bSJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   /*
219c4762a1bSJed Brown      Set array that saves the function norms.  This array is intended
220c4762a1bSJed Brown      when the user wants to save the convergence history for later use
221c4762a1bSJed Brown      rather than just to view the function norms via -snes_monitor.
222c4762a1bSJed Brown   */
223c4762a1bSJed Brown   ierr = SNESSetConvergenceHistory(snes,history,hist_its,50,PETSC_TRUE);CHKERRQ(ierr);
224c4762a1bSJed Brown 
225c4762a1bSJed Brown   /*
226c4762a1bSJed Brown       Add a user provided convergence test; this is to test that SNESNEWTONTR properly calls the
227c4762a1bSJed Brown       user provided test before the specialized test. The convergence context is just an array to
228c4762a1bSJed Brown       test that it gets properly freed at the end
229c4762a1bSJed Brown   */
230c4762a1bSJed Brown   if (use_convergence_test) {
231c4762a1bSJed Brown     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
232c4762a1bSJed Brown     ierr = PetscMalloc1(5,&testarray);CHKERRQ(ierr);
233c4762a1bSJed Brown     ierr = KSPSetConvergenceTest(ksp,ConvergenceTest,testarray,ConvergenceDestroy);CHKERRQ(ierr);
234c4762a1bSJed Brown   }
235c4762a1bSJed Brown 
236c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
237c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
238c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
239c4762a1bSJed Brown   /*
240c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
241c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
242c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
243c4762a1bSJed Brown      this vector to zero by calling VecSet().
244c4762a1bSJed Brown   */
245c4762a1bSJed Brown   ierr = FormInitialGuess(&user,x);CHKERRQ(ierr);
246c4762a1bSJed Brown   ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr);
247c4762a1bSJed Brown   ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
248c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D\n",its);CHKERRQ(ierr);
249c4762a1bSJed Brown 
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   /*
252c4762a1bSJed Brown      Print the convergence history.  This is intended just to demonstrate
253c4762a1bSJed Brown      use of the data attained via SNESSetConvergenceHistory().
254c4762a1bSJed Brown   */
255c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-print_history",&flg);CHKERRQ(ierr);
256c4762a1bSJed Brown   if (flg) {
257c4762a1bSJed Brown     for (i=0; i<its+1; i++) {
258c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_WORLD,"iteration %D: Linear iterations %D Function norm = %g\n",i,hist_its[i],(double)history[i]);CHKERRQ(ierr);
259c4762a1bSJed Brown     }
260c4762a1bSJed Brown   }
261c4762a1bSJed Brown 
262c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
263c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
264c4762a1bSJed Brown      are no longer needed.
265c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
266c4762a1bSJed Brown 
267c4762a1bSJed Brown   if (!matrix_free) {
268c4762a1bSJed Brown     ierr = MatDestroy(&J);CHKERRQ(ierr);
269c4762a1bSJed Brown   }
270c4762a1bSJed Brown   if (fd_coloring) {
271c4762a1bSJed Brown     ierr = MatFDColoringDestroy(&fdcoloring);CHKERRQ(ierr);
272c4762a1bSJed Brown   }
273c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
274c4762a1bSJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
275c4762a1bSJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
276c4762a1bSJed Brown   ierr = PetscFinalize();
277c4762a1bSJed Brown   return ierr;
278c4762a1bSJed Brown }
279c4762a1bSJed Brown /* ------------------------------------------------------------------- */
280c4762a1bSJed Brown /*
281c4762a1bSJed Brown    FormInitialGuess - Forms initial approximation.
282c4762a1bSJed Brown 
283c4762a1bSJed Brown    Input Parameters:
284c4762a1bSJed Brown    user - user-defined application context
285c4762a1bSJed Brown    X - vector
286c4762a1bSJed Brown 
287c4762a1bSJed Brown    Output Parameter:
288c4762a1bSJed Brown    X - vector
289c4762a1bSJed Brown  */
290c4762a1bSJed Brown PetscErrorCode FormInitialGuess(AppCtx *user,Vec X)
291c4762a1bSJed Brown {
292c4762a1bSJed Brown   PetscInt       i,j,row,mx,my;
293c4762a1bSJed Brown   PetscErrorCode ierr;
294c4762a1bSJed Brown   PetscReal      lambda,temp1,temp,hx,hy;
295c4762a1bSJed Brown   PetscScalar    *x;
296c4762a1bSJed Brown 
297c4762a1bSJed Brown   mx     = user->mx;
298c4762a1bSJed Brown   my     = user->my;
299c4762a1bSJed Brown   lambda = user->param;
300c4762a1bSJed Brown 
301c4762a1bSJed Brown   hx = 1.0 / (PetscReal)(mx-1);
302c4762a1bSJed Brown   hy = 1.0 / (PetscReal)(my-1);
303c4762a1bSJed Brown 
304c4762a1bSJed Brown   /*
305c4762a1bSJed Brown      Get a pointer to vector data.
306c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
307c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
308c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
309c4762a1bSJed Brown          the array.
310c4762a1bSJed Brown   */
311c4762a1bSJed Brown   ierr  = VecGetArray(X,&x);CHKERRQ(ierr);
312c4762a1bSJed Brown   temp1 = lambda/(lambda + 1.0);
313c4762a1bSJed Brown   for (j=0; j<my; j++) {
314c4762a1bSJed Brown     temp = (PetscReal)(PetscMin(j,my-j-1))*hy;
315c4762a1bSJed Brown     for (i=0; i<mx; i++) {
316c4762a1bSJed Brown       row = i + j*mx;
317c4762a1bSJed Brown       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
318c4762a1bSJed Brown         x[row] = 0.0;
319c4762a1bSJed Brown         continue;
320c4762a1bSJed Brown       }
321c4762a1bSJed Brown       x[row] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,mx-i-1))*hx,temp));
322c4762a1bSJed Brown     }
323c4762a1bSJed Brown   }
324c4762a1bSJed Brown 
325c4762a1bSJed Brown   /*
326c4762a1bSJed Brown      Restore vector
327c4762a1bSJed Brown   */
328c4762a1bSJed Brown   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
329c4762a1bSJed Brown   return 0;
330c4762a1bSJed Brown }
331c4762a1bSJed Brown /* ------------------------------------------------------------------- */
332c4762a1bSJed Brown /*
333c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
334c4762a1bSJed Brown 
335c4762a1bSJed Brown    Input Parameters:
336c4762a1bSJed Brown .  snes - the SNES context
337c4762a1bSJed Brown .  X - input vector
338c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetFunction()
339c4762a1bSJed Brown 
340c4762a1bSJed Brown    Output Parameter:
341c4762a1bSJed Brown .  F - function vector
342c4762a1bSJed Brown  */
343c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ptr)
344c4762a1bSJed Brown {
345c4762a1bSJed Brown   AppCtx            *user = (AppCtx*)ptr;
346c4762a1bSJed Brown   PetscInt          i,j,row,mx,my;
347c4762a1bSJed Brown   PetscErrorCode    ierr;
348c4762a1bSJed Brown   PetscReal         two = 2.0,one = 1.0,lambda,hx,hy,hxdhy,hydhx;
349c4762a1bSJed Brown   PetscScalar       ut,ub,ul,ur,u,uxx,uyy,sc,*f;
350c4762a1bSJed Brown   const PetscScalar *x;
351c4762a1bSJed Brown 
352c4762a1bSJed Brown   mx     = user->mx;
353c4762a1bSJed Brown   my     = user->my;
354c4762a1bSJed Brown   lambda = user->param;
355c4762a1bSJed Brown   hx     = one / (PetscReal)(mx-1);
356c4762a1bSJed Brown   hy     = one / (PetscReal)(my-1);
357c4762a1bSJed Brown   sc     = hx*hy;
358c4762a1bSJed Brown   hxdhy  = hx/hy;
359c4762a1bSJed Brown   hydhx  = hy/hx;
360c4762a1bSJed Brown 
361c4762a1bSJed Brown   /*
362c4762a1bSJed Brown      Get pointers to vector data
363c4762a1bSJed Brown   */
364c4762a1bSJed Brown   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
365c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
366c4762a1bSJed Brown 
367c4762a1bSJed Brown   /*
368c4762a1bSJed Brown      Compute function
369c4762a1bSJed Brown   */
370c4762a1bSJed Brown   for (j=0; j<my; j++) {
371c4762a1bSJed Brown     for (i=0; i<mx; i++) {
372c4762a1bSJed Brown       row = i + j*mx;
373c4762a1bSJed Brown       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
374c4762a1bSJed Brown         f[row] = x[row];
375c4762a1bSJed Brown         continue;
376c4762a1bSJed Brown       }
377c4762a1bSJed Brown       u      = x[row];
378c4762a1bSJed Brown       ub     = x[row - mx];
379c4762a1bSJed Brown       ul     = x[row - 1];
380c4762a1bSJed Brown       ut     = x[row + mx];
381c4762a1bSJed Brown       ur     = x[row + 1];
382c4762a1bSJed Brown       uxx    = (-ur + two*u - ul)*hydhx;
383c4762a1bSJed Brown       uyy    = (-ut + two*u - ub)*hxdhy;
384c4762a1bSJed Brown       f[row] = uxx + uyy - sc*lambda*PetscExpScalar(u);
385c4762a1bSJed Brown     }
386c4762a1bSJed Brown   }
387c4762a1bSJed Brown 
388c4762a1bSJed Brown   /*
389c4762a1bSJed Brown      Restore vectors
390c4762a1bSJed Brown   */
391c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
392c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
393c4762a1bSJed Brown   return 0;
394c4762a1bSJed Brown }
395c4762a1bSJed Brown /* ------------------------------------------------------------------- */
396c4762a1bSJed Brown /*
397c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
398c4762a1bSJed Brown 
399c4762a1bSJed Brown    Input Parameters:
400c4762a1bSJed Brown .  snes - the SNES context
401c4762a1bSJed Brown .  x - input vector
402c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetJacobian()
403c4762a1bSJed Brown 
404c4762a1bSJed Brown    Output Parameters:
405c4762a1bSJed Brown .  A - Jacobian matrix
406c4762a1bSJed Brown .  B - optionally different preconditioning matrix
407c4762a1bSJed Brown .  flag - flag indicating matrix structure
408c4762a1bSJed Brown */
409c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec X,Mat J,Mat jac,void *ptr)
410c4762a1bSJed Brown {
411c4762a1bSJed Brown   AppCtx            *user = (AppCtx*)ptr;   /* user-defined applicatin context */
412c4762a1bSJed Brown   PetscInt          i,j,row,mx,my,col[5];
413c4762a1bSJed Brown   PetscErrorCode    ierr;
414c4762a1bSJed Brown   PetscScalar       two = 2.0,one = 1.0,lambda,v[5],sc;
415c4762a1bSJed Brown   const PetscScalar *x;
416c4762a1bSJed Brown   PetscReal         hx,hy,hxdhy,hydhx;
417c4762a1bSJed Brown 
418c4762a1bSJed Brown   mx     = user->mx;
419c4762a1bSJed Brown   my     = user->my;
420c4762a1bSJed Brown   lambda = user->param;
421c4762a1bSJed Brown   hx     = 1.0 / (PetscReal)(mx-1);
422c4762a1bSJed Brown   hy     = 1.0 / (PetscReal)(my-1);
423c4762a1bSJed Brown   sc     = hx*hy;
424c4762a1bSJed Brown   hxdhy  = hx/hy;
425c4762a1bSJed Brown   hydhx  = hy/hx;
426c4762a1bSJed Brown 
427c4762a1bSJed Brown   /*
428c4762a1bSJed Brown      Get pointer to vector data
429c4762a1bSJed Brown   */
430c4762a1bSJed Brown   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
431c4762a1bSJed Brown 
432c4762a1bSJed Brown   /*
433c4762a1bSJed Brown      Compute entries of the Jacobian
434c4762a1bSJed Brown   */
435c4762a1bSJed Brown   for (j=0; j<my; j++) {
436c4762a1bSJed Brown     for (i=0; i<mx; i++) {
437c4762a1bSJed Brown       row = i + j*mx;
438c4762a1bSJed Brown       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
439c4762a1bSJed Brown         ierr = MatSetValues(jac,1,&row,1,&row,&one,INSERT_VALUES);CHKERRQ(ierr);
440c4762a1bSJed Brown         continue;
441c4762a1bSJed Brown       }
442c4762a1bSJed Brown       v[0] = -hxdhy; col[0] = row - mx;
443c4762a1bSJed Brown       v[1] = -hydhx; col[1] = row - 1;
444c4762a1bSJed Brown       v[2] = two*(hydhx + hxdhy) - sc*lambda*PetscExpScalar(x[row]); col[2] = row;
445c4762a1bSJed Brown       v[3] = -hydhx; col[3] = row + 1;
446c4762a1bSJed Brown       v[4] = -hxdhy; col[4] = row + mx;
447c4762a1bSJed Brown       ierr = MatSetValues(jac,1,&row,5,col,v,INSERT_VALUES);CHKERRQ(ierr);
448c4762a1bSJed Brown     }
449c4762a1bSJed Brown   }
450c4762a1bSJed Brown 
451c4762a1bSJed Brown   /*
452c4762a1bSJed Brown      Restore vector
453c4762a1bSJed Brown   */
454c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
455c4762a1bSJed Brown 
456c4762a1bSJed Brown   /*
457c4762a1bSJed Brown      Assemble matrix
458c4762a1bSJed Brown   */
459c4762a1bSJed Brown   ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
460c4762a1bSJed Brown   ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
461c4762a1bSJed Brown 
462c4762a1bSJed Brown   if (jac != J) {
463c4762a1bSJed Brown     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
464c4762a1bSJed Brown     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
465c4762a1bSJed Brown   }
466c4762a1bSJed Brown 
467c4762a1bSJed Brown   return 0;
468c4762a1bSJed Brown }
469c4762a1bSJed Brown 
470c4762a1bSJed Brown PetscErrorCode ConvergenceTest(KSP ksp,PetscInt it,PetscReal nrm,KSPConvergedReason *reason,void *ctx)
471c4762a1bSJed Brown {
472c4762a1bSJed Brown   PetscErrorCode ierr;
473c4762a1bSJed Brown 
474c4762a1bSJed Brown   PetscFunctionBegin;
475c4762a1bSJed Brown   *reason = KSP_CONVERGED_ITERATING;
476c4762a1bSJed Brown   if (it > 1) {
477c4762a1bSJed Brown     *reason = KSP_CONVERGED_ITS;
478c4762a1bSJed Brown     ierr = PetscInfo(NULL,"User provided convergence test returning after 2 iterations\n");CHKERRQ(ierr);
479c4762a1bSJed Brown   }
480c4762a1bSJed Brown   PetscFunctionReturn(0);
481c4762a1bSJed Brown }
482c4762a1bSJed Brown 
483c4762a1bSJed Brown PetscErrorCode ConvergenceDestroy(void* ctx)
484c4762a1bSJed Brown {
485c4762a1bSJed Brown   PetscErrorCode ierr;
486c4762a1bSJed Brown 
487c4762a1bSJed Brown   PetscFunctionBegin;
488c4762a1bSJed Brown   ierr = PetscInfo(NULL,"User provided convergence destroy called\n");CHKERRQ(ierr);
489c4762a1bSJed Brown   ierr = PetscFree(ctx);CHKERRQ(ierr);
490c4762a1bSJed Brown   PetscFunctionReturn(0);
491c4762a1bSJed Brown }
492c4762a1bSJed Brown 
493c4762a1bSJed Brown PetscErrorCode postcheck(SNES snes,Vec x,Vec y,Vec w,PetscBool *changed_y,PetscBool *changed_w,void *ctx)
494c4762a1bSJed Brown {
495c4762a1bSJed Brown   PetscErrorCode ierr;
496c4762a1bSJed Brown   PetscReal      norm;
497c4762a1bSJed Brown   Vec            tmp;
498c4762a1bSJed Brown 
499c4762a1bSJed Brown   PetscFunctionBegin;
500c4762a1bSJed Brown   ierr = VecDuplicate(x,&tmp);CHKERRQ(ierr);
501c4762a1bSJed Brown   ierr = VecWAXPY(tmp,-1.0,x,w);CHKERRQ(ierr);
502c4762a1bSJed Brown   ierr = VecNorm(tmp,NORM_2,&norm);CHKERRQ(ierr);
503c4762a1bSJed Brown   ierr = VecDestroy(&tmp);CHKERRQ(ierr);
504c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Norm of search step %g\n",(double)norm);CHKERRQ(ierr);
505c4762a1bSJed Brown   PetscFunctionReturn(0);
506c4762a1bSJed Brown }
507c4762a1bSJed Brown 
508c4762a1bSJed Brown 
509c4762a1bSJed Brown /*TEST
510c4762a1bSJed Brown 
511c4762a1bSJed Brown    build:
512c4762a1bSJed Brown       requires: !single
513c4762a1bSJed Brown 
514c4762a1bSJed Brown    test:
515c4762a1bSJed Brown       args: -ksp_gmres_cgs_refinement_type refine_always
516c4762a1bSJed Brown 
517c4762a1bSJed Brown    test:
518c4762a1bSJed Brown       suffix: 2
519c4762a1bSJed Brown       args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always
520c4762a1bSJed Brown 
521c4762a1bSJed Brown    test:
522c4762a1bSJed Brown       suffix: 2a
523c4762a1bSJed Brown       filter: grep -i KSPConvergedDefault > /dev/null && echo "Found KSPConvergedDefault"
524c4762a1bSJed Brown       args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always -info
525956f8c0dSBarry Smith       requires: define(PETSC_USE_INFO)
526c4762a1bSJed Brown 
527c4762a1bSJed Brown    test:
528c4762a1bSJed Brown       suffix: 2b
529c4762a1bSJed Brown       filter: grep -i  "User provided convergence test" > /dev/null  && echo "Found User provided convergence test"
530c4762a1bSJed Brown       args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always -use_convergence_test -info
531956f8c0dSBarry Smith       requires: define(PETSC_USE_INFO)
532c4762a1bSJed Brown 
533c4762a1bSJed Brown    test:
534c4762a1bSJed Brown       suffix: 3
535c4762a1bSJed Brown       args: -snes_monitor_short -mat_coloring_type sl -snes_fd_coloring -mx 8 -my 11 -ksp_gmres_cgs_refinement_type refine_always
536c4762a1bSJed Brown 
537c4762a1bSJed Brown    test:
538c4762a1bSJed Brown       suffix: 4
539c4762a1bSJed Brown       args: -pc -par 6.807 -snes_monitor -snes_converged_reason
540c4762a1bSJed Brown 
541c4762a1bSJed Brown TEST*/
542c4762a1bSJed Brown 
543c4762a1bSJed Brown 
544