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