xref: /petsc/src/snes/tests/ex1.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 /* ------------------------------------------------------------------------
10c4762a1bSJed Brown 
11c4762a1bSJed Brown     Solid Fuel Ignition (SFI) problem.  This problem is modeled by
12c4762a1bSJed Brown     the partial differential equation
13c4762a1bSJed Brown 
14c4762a1bSJed Brown             -Laplacian u - lambda*exp(u) = 0,  0 < x,y < 1,
15c4762a1bSJed Brown 
16c4762a1bSJed Brown     with boundary conditions
17c4762a1bSJed Brown 
18c4762a1bSJed Brown              u = 0  for  x = 0, x = 1, y = 0, y = 1.
19c4762a1bSJed Brown 
20c4762a1bSJed Brown     A finite difference approximation with the usual 5-point stencil
21c4762a1bSJed Brown     is used to discretize the boundary value problem to obtain a nonlinear
22c4762a1bSJed Brown     system of equations.
23c4762a1bSJed Brown 
24c4762a1bSJed Brown     The parallel version of this code is snes/tutorials/ex5.c
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   ------------------------------------------------------------------------- */
27c4762a1bSJed Brown 
28c4762a1bSJed Brown /*
29c4762a1bSJed Brown    Include "petscsnes.h" so that we can use SNES solvers.  Note that
30c4762a1bSJed Brown    this file automatically includes:
31c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
32c4762a1bSJed Brown      petscmat.h - matrices
33c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
34c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
35c4762a1bSJed Brown      petscksp.h   - linear solvers
36c4762a1bSJed Brown */
37c4762a1bSJed Brown 
38c4762a1bSJed Brown #include <petscsnes.h>
39c4762a1bSJed Brown 
40c4762a1bSJed Brown /*
41c4762a1bSJed Brown    User-defined application context - contains data needed by the
42c4762a1bSJed Brown    application-provided call-back routines, FormJacobian() and
43c4762a1bSJed Brown    FormFunction().
44c4762a1bSJed Brown */
45c4762a1bSJed Brown typedef struct {
46c4762a1bSJed Brown   PetscReal param; /* test problem parameter */
47c4762a1bSJed Brown   PetscInt  mx;    /* Discretization in x-direction */
48c4762a1bSJed Brown   PetscInt  my;    /* Discretization in y-direction */
49c4762a1bSJed Brown } AppCtx;
50c4762a1bSJed Brown 
51c4762a1bSJed Brown /*
52c4762a1bSJed Brown    User-defined routines
53c4762a1bSJed Brown */
54c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *);
55c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES, Vec, Vec, void *);
56c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(AppCtx *, Vec);
57c4762a1bSJed Brown extern PetscErrorCode ConvergenceTest(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
58c4762a1bSJed Brown extern PetscErrorCode ConvergenceDestroy(void *);
59c4762a1bSJed Brown extern PetscErrorCode postcheck(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *);
60c4762a1bSJed Brown 
619371c9d4SSatish Balay int main(int argc, char **argv) {
62c4762a1bSJed Brown   SNES          snes; /* nonlinear solver context */
63c4762a1bSJed Brown   Vec           x, r; /* solution, residual vectors */
64c4762a1bSJed Brown   Mat           J;    /* Jacobian matrix */
65c4762a1bSJed Brown   AppCtx        user; /* user-defined application context */
66c4762a1bSJed Brown   PetscInt      i, its, N, hist_its[50];
67c4762a1bSJed Brown   PetscMPIInt   size;
68c4762a1bSJed Brown   PetscReal     bratu_lambda_max = 6.81, bratu_lambda_min = 0., history[50];
69c4762a1bSJed Brown   MatFDColoring fdcoloring;
70c4762a1bSJed Brown   PetscBool     matrix_free = PETSC_FALSE, flg, fd_coloring = PETSC_FALSE, use_convergence_test = PETSC_FALSE, pc = PETSC_FALSE;
71c4762a1bSJed Brown   KSP           ksp;
72c4762a1bSJed Brown   PetscInt     *testarray;
73c4762a1bSJed Brown 
74327415f7SBarry Smith   PetscFunctionBeginUser;
759566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
769566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
77be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   /*
80c4762a1bSJed Brown      Initialize problem parameters
81c4762a1bSJed Brown   */
829371c9d4SSatish Balay   user.mx    = 4;
839371c9d4SSatish Balay   user.my    = 4;
849371c9d4SSatish Balay   user.param = 6.0;
859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &user.mx, NULL));
869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &user.my, NULL));
879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-par", &user.param, NULL));
889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-pc", &pc, NULL));
89e00437b9SBarry Smith   PetscCheck(user.param < bratu_lambda_max && user.param > bratu_lambda_min, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Lambda is out of range");
90c4762a1bSJed Brown   N = user.mx * user.my;
919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_convergence_test", &use_convergence_test, NULL));
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94c4762a1bSJed Brown      Create nonlinear solver context
95c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96c4762a1bSJed Brown 
979566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   if (pc) {
1009566063dSJacob Faibussowitsch     PetscCall(SNESSetType(snes, SNESNEWTONTR));
1019566063dSJacob Faibussowitsch     PetscCall(SNESNewtonTRSetPostCheck(snes, postcheck, NULL));
102c4762a1bSJed Brown   }
103c4762a1bSJed Brown 
104c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105c4762a1bSJed Brown      Create vector data structures; set function evaluation routine
106c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
107c4762a1bSJed Brown 
1089566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
1099566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, PETSC_DECIDE, N));
1109566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
1119566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   /*
114c4762a1bSJed Brown      Set function evaluation routine and vector.  Whenever the nonlinear
115c4762a1bSJed Brown      solver needs to evaluate the nonlinear function, it will call this
116c4762a1bSJed Brown      routine.
117c4762a1bSJed Brown       - Note that the final routine argument is the user-defined
118c4762a1bSJed Brown         context that provides application-specific data for the
119c4762a1bSJed Brown         function evaluation routine.
120c4762a1bSJed Brown   */
1219566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, FormFunction, (void *)&user));
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine
125c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   /*
128c4762a1bSJed Brown      Create matrix. Here we only approximately preallocate storage space
129c4762a1bSJed Brown      for the Jacobian.  See the users manual for a discussion of better
130c4762a1bSJed Brown      techniques for preallocating matrix memory.
131c4762a1bSJed Brown   */
1329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_mf", &matrix_free, NULL));
133c4762a1bSJed Brown   if (!matrix_free) {
134c4762a1bSJed Brown     PetscBool matrix_free_operator = PETSC_FALSE;
1359566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_mf_operator", &matrix_free_operator, NULL));
136c4762a1bSJed Brown     if (matrix_free_operator) matrix_free = PETSC_FALSE;
137c4762a1bSJed Brown   }
138*48a46eb9SPierre Jolivet   if (!matrix_free) PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, N, N, 5, NULL, &J));
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   /*
141c4762a1bSJed Brown      This option will cause the Jacobian to be computed via finite differences
142c4762a1bSJed Brown     efficiently using a coloring of the columns of the matrix.
143c4762a1bSJed Brown   */
1449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_fd_coloring", &fd_coloring, NULL));
145e00437b9SBarry Smith   PetscCheck(!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");
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   if (fd_coloring) {
148c4762a1bSJed Brown     ISColoring  iscoloring;
149c4762a1bSJed Brown     MatColoring mc;
150c4762a1bSJed Brown 
151c4762a1bSJed Brown     /*
152c4762a1bSJed Brown       This initializes the nonzero structure of the Jacobian. This is artificial
153c4762a1bSJed Brown       because clearly if we had a routine to compute the Jacobian we won't need
154c4762a1bSJed Brown       to use finite differences.
155c4762a1bSJed Brown     */
1569566063dSJacob Faibussowitsch     PetscCall(FormJacobian(snes, x, J, J, &user));
157c4762a1bSJed Brown 
158c4762a1bSJed Brown     /*
159c4762a1bSJed Brown        Color the matrix, i.e. determine groups of columns that share no common
160a5b23f4aSJose E. Roman       rows. These columns in the Jacobian can all be computed simultaneously.
161c4762a1bSJed Brown     */
1629566063dSJacob Faibussowitsch     PetscCall(MatColoringCreate(J, &mc));
1639566063dSJacob Faibussowitsch     PetscCall(MatColoringSetType(mc, MATCOLORINGSL));
1649566063dSJacob Faibussowitsch     PetscCall(MatColoringSetFromOptions(mc));
1659566063dSJacob Faibussowitsch     PetscCall(MatColoringApply(mc, &iscoloring));
1669566063dSJacob Faibussowitsch     PetscCall(MatColoringDestroy(&mc));
167c4762a1bSJed Brown     /*
168c4762a1bSJed Brown        Create the data structure that SNESComputeJacobianDefaultColor() uses
169c4762a1bSJed Brown        to compute the actual Jacobians via finite differences.
170c4762a1bSJed Brown     */
1719566063dSJacob Faibussowitsch     PetscCall(MatFDColoringCreate(J, iscoloring, &fdcoloring));
1729566063dSJacob Faibussowitsch     PetscCall(MatFDColoringSetFunction(fdcoloring, (PetscErrorCode(*)(void))FormFunction, &user));
1739566063dSJacob Faibussowitsch     PetscCall(MatFDColoringSetFromOptions(fdcoloring));
1749566063dSJacob Faibussowitsch     PetscCall(MatFDColoringSetUp(J, iscoloring, fdcoloring));
175c4762a1bSJed Brown     /*
176c4762a1bSJed Brown         Tell SNES to use the routine SNESComputeJacobianDefaultColor()
177c4762a1bSJed Brown       to compute Jacobians.
178c4762a1bSJed Brown     */
1799566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, J, J, SNESComputeJacobianDefaultColor, fdcoloring));
1809566063dSJacob Faibussowitsch     PetscCall(ISColoringDestroy(&iscoloring));
181c4762a1bSJed Brown   }
182c4762a1bSJed Brown   /*
183c4762a1bSJed Brown      Set Jacobian matrix data structure and default Jacobian evaluation
184c4762a1bSJed Brown      routine.  Whenever the nonlinear solver needs to compute the
185c4762a1bSJed Brown      Jacobian matrix, it will call this routine.
186c4762a1bSJed Brown       - Note that the final routine argument is the user-defined
187c4762a1bSJed Brown         context that provides application-specific data for the
188c4762a1bSJed Brown         Jacobian evaluation routine.
189c4762a1bSJed Brown       - The user can override with:
190c4762a1bSJed Brown          -snes_fd : default finite differencing approximation of Jacobian
191c4762a1bSJed Brown          -snes_mf : matrix-free Newton-Krylov method with no preconditioning
192c4762a1bSJed Brown                     (unless user explicitly sets preconditioner)
193c4762a1bSJed Brown          -snes_mf_operator : form preconditioning matrix as set by the user,
194c4762a1bSJed Brown                              but use matrix-free approx for Jacobian-vector
195c4762a1bSJed Brown                              products within Newton-Krylov method
196c4762a1bSJed Brown   */
197c4762a1bSJed Brown   else if (!matrix_free) {
1989566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, (void *)&user));
199c4762a1bSJed Brown   }
200c4762a1bSJed Brown 
201c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
202c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
203c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
204c4762a1bSJed Brown 
205c4762a1bSJed Brown   /*
206c4762a1bSJed Brown      Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
207c4762a1bSJed Brown   */
2089566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
209c4762a1bSJed Brown 
210c4762a1bSJed Brown   /*
211c4762a1bSJed Brown      Set array that saves the function norms.  This array is intended
212c4762a1bSJed Brown      when the user wants to save the convergence history for later use
213c4762a1bSJed Brown      rather than just to view the function norms via -snes_monitor.
214c4762a1bSJed Brown   */
2159566063dSJacob Faibussowitsch   PetscCall(SNESSetConvergenceHistory(snes, history, hist_its, 50, PETSC_TRUE));
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   /*
218c4762a1bSJed Brown       Add a user provided convergence test; this is to test that SNESNEWTONTR properly calls the
219c4762a1bSJed Brown       user provided test before the specialized test. The convergence context is just an array to
220c4762a1bSJed Brown       test that it gets properly freed at the end
221c4762a1bSJed Brown   */
222c4762a1bSJed Brown   if (use_convergence_test) {
2239566063dSJacob Faibussowitsch     PetscCall(SNESGetKSP(snes, &ksp));
2249566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5, &testarray));
2259566063dSJacob Faibussowitsch     PetscCall(KSPSetConvergenceTest(ksp, ConvergenceTest, testarray, ConvergenceDestroy));
226c4762a1bSJed Brown   }
227c4762a1bSJed Brown 
228c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
229c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
230c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
231c4762a1bSJed Brown   /*
232c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
233c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
234c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
235c4762a1bSJed Brown      this vector to zero by calling VecSet().
236c4762a1bSJed Brown   */
2379566063dSJacob Faibussowitsch   PetscCall(FormInitialGuess(&user, x));
2389566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, x));
2399566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
24063a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its));
241c4762a1bSJed Brown 
242c4762a1bSJed Brown   /*
243c4762a1bSJed Brown      Print the convergence history.  This is intended just to demonstrate
244c4762a1bSJed Brown      use of the data attained via SNESSetConvergenceHistory().
245c4762a1bSJed Brown   */
2469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-print_history", &flg));
247c4762a1bSJed Brown   if (flg) {
248*48a46eb9SPierre Jolivet     for (i = 0; i < its + 1; i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "iteration %" PetscInt_FMT ": Linear iterations %" PetscInt_FMT " Function norm = %g\n", i, hist_its[i], (double)history[i]));
249c4762a1bSJed Brown   }
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
252c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
253c4762a1bSJed Brown      are no longer needed.
254c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
255c4762a1bSJed Brown 
256*48a46eb9SPierre Jolivet   if (!matrix_free) PetscCall(MatDestroy(&J));
257*48a46eb9SPierre Jolivet   if (fd_coloring) PetscCall(MatFDColoringDestroy(&fdcoloring));
2589566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
2599566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
2609566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
2619566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
262b122ec5aSJacob Faibussowitsch   return 0;
263c4762a1bSJed Brown }
264c4762a1bSJed Brown /* ------------------------------------------------------------------- */
265c4762a1bSJed Brown /*
266c4762a1bSJed Brown    FormInitialGuess - Forms initial approximation.
267c4762a1bSJed Brown 
268c4762a1bSJed Brown    Input Parameters:
269c4762a1bSJed Brown    user - user-defined application context
270c4762a1bSJed Brown    X - vector
271c4762a1bSJed Brown 
272c4762a1bSJed Brown    Output Parameter:
273c4762a1bSJed Brown    X - vector
274c4762a1bSJed Brown  */
2759371c9d4SSatish Balay PetscErrorCode FormInitialGuess(AppCtx *user, Vec X) {
276c4762a1bSJed Brown   PetscInt     i, j, row, mx, my;
277c4762a1bSJed Brown   PetscReal    lambda, temp1, temp, hx, hy;
278c4762a1bSJed Brown   PetscScalar *x;
279c4762a1bSJed Brown 
280c4762a1bSJed Brown   mx     = user->mx;
281c4762a1bSJed Brown   my     = user->my;
282c4762a1bSJed Brown   lambda = user->param;
283c4762a1bSJed Brown 
284c4762a1bSJed Brown   hx = 1.0 / (PetscReal)(mx - 1);
285c4762a1bSJed Brown   hy = 1.0 / (PetscReal)(my - 1);
286c4762a1bSJed Brown 
287c4762a1bSJed Brown   /*
288c4762a1bSJed Brown      Get a pointer to vector data.
289c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
290c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
291c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
292c4762a1bSJed Brown          the array.
293c4762a1bSJed Brown   */
2949566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
295c4762a1bSJed Brown   temp1 = lambda / (lambda + 1.0);
296c4762a1bSJed Brown   for (j = 0; j < my; j++) {
297c4762a1bSJed Brown     temp = (PetscReal)(PetscMin(j, my - j - 1)) * hy;
298c4762a1bSJed Brown     for (i = 0; i < mx; i++) {
299c4762a1bSJed Brown       row = i + j * mx;
300c4762a1bSJed Brown       if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) {
301c4762a1bSJed Brown         x[row] = 0.0;
302c4762a1bSJed Brown         continue;
303c4762a1bSJed Brown       }
304c4762a1bSJed Brown       x[row] = temp1 * PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i, mx - i - 1)) * hx, temp));
305c4762a1bSJed Brown     }
306c4762a1bSJed Brown   }
307c4762a1bSJed Brown 
308c4762a1bSJed Brown   /*
309c4762a1bSJed Brown      Restore vector
310c4762a1bSJed Brown   */
3119566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, &x));
312c4762a1bSJed Brown   return 0;
313c4762a1bSJed Brown }
314c4762a1bSJed Brown /* ------------------------------------------------------------------- */
315c4762a1bSJed Brown /*
316c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
317c4762a1bSJed Brown 
318c4762a1bSJed Brown    Input Parameters:
319c4762a1bSJed Brown .  snes - the SNES context
320c4762a1bSJed Brown .  X - input vector
321c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetFunction()
322c4762a1bSJed Brown 
323c4762a1bSJed Brown    Output Parameter:
324c4762a1bSJed Brown .  F - function vector
325c4762a1bSJed Brown  */
3269371c9d4SSatish Balay PetscErrorCode FormFunction(SNES snes, Vec X, Vec F, void *ptr) {
327c4762a1bSJed Brown   AppCtx            *user = (AppCtx *)ptr;
328c4762a1bSJed Brown   PetscInt           i, j, row, mx, my;
329c4762a1bSJed Brown   PetscReal          two = 2.0, one = 1.0, lambda, hx, hy, hxdhy, hydhx;
330c4762a1bSJed Brown   PetscScalar        ut, ub, ul, ur, u, uxx, uyy, sc, *f;
331c4762a1bSJed Brown   const PetscScalar *x;
332c4762a1bSJed Brown 
333c4762a1bSJed Brown   mx     = user->mx;
334c4762a1bSJed Brown   my     = user->my;
335c4762a1bSJed Brown   lambda = user->param;
336c4762a1bSJed Brown   hx     = one / (PetscReal)(mx - 1);
337c4762a1bSJed Brown   hy     = one / (PetscReal)(my - 1);
338c4762a1bSJed Brown   sc     = hx * hy;
339c4762a1bSJed Brown   hxdhy  = hx / hy;
340c4762a1bSJed Brown   hydhx  = hy / hx;
341c4762a1bSJed Brown 
342c4762a1bSJed Brown   /*
343c4762a1bSJed Brown      Get pointers to vector data
344c4762a1bSJed Brown   */
3459566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
3469566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
347c4762a1bSJed Brown 
348c4762a1bSJed Brown   /*
349c4762a1bSJed Brown      Compute function
350c4762a1bSJed Brown   */
351c4762a1bSJed Brown   for (j = 0; j < my; j++) {
352c4762a1bSJed Brown     for (i = 0; i < mx; i++) {
353c4762a1bSJed Brown       row = i + j * mx;
354c4762a1bSJed Brown       if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) {
355c4762a1bSJed Brown         f[row] = x[row];
356c4762a1bSJed Brown         continue;
357c4762a1bSJed Brown       }
358c4762a1bSJed Brown       u      = x[row];
359c4762a1bSJed Brown       ub     = x[row - mx];
360c4762a1bSJed Brown       ul     = x[row - 1];
361c4762a1bSJed Brown       ut     = x[row + mx];
362c4762a1bSJed Brown       ur     = x[row + 1];
363c4762a1bSJed Brown       uxx    = (-ur + two * u - ul) * hydhx;
364c4762a1bSJed Brown       uyy    = (-ut + two * u - ub) * hxdhy;
365c4762a1bSJed Brown       f[row] = uxx + uyy - sc * lambda * PetscExpScalar(u);
366c4762a1bSJed Brown     }
367c4762a1bSJed Brown   }
368c4762a1bSJed Brown 
369c4762a1bSJed Brown   /*
370c4762a1bSJed Brown      Restore vectors
371c4762a1bSJed Brown   */
3729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
3739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
374c4762a1bSJed Brown   return 0;
375c4762a1bSJed Brown }
376c4762a1bSJed Brown /* ------------------------------------------------------------------- */
377c4762a1bSJed Brown /*
378c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
379c4762a1bSJed Brown 
380c4762a1bSJed Brown    Input Parameters:
381c4762a1bSJed Brown .  snes - the SNES context
382c4762a1bSJed Brown .  x - input vector
383c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetJacobian()
384c4762a1bSJed Brown 
385c4762a1bSJed Brown    Output Parameters:
386c4762a1bSJed Brown .  A - Jacobian matrix
387c4762a1bSJed Brown .  B - optionally different preconditioning matrix
388c4762a1bSJed Brown .  flag - flag indicating matrix structure
389c4762a1bSJed Brown */
3909371c9d4SSatish Balay PetscErrorCode FormJacobian(SNES snes, Vec X, Mat J, Mat jac, void *ptr) {
391c4762a1bSJed Brown   AppCtx            *user = (AppCtx *)ptr; /* user-defined applicatin context */
392c4762a1bSJed Brown   PetscInt           i, j, row, mx, my, col[5];
393c4762a1bSJed Brown   PetscScalar        two = 2.0, one = 1.0, lambda, v[5], sc;
394c4762a1bSJed Brown   const PetscScalar *x;
395c4762a1bSJed Brown   PetscReal          hx, hy, hxdhy, hydhx;
396c4762a1bSJed Brown 
397c4762a1bSJed Brown   mx     = user->mx;
398c4762a1bSJed Brown   my     = user->my;
399c4762a1bSJed Brown   lambda = user->param;
400c4762a1bSJed Brown   hx     = 1.0 / (PetscReal)(mx - 1);
401c4762a1bSJed Brown   hy     = 1.0 / (PetscReal)(my - 1);
402c4762a1bSJed Brown   sc     = hx * hy;
403c4762a1bSJed Brown   hxdhy  = hx / hy;
404c4762a1bSJed Brown   hydhx  = hy / hx;
405c4762a1bSJed Brown 
406c4762a1bSJed Brown   /*
407c4762a1bSJed Brown      Get pointer to vector data
408c4762a1bSJed Brown   */
4099566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
410c4762a1bSJed Brown 
411c4762a1bSJed Brown   /*
412c4762a1bSJed Brown      Compute entries of the Jacobian
413c4762a1bSJed Brown   */
414c4762a1bSJed Brown   for (j = 0; j < my; j++) {
415c4762a1bSJed Brown     for (i = 0; i < mx; i++) {
416c4762a1bSJed Brown       row = i + j * mx;
417c4762a1bSJed Brown       if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) {
4189566063dSJacob Faibussowitsch         PetscCall(MatSetValues(jac, 1, &row, 1, &row, &one, INSERT_VALUES));
419c4762a1bSJed Brown         continue;
420c4762a1bSJed Brown       }
4219371c9d4SSatish Balay       v[0]   = -hxdhy;
4229371c9d4SSatish Balay       col[0] = row - mx;
4239371c9d4SSatish Balay       v[1]   = -hydhx;
4249371c9d4SSatish Balay       col[1] = row - 1;
4259371c9d4SSatish Balay       v[2]   = two * (hydhx + hxdhy) - sc * lambda * PetscExpScalar(x[row]);
4269371c9d4SSatish Balay       col[2] = row;
4279371c9d4SSatish Balay       v[3]   = -hydhx;
4289371c9d4SSatish Balay       col[3] = row + 1;
4299371c9d4SSatish Balay       v[4]   = -hxdhy;
4309371c9d4SSatish Balay       col[4] = row + mx;
4319566063dSJacob Faibussowitsch       PetscCall(MatSetValues(jac, 1, &row, 5, col, v, INSERT_VALUES));
432c4762a1bSJed Brown     }
433c4762a1bSJed Brown   }
434c4762a1bSJed Brown 
435c4762a1bSJed Brown   /*
436c4762a1bSJed Brown      Restore vector
437c4762a1bSJed Brown   */
4389566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
439c4762a1bSJed Brown 
440c4762a1bSJed Brown   /*
441c4762a1bSJed Brown      Assemble matrix
442c4762a1bSJed Brown   */
4439566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
4449566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
445c4762a1bSJed Brown 
446c4762a1bSJed Brown   if (jac != J) {
4479566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
4489566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
449c4762a1bSJed Brown   }
450c4762a1bSJed Brown 
451c4762a1bSJed Brown   return 0;
452c4762a1bSJed Brown }
453c4762a1bSJed Brown 
4549371c9d4SSatish Balay PetscErrorCode ConvergenceTest(KSP ksp, PetscInt it, PetscReal nrm, KSPConvergedReason *reason, void *ctx) {
455c4762a1bSJed Brown   PetscFunctionBegin;
456c4762a1bSJed Brown   *reason = KSP_CONVERGED_ITERATING;
457c4762a1bSJed Brown   if (it > 1) {
458c4762a1bSJed Brown     *reason = KSP_CONVERGED_ITS;
4599566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "User provided convergence test returning after 2 iterations\n"));
460c4762a1bSJed Brown   }
461c4762a1bSJed Brown   PetscFunctionReturn(0);
462c4762a1bSJed Brown }
463c4762a1bSJed Brown 
4649371c9d4SSatish Balay PetscErrorCode ConvergenceDestroy(void *ctx) {
465c4762a1bSJed Brown   PetscFunctionBegin;
4669566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "User provided convergence destroy called\n"));
4679566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
468c4762a1bSJed Brown   PetscFunctionReturn(0);
469c4762a1bSJed Brown }
470c4762a1bSJed Brown 
4719371c9d4SSatish Balay PetscErrorCode postcheck(SNES snes, Vec x, Vec y, Vec w, PetscBool *changed_y, PetscBool *changed_w, void *ctx) {
472c4762a1bSJed Brown   PetscReal norm;
473c4762a1bSJed Brown   Vec       tmp;
474c4762a1bSJed Brown 
475c4762a1bSJed Brown   PetscFunctionBegin;
4769566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &tmp));
4779566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(tmp, -1.0, x, w));
4789566063dSJacob Faibussowitsch   PetscCall(VecNorm(tmp, NORM_2, &norm));
4799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp));
4809566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of search step %g\n", (double)norm));
481c4762a1bSJed Brown   PetscFunctionReturn(0);
482c4762a1bSJed Brown }
483c4762a1bSJed Brown 
484c4762a1bSJed Brown /*TEST
485c4762a1bSJed Brown 
486c4762a1bSJed Brown    build:
487c4762a1bSJed Brown       requires: !single
488c4762a1bSJed Brown 
489c4762a1bSJed Brown    test:
490c4762a1bSJed Brown       args: -ksp_gmres_cgs_refinement_type refine_always
491c4762a1bSJed Brown 
492c4762a1bSJed Brown    test:
493c4762a1bSJed Brown       suffix: 2
494c4762a1bSJed Brown       args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always
495c4762a1bSJed Brown 
496c4762a1bSJed Brown    test:
497c4762a1bSJed Brown       suffix: 2a
498c4762a1bSJed Brown       filter: grep -i KSPConvergedDefault > /dev/null && echo "Found KSPConvergedDefault"
499c4762a1bSJed Brown       args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always -info
500dfd57a17SPierre Jolivet       requires: defined(PETSC_USE_INFO)
501c4762a1bSJed Brown 
502c4762a1bSJed Brown    test:
503c4762a1bSJed Brown       suffix: 2b
504c4762a1bSJed Brown       filter: grep -i  "User provided convergence test" > /dev/null  && echo "Found User provided convergence test"
505c4762a1bSJed Brown       args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always -use_convergence_test -info
506dfd57a17SPierre Jolivet       requires: defined(PETSC_USE_INFO)
507c4762a1bSJed Brown 
508c4762a1bSJed Brown    test:
509c4762a1bSJed Brown       suffix: 3
510c4762a1bSJed Brown       args: -snes_monitor_short -mat_coloring_type sl -snes_fd_coloring -mx 8 -my 11 -ksp_gmres_cgs_refinement_type refine_always
511c4762a1bSJed Brown 
512c4762a1bSJed Brown    test:
513c4762a1bSJed Brown       suffix: 4
514c4762a1bSJed Brown       args: -pc -par 6.807 -snes_monitor -snes_converged_reason
515c4762a1bSJed Brown 
516c4762a1bSJed Brown TEST*/
517