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