xref: /petsc/src/snes/tests/ex1.c (revision 9be84c52375fa69a6481fcef8556b9a7b38e14a0)
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 
9f6dfbefdSBarry Smith /*
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 
26f6dfbefdSBarry Smith */
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 
61d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
62d71ae5a4SJacob Faibussowitsch {
63c4762a1bSJed Brown   SNES          snes; /* nonlinear solver context */
64c4762a1bSJed Brown   Vec           x, r; /* solution, residual vectors */
65c4762a1bSJed Brown   Mat           J;    /* Jacobian matrix */
66c4762a1bSJed Brown   AppCtx        user; /* user-defined application context */
67c4762a1bSJed Brown   PetscInt      i, its, N, hist_its[50];
68c4762a1bSJed Brown   PetscMPIInt   size;
69c4762a1bSJed Brown   PetscReal     bratu_lambda_max = 6.81, bratu_lambda_min = 0., history[50];
70c4762a1bSJed Brown   MatFDColoring fdcoloring;
71*9be84c52SStefano Zampini   PetscBool     matrix_free = PETSC_FALSE, flg, fd_coloring = PETSC_FALSE, use_convergence_test = PETSC_FALSE, pc = PETSC_FALSE, prunejacobian = PETSC_FALSE, null_appctx = PETSC_TRUE;
72c4762a1bSJed Brown   KSP           ksp;
73c4762a1bSJed Brown   PetscInt     *testarray;
74c4762a1bSJed Brown 
75327415f7SBarry Smith   PetscFunctionBeginUser;
769566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
779566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
78be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /*
81c4762a1bSJed Brown      Initialize problem parameters
82c4762a1bSJed Brown   */
839371c9d4SSatish Balay   user.mx    = 4;
849371c9d4SSatish Balay   user.my    = 4;
859371c9d4SSatish Balay   user.param = 6.0;
869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &user.mx, NULL));
879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &user.my, NULL));
889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-par", &user.param, NULL));
899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-pc", &pc, NULL));
90e00437b9SBarry Smith   PetscCheck(user.param < bratu_lambda_max && user.param > bratu_lambda_min, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Lambda is out of range");
91c4762a1bSJed Brown   N = user.mx * user.my;
929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_convergence_test", &use_convergence_test, NULL));
9378e7fe0eSHong Zhang   PetscCall(PetscOptionsGetBool(NULL, NULL, "-prune_jacobian", &prunejacobian, NULL));
94*9be84c52SStefano Zampini   PetscCall(PetscOptionsGetBool(NULL, NULL, "-null_appctx", &null_appctx, NULL));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97c4762a1bSJed Brown      Create nonlinear solver context
98c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99c4762a1bSJed Brown 
1009566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   if (pc) {
1039566063dSJacob Faibussowitsch     PetscCall(SNESSetType(snes, SNESNEWTONTR));
1049566063dSJacob Faibussowitsch     PetscCall(SNESNewtonTRSetPostCheck(snes, postcheck, NULL));
105c4762a1bSJed Brown   }
106c4762a1bSJed Brown 
107*9be84c52SStefano Zampini   /* Test application context handling from Python */
108*9be84c52SStefano Zampini   if (!null_appctx) { PetscCall(SNESSetApplicationContext(snes, (void *)&user)); }
109*9be84c52SStefano Zampini 
110c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111c4762a1bSJed Brown      Create vector data structures; set function evaluation routine
112c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
113c4762a1bSJed Brown 
1149566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
1159566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, PETSC_DECIDE, N));
1169566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
1179566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
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   */
1279566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, FormFunction, (void *)&user));
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   */
1389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_mf", &matrix_free, NULL));
139c4762a1bSJed Brown   if (!matrix_free) {
140c4762a1bSJed Brown     PetscBool matrix_free_operator = PETSC_FALSE;
1419566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_mf_operator", &matrix_free_operator, NULL));
142c4762a1bSJed Brown     if (matrix_free_operator) matrix_free = PETSC_FALSE;
143c4762a1bSJed Brown   }
14448a46eb9SPierre Jolivet   if (!matrix_free) PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, N, N, 5, NULL, &J));
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   /*
147c4762a1bSJed Brown      This option will cause the Jacobian to be computed via finite differences
148c4762a1bSJed Brown     efficiently using a coloring of the columns of the matrix.
149c4762a1bSJed Brown   */
1509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_fd_coloring", &fd_coloring, NULL));
151e00437b9SBarry 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");
152c4762a1bSJed Brown 
153c4762a1bSJed Brown   if (fd_coloring) {
154c4762a1bSJed Brown     ISColoring  iscoloring;
155c4762a1bSJed Brown     MatColoring mc;
15678e7fe0eSHong Zhang     if (prunejacobian) {
15778e7fe0eSHong Zhang       /* Initialize x with random nonzero values so that the nonzeros in the Jacobian
15878e7fe0eSHong Zhang          can better reflect the sparsity structure of the Jacobian. */
15978e7fe0eSHong Zhang       PetscRandom rctx;
16078e7fe0eSHong Zhang       PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rctx));
16178e7fe0eSHong Zhang       PetscCall(PetscRandomSetInterval(rctx, 1.0, 2.0));
16278e7fe0eSHong Zhang       PetscCall(VecSetRandom(x, rctx));
16378e7fe0eSHong Zhang       PetscCall(PetscRandomDestroy(&rctx));
16478e7fe0eSHong Zhang     }
165c4762a1bSJed Brown 
166c4762a1bSJed Brown     /*
167c4762a1bSJed Brown       This initializes the nonzero structure of the Jacobian. This is artificial
168c4762a1bSJed Brown       because clearly if we had a routine to compute the Jacobian we won't need
169c4762a1bSJed Brown       to use finite differences.
170c4762a1bSJed Brown     */
1719566063dSJacob Faibussowitsch     PetscCall(FormJacobian(snes, x, J, J, &user));
172c4762a1bSJed Brown 
173c4762a1bSJed Brown     /*
174c4762a1bSJed Brown        Color the matrix, i.e. determine groups of columns that share no common
175a5b23f4aSJose E. Roman       rows. These columns in the Jacobian can all be computed simultaneously.
176c4762a1bSJed Brown     */
1779566063dSJacob Faibussowitsch     PetscCall(MatColoringCreate(J, &mc));
1789566063dSJacob Faibussowitsch     PetscCall(MatColoringSetType(mc, MATCOLORINGSL));
1799566063dSJacob Faibussowitsch     PetscCall(MatColoringSetFromOptions(mc));
1809566063dSJacob Faibussowitsch     PetscCall(MatColoringApply(mc, &iscoloring));
1819566063dSJacob Faibussowitsch     PetscCall(MatColoringDestroy(&mc));
182c4762a1bSJed Brown     /*
183c4762a1bSJed Brown        Create the data structure that SNESComputeJacobianDefaultColor() uses
184c4762a1bSJed Brown        to compute the actual Jacobians via finite differences.
185c4762a1bSJed Brown     */
1869566063dSJacob Faibussowitsch     PetscCall(MatFDColoringCreate(J, iscoloring, &fdcoloring));
1879566063dSJacob Faibussowitsch     PetscCall(MatFDColoringSetFunction(fdcoloring, (PetscErrorCode(*)(void))FormFunction, &user));
1889566063dSJacob Faibussowitsch     PetscCall(MatFDColoringSetFromOptions(fdcoloring));
1899566063dSJacob Faibussowitsch     PetscCall(MatFDColoringSetUp(J, iscoloring, fdcoloring));
190c4762a1bSJed Brown     /*
191c4762a1bSJed Brown         Tell SNES to use the routine SNESComputeJacobianDefaultColor()
192c4762a1bSJed Brown       to compute Jacobians.
193c4762a1bSJed Brown     */
1949566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, J, J, SNESComputeJacobianDefaultColor, fdcoloring));
1959566063dSJacob Faibussowitsch     PetscCall(ISColoringDestroy(&iscoloring));
19678e7fe0eSHong Zhang     if (prunejacobian) PetscCall(SNESPruneJacobianColor(snes, J, J));
197c4762a1bSJed Brown   }
198c4762a1bSJed Brown   /*
199c4762a1bSJed Brown      Set Jacobian matrix data structure and default Jacobian evaluation
200c4762a1bSJed Brown      routine.  Whenever the nonlinear solver needs to compute the
201c4762a1bSJed Brown      Jacobian matrix, it will call this routine.
202c4762a1bSJed Brown       - Note that the final routine argument is the user-defined
203c4762a1bSJed Brown         context that provides application-specific data for the
204c4762a1bSJed Brown         Jacobian evaluation routine.
205c4762a1bSJed Brown       - The user can override with:
206c4762a1bSJed Brown          -snes_fd : default finite differencing approximation of Jacobian
207c4762a1bSJed Brown          -snes_mf : matrix-free Newton-Krylov method with no preconditioning
208c4762a1bSJed Brown                     (unless user explicitly sets preconditioner)
209c4762a1bSJed Brown          -snes_mf_operator : form preconditioning matrix as set by the user,
210c4762a1bSJed Brown                              but use matrix-free approx for Jacobian-vector
211c4762a1bSJed Brown                              products within Newton-Krylov method
212c4762a1bSJed Brown   */
213c4762a1bSJed Brown   else if (!matrix_free) {
2149566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, (void *)&user));
215c4762a1bSJed Brown   }
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
218c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
219c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   /*
222c4762a1bSJed Brown      Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
223c4762a1bSJed Brown   */
2249566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
225c4762a1bSJed Brown 
226c4762a1bSJed Brown   /*
227c4762a1bSJed Brown      Set array that saves the function norms.  This array is intended
228c4762a1bSJed Brown      when the user wants to save the convergence history for later use
229c4762a1bSJed Brown      rather than just to view the function norms via -snes_monitor.
230c4762a1bSJed Brown   */
2319566063dSJacob Faibussowitsch   PetscCall(SNESSetConvergenceHistory(snes, history, hist_its, 50, PETSC_TRUE));
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   /*
234c4762a1bSJed Brown       Add a user provided convergence test; this is to test that SNESNEWTONTR properly calls the
235c4762a1bSJed Brown       user provided test before the specialized test. The convergence context is just an array to
236c4762a1bSJed Brown       test that it gets properly freed at the end
237c4762a1bSJed Brown   */
238c4762a1bSJed Brown   if (use_convergence_test) {
2399566063dSJacob Faibussowitsch     PetscCall(SNESGetKSP(snes, &ksp));
2409566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5, &testarray));
2419566063dSJacob Faibussowitsch     PetscCall(KSPSetConvergenceTest(ksp, ConvergenceTest, testarray, ConvergenceDestroy));
242c4762a1bSJed Brown   }
243c4762a1bSJed Brown 
244c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
245c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
246c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
247c4762a1bSJed Brown   /*
248c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
249c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
250c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
251c4762a1bSJed Brown      this vector to zero by calling VecSet().
252c4762a1bSJed Brown   */
2539566063dSJacob Faibussowitsch   PetscCall(FormInitialGuess(&user, x));
2549566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, x));
2559566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
25663a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its));
257c4762a1bSJed Brown 
258c4762a1bSJed Brown   /*
259c4762a1bSJed Brown      Print the convergence history.  This is intended just to demonstrate
260c4762a1bSJed Brown      use of the data attained via SNESSetConvergenceHistory().
261c4762a1bSJed Brown   */
2629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-print_history", &flg));
263c4762a1bSJed Brown   if (flg) {
26448a46eb9SPierre 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]));
265c4762a1bSJed Brown   }
266c4762a1bSJed Brown 
267c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
268c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
269c4762a1bSJed Brown      are no longer needed.
270c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
271c4762a1bSJed Brown 
27248a46eb9SPierre Jolivet   if (!matrix_free) PetscCall(MatDestroy(&J));
27348a46eb9SPierre Jolivet   if (fd_coloring) PetscCall(MatFDColoringDestroy(&fdcoloring));
2749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
2759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
2769566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
2779566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
278b122ec5aSJacob Faibussowitsch   return 0;
279c4762a1bSJed Brown }
280f6dfbefdSBarry Smith 
281c4762a1bSJed Brown /*
282c4762a1bSJed Brown    FormInitialGuess - Forms initial approximation.
283c4762a1bSJed Brown 
284c4762a1bSJed Brown    Input Parameters:
285c4762a1bSJed Brown    user - user-defined application context
286c4762a1bSJed Brown    X - vector
287c4762a1bSJed Brown 
288c4762a1bSJed Brown    Output Parameter:
289c4762a1bSJed Brown    X - vector
290c4762a1bSJed Brown  */
291d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(AppCtx *user, Vec X)
292d71ae5a4SJacob Faibussowitsch {
293c4762a1bSJed Brown   PetscInt     i, j, row, mx, my;
294c4762a1bSJed Brown   PetscReal    lambda, temp1, temp, hx, hy;
295c4762a1bSJed Brown   PetscScalar *x;
296c4762a1bSJed Brown 
2973ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
298c4762a1bSJed Brown   mx     = user->mx;
299c4762a1bSJed Brown   my     = user->my;
300c4762a1bSJed Brown   lambda = user->param;
301c4762a1bSJed Brown 
302c4762a1bSJed Brown   hx = 1.0 / (PetscReal)(mx - 1);
303c4762a1bSJed Brown   hy = 1.0 / (PetscReal)(my - 1);
304c4762a1bSJed Brown 
305c4762a1bSJed Brown   /*
306c4762a1bSJed Brown      Get a pointer to vector data.
307c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
308c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
309c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
310c4762a1bSJed Brown          the array.
311c4762a1bSJed Brown   */
3129566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
313c4762a1bSJed Brown   temp1 = lambda / (lambda + 1.0);
314c4762a1bSJed Brown   for (j = 0; j < my; j++) {
315c4762a1bSJed Brown     temp = (PetscReal)(PetscMin(j, my - j - 1)) * hy;
316c4762a1bSJed Brown     for (i = 0; i < mx; i++) {
317c4762a1bSJed Brown       row = i + j * mx;
318c4762a1bSJed Brown       if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) {
319c4762a1bSJed Brown         x[row] = 0.0;
320c4762a1bSJed Brown         continue;
321c4762a1bSJed Brown       }
322c4762a1bSJed Brown       x[row] = temp1 * PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i, mx - i - 1)) * hx, temp));
323c4762a1bSJed Brown     }
324c4762a1bSJed Brown   }
325c4762a1bSJed Brown 
326c4762a1bSJed Brown   /*
327c4762a1bSJed Brown      Restore vector
328c4762a1bSJed Brown   */
3299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, &x));
3303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
331c4762a1bSJed Brown }
332f6dfbefdSBarry Smith 
333c4762a1bSJed Brown /*
334c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
335c4762a1bSJed Brown 
336c4762a1bSJed Brown    Input Parameters:
337c4762a1bSJed Brown .  snes - the SNES context
338c4762a1bSJed Brown .  X - input vector
339c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetFunction()
340c4762a1bSJed Brown 
341c4762a1bSJed Brown    Output Parameter:
342c4762a1bSJed Brown .  F - function vector
343c4762a1bSJed Brown  */
344d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(SNES snes, Vec X, Vec F, void *ptr)
345d71ae5a4SJacob Faibussowitsch {
346c4762a1bSJed Brown   AppCtx            *user = (AppCtx *)ptr;
347c4762a1bSJed Brown   PetscInt           i, j, row, mx, my;
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 
3523ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
353c4762a1bSJed Brown   mx     = user->mx;
354c4762a1bSJed Brown   my     = user->my;
355c4762a1bSJed Brown   lambda = user->param;
356c4762a1bSJed Brown   hx     = one / (PetscReal)(mx - 1);
357c4762a1bSJed Brown   hy     = one / (PetscReal)(my - 1);
358c4762a1bSJed Brown   sc     = hx * hy;
359c4762a1bSJed Brown   hxdhy  = hx / hy;
360c4762a1bSJed Brown   hydhx  = hy / hx;
361c4762a1bSJed Brown 
362c4762a1bSJed Brown   /*
363c4762a1bSJed Brown      Get pointers to vector data
364c4762a1bSJed Brown   */
3659566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
3669566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
367c4762a1bSJed Brown 
368c4762a1bSJed Brown   /*
369c4762a1bSJed Brown      Compute function
370c4762a1bSJed Brown   */
371c4762a1bSJed Brown   for (j = 0; j < my; j++) {
372c4762a1bSJed Brown     for (i = 0; i < mx; i++) {
373c4762a1bSJed Brown       row = i + j * mx;
374c4762a1bSJed Brown       if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) {
375c4762a1bSJed Brown         f[row] = x[row];
376c4762a1bSJed Brown         continue;
377c4762a1bSJed Brown       }
378c4762a1bSJed Brown       u      = x[row];
379c4762a1bSJed Brown       ub     = x[row - mx];
380c4762a1bSJed Brown       ul     = x[row - 1];
381c4762a1bSJed Brown       ut     = x[row + mx];
382c4762a1bSJed Brown       ur     = x[row + 1];
383c4762a1bSJed Brown       uxx    = (-ur + two * u - ul) * hydhx;
384c4762a1bSJed Brown       uyy    = (-ut + two * u - ub) * hxdhy;
385c4762a1bSJed Brown       f[row] = uxx + uyy - sc * lambda * PetscExpScalar(u);
386c4762a1bSJed Brown     }
387c4762a1bSJed Brown   }
388c4762a1bSJed Brown 
389c4762a1bSJed Brown   /*
390c4762a1bSJed Brown      Restore vectors
391c4762a1bSJed Brown   */
3929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
3939566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
3943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
395c4762a1bSJed Brown }
396f6dfbefdSBarry Smith 
397c4762a1bSJed Brown /*
398c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
399c4762a1bSJed Brown 
400c4762a1bSJed Brown    Input Parameters:
401c4762a1bSJed Brown .  snes - the SNES context
402c4762a1bSJed Brown .  x - input vector
403c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetJacobian()
404c4762a1bSJed Brown 
405c4762a1bSJed Brown    Output Parameters:
406c4762a1bSJed Brown .  A - Jacobian matrix
407c4762a1bSJed Brown .  B - optionally different preconditioning matrix
408c4762a1bSJed Brown .  flag - flag indicating matrix structure
409c4762a1bSJed Brown */
410d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec X, Mat J, Mat jac, void *ptr)
411d71ae5a4SJacob Faibussowitsch {
412da81f932SPierre Jolivet   AppCtx            *user = (AppCtx *)ptr; /* user-defined application context */
413c4762a1bSJed Brown   PetscInt           i, j, row, mx, my, col[5];
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 
4183ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
419c4762a1bSJed Brown   mx     = user->mx;
420c4762a1bSJed Brown   my     = user->my;
421c4762a1bSJed Brown   lambda = user->param;
422c4762a1bSJed Brown   hx     = 1.0 / (PetscReal)(mx - 1);
423c4762a1bSJed Brown   hy     = 1.0 / (PetscReal)(my - 1);
424c4762a1bSJed Brown   sc     = hx * hy;
425c4762a1bSJed Brown   hxdhy  = hx / hy;
426c4762a1bSJed Brown   hydhx  = hy / hx;
427c4762a1bSJed Brown 
428c4762a1bSJed Brown   /*
429c4762a1bSJed Brown      Get pointer to vector data
430c4762a1bSJed Brown   */
4319566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
432c4762a1bSJed Brown 
433c4762a1bSJed Brown   /*
434c4762a1bSJed Brown      Compute entries of the Jacobian
435c4762a1bSJed Brown   */
436c4762a1bSJed Brown   for (j = 0; j < my; j++) {
437c4762a1bSJed Brown     for (i = 0; i < mx; i++) {
438c4762a1bSJed Brown       row = i + j * mx;
439c4762a1bSJed Brown       if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) {
4409566063dSJacob Faibussowitsch         PetscCall(MatSetValues(jac, 1, &row, 1, &row, &one, INSERT_VALUES));
441c4762a1bSJed Brown         continue;
442c4762a1bSJed Brown       }
4439371c9d4SSatish Balay       v[0]   = -hxdhy;
4449371c9d4SSatish Balay       col[0] = row - mx;
4459371c9d4SSatish Balay       v[1]   = -hydhx;
4469371c9d4SSatish Balay       col[1] = row - 1;
4479371c9d4SSatish Balay       v[2]   = two * (hydhx + hxdhy) - sc * lambda * PetscExpScalar(x[row]);
4489371c9d4SSatish Balay       col[2] = row;
4499371c9d4SSatish Balay       v[3]   = -hydhx;
4509371c9d4SSatish Balay       col[3] = row + 1;
4519371c9d4SSatish Balay       v[4]   = -hxdhy;
4529371c9d4SSatish Balay       col[4] = row + mx;
4539566063dSJacob Faibussowitsch       PetscCall(MatSetValues(jac, 1, &row, 5, col, v, INSERT_VALUES));
454c4762a1bSJed Brown     }
455c4762a1bSJed Brown   }
456c4762a1bSJed Brown 
457c4762a1bSJed Brown   /*
458c4762a1bSJed Brown      Restore vector
459c4762a1bSJed Brown   */
4609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
461c4762a1bSJed Brown 
462c4762a1bSJed Brown   /*
463c4762a1bSJed Brown      Assemble matrix
464c4762a1bSJed Brown   */
4659566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
4669566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
467c4762a1bSJed Brown 
468c4762a1bSJed Brown   if (jac != J) {
4699566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
4709566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
471c4762a1bSJed Brown   }
472c4762a1bSJed Brown 
4733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
474c4762a1bSJed Brown }
475c4762a1bSJed Brown 
476d71ae5a4SJacob Faibussowitsch PetscErrorCode ConvergenceTest(KSP ksp, PetscInt it, PetscReal nrm, KSPConvergedReason *reason, void *ctx)
477d71ae5a4SJacob Faibussowitsch {
4783ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
479c4762a1bSJed Brown   *reason = KSP_CONVERGED_ITERATING;
480c4762a1bSJed Brown   if (it > 1) {
481c4762a1bSJed Brown     *reason = KSP_CONVERGED_ITS;
4829566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "User provided convergence test returning after 2 iterations\n"));
483c4762a1bSJed Brown   }
4843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
485c4762a1bSJed Brown }
486c4762a1bSJed Brown 
487d71ae5a4SJacob Faibussowitsch PetscErrorCode ConvergenceDestroy(void *ctx)
488d71ae5a4SJacob Faibussowitsch {
4893ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
4909566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "User provided convergence destroy called\n"));
4919566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
4923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
493c4762a1bSJed Brown }
494c4762a1bSJed Brown 
495d71ae5a4SJacob Faibussowitsch PetscErrorCode postcheck(SNES snes, Vec x, Vec y, Vec w, PetscBool *changed_y, PetscBool *changed_w, void *ctx)
496d71ae5a4SJacob Faibussowitsch {
497c4762a1bSJed Brown   PetscReal norm;
498c4762a1bSJed Brown   Vec       tmp;
499c4762a1bSJed Brown 
5003ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
5019566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &tmp));
5029566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(tmp, -1.0, x, w));
5039566063dSJacob Faibussowitsch   PetscCall(VecNorm(tmp, NORM_2, &norm));
5049566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp));
5059566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of search step %g\n", (double)norm));
5063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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
525dfd57a17SPierre Jolivet       requires: defined(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
531dfd57a17SPierre Jolivet       requires: defined(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 
54178e7fe0eSHong Zhang    test:
54278e7fe0eSHong Zhang       suffix: 5
54378e7fe0eSHong Zhang       args: -snes_monitor_short -mat_coloring_type sl -snes_fd_coloring -mx 8 -my 11 -ksp_gmres_cgs_refinement_type refine_always -prune_jacobian
54478e7fe0eSHong Zhang       output_file: output/ex1_3.out
545ccb5f961SBarry Smith 
546ccb5f961SBarry Smith    test:
547ccb5f961SBarry Smith       suffix: 6
548ccb5f961SBarry Smith       args: -snes_monitor draw:image:testfile -viewer_view
549ccb5f961SBarry Smith 
550*9be84c52SStefano Zampini    test:
551*9be84c52SStefano Zampini       suffix: python
552*9be84c52SStefano Zampini       requires: petsc4py
553*9be84c52SStefano Zampini       args: -python -snes_type python -snes_python_type ex1.py:MySNES -snes_view -null_appctx {{0 1}separate output}
554*9be84c52SStefano Zampini       localrunfiles: ex1.py
555*9be84c52SStefano Zampini 
556c4762a1bSJed Brown TEST*/
557