xref: /petsc/src/snes/tutorials/ex6.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4421ceaSFande Kong 
2c4421ceaSFande Kong static char help[] = "Newton method to solve u'' + u^{2} = f, sequentially.\n\
3c4421ceaSFande Kong This example employs a user-defined reasonview routine.\n\n";
4c4421ceaSFande Kong 
5c4421ceaSFande Kong #include <petscsnes.h>
6c4421ceaSFande Kong 
7c4421ceaSFande Kong /*
8c4421ceaSFande Kong    User-defined routines
9c4421ceaSFande Kong */
10c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);
11c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
12c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode FormInitialGuess(Vec);
13c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode MySNESConvergedReasonView(SNES,void*);
14c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode MyKSPConvergedReasonView(KSP,void*);
15c4421ceaSFande Kong 
16c4421ceaSFande Kong /*
17c4421ceaSFande Kong    User-defined context for monitoring
18c4421ceaSFande Kong */
19c4421ceaSFande Kong typedef struct {
20c4421ceaSFande Kong   PetscViewer viewer;
21c4421ceaSFande Kong } ReasonViewCtx;
22c4421ceaSFande Kong 
23c4421ceaSFande Kong int main(int argc,char **argv)
24c4421ceaSFande Kong {
25c4421ceaSFande Kong   SNES           snes;                /* SNES context */
26c4421ceaSFande Kong   KSP            ksp;                 /* KSP context */
27c4421ceaSFande Kong   Vec            x,r,F,U;             /* vectors */
28c4421ceaSFande Kong   Mat            J;                   /* Jacobian matrix */
29c4421ceaSFande Kong   ReasonViewCtx     monP;             /* monitoring context */
30c4421ceaSFande Kong   PetscInt       its,n = 5,i;
31c4421ceaSFande Kong   PetscMPIInt    size;
32c4421ceaSFande Kong   PetscScalar    h,xp,v;
33c4421ceaSFande Kong   MPI_Comm       comm;
34c4421ceaSFande Kong 
35*327415f7SBarry Smith   PetscFunctionBeginUser;
369566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
379566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
38be096a46SBarry Smith   PetscCheck(size == 1,PETSC_COMM_SELF,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
40c4421ceaSFande Kong   h    = 1.0/(n-1);
41c4421ceaSFande Kong   comm = PETSC_COMM_WORLD;
42c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
43c4421ceaSFande Kong      Create nonlinear solver context
44c4421ceaSFande Kong      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
45c4421ceaSFande Kong 
469566063dSJacob Faibussowitsch   PetscCall(SNESCreate(comm,&snes));
47c4421ceaSFande Kong 
48c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
49c4421ceaSFande Kong      Create vector data structures; set function evaluation routine
50c4421ceaSFande Kong      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
51c4421ceaSFande Kong   /*
52c4421ceaSFande Kong      Note that we form 1 vector from scratch and then duplicate as needed.
53c4421ceaSFande Kong   */
549566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm,&x));
559566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x,PETSC_DECIDE,n));
569566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
579566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&r));
589566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&F));
599566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&U));
60c4421ceaSFande Kong 
61c4421ceaSFande Kong   /*
62c4421ceaSFande Kong      Set function evaluation routine and vector
63c4421ceaSFande Kong   */
649566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes,r,FormFunction,(void*)F));
65c4421ceaSFande Kong 
66c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67c4421ceaSFande Kong      Create matrix data structure; set Jacobian evaluation routine
68c4421ceaSFande Kong      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
69c4421ceaSFande Kong 
709566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm,&J));
719566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n));
729566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
739566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(J,3,NULL));
74c4421ceaSFande Kong 
75c4421ceaSFande Kong   /*
76c4421ceaSFande Kong      Set Jacobian matrix data structure and default Jacobian evaluation
77c4421ceaSFande Kong      routine. User can override with:
78c4421ceaSFande Kong      -snes_fd : default finite differencing approximation of Jacobian
79c4421ceaSFande Kong      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
80c4421ceaSFande Kong                 (unless user explicitly sets preconditioner)
81c4421ceaSFande Kong      -snes_mf_operator : form preconditioning matrix as set by the user,
82c4421ceaSFande Kong                          but use matrix-free approx for Jacobian-vector
83c4421ceaSFande Kong                          products within Newton-Krylov method
84c4421ceaSFande Kong   */
85c4421ceaSFande Kong 
869566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes,J,J,FormJacobian,NULL));
87c4421ceaSFande Kong 
88c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
89c4421ceaSFande Kong      Customize nonlinear solver; set runtime options
90c4421ceaSFande Kong    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
91c4421ceaSFande Kong 
92c4421ceaSFande Kong   /*
93c4421ceaSFande Kong      Set an optional user-defined reasonview routine
94c4421ceaSFande Kong   */
959566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIGetStdout(comm,&monP.viewer));
96c4421ceaSFande Kong   /* Just make sure we can not repeat addding the same function
97c4421ceaSFande Kong    * PETSc will be able to igore the repeated function
98c4421ceaSFande Kong    */
99c4421ceaSFande Kong   for (i=0; i<4; i++) {
1009566063dSJacob Faibussowitsch     PetscCall(SNESConvergedReasonViewSet(snes,MySNESConvergedReasonView,&monP,0));
101c4421ceaSFande Kong   }
1029566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes,&ksp));
103c4421ceaSFande Kong   /* Just make sure we can not repeat addding the same function
104c4421ceaSFande Kong    * PETSc will be able to igore the repeated function
105c4421ceaSFande Kong    */
106c4421ceaSFande Kong   for (i=0; i<4; i++) {
1079566063dSJacob Faibussowitsch     PetscCall(KSPConvergedReasonViewSet(ksp,MyKSPConvergedReasonView,&monP,0));
108c4421ceaSFande Kong   }
109c4421ceaSFande Kong   /*
110c4421ceaSFande Kong      Set SNES/KSP/KSP/PC runtime options, e.g.,
111c4421ceaSFande Kong          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
112c4421ceaSFande Kong   */
1139566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
114c4421ceaSFande Kong 
115c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116c4421ceaSFande Kong      Initialize application:
117c4421ceaSFande Kong      Store right-hand-side of PDE and exact solution
118c4421ceaSFande Kong    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119c4421ceaSFande Kong 
120c4421ceaSFande Kong   xp = 0.0;
121c4421ceaSFande Kong   for (i=0; i<n; i++) {
122c4421ceaSFande Kong     v    = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */
1239566063dSJacob Faibussowitsch     PetscCall(VecSetValues(F,1,&i,&v,INSERT_VALUES));
124c4421ceaSFande Kong     v    = xp*xp*xp;
1259566063dSJacob Faibussowitsch     PetscCall(VecSetValues(U,1,&i,&v,INSERT_VALUES));
126c4421ceaSFande Kong     xp  += h;
127c4421ceaSFande Kong   }
128c4421ceaSFande Kong 
129c4421ceaSFande Kong   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130c4421ceaSFande Kong      Evaluate initial guess; then solve nonlinear system
131c4421ceaSFande Kong    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
132c4421ceaSFande Kong   /*
133c4421ceaSFande Kong      Note: The user should initialize the vector, x, with the initial guess
134c4421ceaSFande Kong      for the nonlinear solver prior to calling SNESSolve().  In particular,
135c4421ceaSFande Kong      to employ an initial guess of zero, the user should explicitly set
136c4421ceaSFande Kong      this vector to zero by calling VecSet().
137c4421ceaSFande Kong   */
1389566063dSJacob Faibussowitsch   PetscCall(FormInitialGuess(x));
1399566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes,NULL,x));
1409566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes,&its));
141c4421ceaSFande Kong 
142c4421ceaSFande Kong   /*
143c4421ceaSFande Kong      Free work space.  All PETSc objects should be destroyed when they
144c4421ceaSFande Kong      are no longer needed.
145c4421ceaSFande Kong   */
1469566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));  PetscCall(VecDestroy(&r));
1479566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));  PetscCall(VecDestroy(&F));
1489566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));  PetscCall(SNESDestroy(&snes));
1499566063dSJacob Faibussowitsch   /*PetscCall(PetscViewerDestroy(&monP.viewer));*/
1509566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
151b122ec5aSJacob Faibussowitsch   return 0;
152c4421ceaSFande Kong }
153c4421ceaSFande Kong /* ------------------------------------------------------------------- */
154c4421ceaSFande Kong /*
155c4421ceaSFande Kong    FormInitialGuess - Computes initial guess.
156c4421ceaSFande Kong 
157c4421ceaSFande Kong    Input/Output Parameter:
158c4421ceaSFande Kong .  x - the solution vector
159c4421ceaSFande Kong */
160c4421ceaSFande Kong PetscErrorCode FormInitialGuess(Vec x)
161c4421ceaSFande Kong {
162c4421ceaSFande Kong   PetscScalar    pfive = .50;
1639566063dSJacob Faibussowitsch   PetscCall(VecSet(x,pfive));
164c4421ceaSFande Kong   return 0;
165c4421ceaSFande Kong }
166c4421ceaSFande Kong /* ------------------------------------------------------------------- */
167c4421ceaSFande Kong /*
168c4421ceaSFande Kong    FormFunction - Evaluates nonlinear function, F(x).
169c4421ceaSFande Kong 
170c4421ceaSFande Kong    Input Parameters:
171c4421ceaSFande Kong .  snes - the SNES context
172c4421ceaSFande Kong .  x - input vector
173c4421ceaSFande Kong .  ctx - optional user-defined context, as set by SNESSetFunction()
174c4421ceaSFande Kong 
175c4421ceaSFande Kong    Output Parameter:
176c4421ceaSFande Kong .  f - function vector
177c4421ceaSFande Kong 
178c4421ceaSFande Kong    Note:
179c4421ceaSFande Kong    The user-defined context can contain any application-specific data
180c4421ceaSFande Kong    needed for the function evaluation (such as various parameters, work
181c4421ceaSFande Kong    vectors, and grid information).  In this program the context is just
182c4421ceaSFande Kong    a vector containing the right-hand-side of the discretized PDE.
183c4421ceaSFande Kong  */
184c4421ceaSFande Kong 
185c4421ceaSFande Kong PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx)
186c4421ceaSFande Kong {
187c4421ceaSFande Kong   Vec               g = (Vec)ctx;
188c4421ceaSFande Kong   const PetscScalar *xx,*gg;
189c4421ceaSFande Kong   PetscScalar       *ff,d;
190c4421ceaSFande Kong   PetscInt          i,n;
191c4421ceaSFande Kong 
192c4421ceaSFande Kong   /*
193c4421ceaSFande Kong      Get pointers to vector data.
194c4421ceaSFande Kong        - For default PETSc vectors, VecGetArray() returns a pointer to
195c4421ceaSFande Kong          the data array.  Otherwise, the routine is implementation dependent.
196c4421ceaSFande Kong        - You MUST call VecRestoreArray() when you no longer need access to
197c4421ceaSFande Kong          the array.
198c4421ceaSFande Kong   */
1999566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x,&xx));
2009566063dSJacob Faibussowitsch   PetscCall(VecGetArray(f,&ff));
2019566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(g,&gg));
202c4421ceaSFande Kong 
203c4421ceaSFande Kong   /*
204c4421ceaSFande Kong      Compute function
205c4421ceaSFande Kong   */
2069566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x,&n));
207c4421ceaSFande Kong   d     = (PetscReal)(n - 1); d = d*d;
208c4421ceaSFande Kong   ff[0] = xx[0];
209c4421ceaSFande Kong   for (i=1; i<n-1; i++) ff[i] = d*(xx[i-1] - 2.0*xx[i] + xx[i+1]) + xx[i]*xx[i] - gg[i];
210c4421ceaSFande Kong   ff[n-1] = xx[n-1] - 1.0;
211c4421ceaSFande Kong 
212c4421ceaSFande Kong   /*
213c4421ceaSFande Kong      Restore vectors
214c4421ceaSFande Kong   */
2159566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x,&xx));
2169566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(f,&ff));
2179566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(g,&gg));
218c4421ceaSFande Kong   return 0;
219c4421ceaSFande Kong }
220c4421ceaSFande Kong /* ------------------------------------------------------------------- */
221c4421ceaSFande Kong /*
222c4421ceaSFande Kong    FormJacobian - Evaluates Jacobian matrix.
223c4421ceaSFande Kong 
224c4421ceaSFande Kong    Input Parameters:
225c4421ceaSFande Kong .  snes - the SNES context
226c4421ceaSFande Kong .  x - input vector
227c4421ceaSFande Kong .  dummy - optional user-defined context (not used here)
228c4421ceaSFande Kong 
229c4421ceaSFande Kong    Output Parameters:
230c4421ceaSFande Kong .  jac - Jacobian matrix
231c4421ceaSFande Kong .  B - optionally different preconditioning matrix
232c4421ceaSFande Kong 
233c4421ceaSFande Kong */
234c4421ceaSFande Kong 
235c4421ceaSFande Kong PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *dummy)
236c4421ceaSFande Kong {
237c4421ceaSFande Kong   const PetscScalar *xx;
238c4421ceaSFande Kong   PetscScalar       A[3],d;
239c4421ceaSFande Kong   PetscInt          i,n,j[3];
240c4421ceaSFande Kong 
241c4421ceaSFande Kong   /*
242c4421ceaSFande Kong      Get pointer to vector data
243c4421ceaSFande Kong   */
2449566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x,&xx));
245c4421ceaSFande Kong 
246c4421ceaSFande Kong   /*
247c4421ceaSFande Kong      Compute Jacobian entries and insert into matrix.
248c4421ceaSFande Kong       - Note that in this case we set all elements for a particular
249c4421ceaSFande Kong         row at once.
250c4421ceaSFande Kong   */
2519566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x,&n));
252c4421ceaSFande Kong   d    = (PetscReal)(n - 1); d = d*d;
253c4421ceaSFande Kong 
254c4421ceaSFande Kong   /*
255c4421ceaSFande Kong      Interior grid points
256c4421ceaSFande Kong   */
257c4421ceaSFande Kong   for (i=1; i<n-1; i++) {
258c4421ceaSFande Kong     j[0] = i - 1; j[1] = i; j[2] = i + 1;
259c4421ceaSFande Kong     A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i];
2609566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B,1,&i,3,j,A,INSERT_VALUES));
261c4421ceaSFande Kong   }
262c4421ceaSFande Kong 
263c4421ceaSFande Kong   /*
264c4421ceaSFande Kong      Boundary points
265c4421ceaSFande Kong   */
266c4421ceaSFande Kong   i = 0;   A[0] = 1.0;
267c4421ceaSFande Kong 
2689566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES));
269c4421ceaSFande Kong 
270c4421ceaSFande Kong   i = n-1; A[0] = 1.0;
271c4421ceaSFande Kong 
2729566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES));
273c4421ceaSFande Kong 
274c4421ceaSFande Kong   /*
275c4421ceaSFande Kong      Restore vector
276c4421ceaSFande Kong   */
2779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x,&xx));
278c4421ceaSFande Kong 
279c4421ceaSFande Kong   /*
280c4421ceaSFande Kong      Assemble matrix
281c4421ceaSFande Kong   */
2829566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
2839566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
284c4421ceaSFande Kong   if (jac != B) {
2859566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
2869566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
287c4421ceaSFande Kong   }
288c4421ceaSFande Kong   return 0;
289c4421ceaSFande Kong }
290c4421ceaSFande Kong 
291c4421ceaSFande Kong PetscErrorCode MySNESConvergedReasonView(SNES snes,void *ctx)
292c4421ceaSFande Kong {
293c4421ceaSFande Kong   ReasonViewCtx         *monP = (ReasonViewCtx*) ctx;
294c4421ceaSFande Kong   PetscViewer           viewer = monP->viewer;
295c4421ceaSFande Kong   SNESConvergedReason   reason;
296c4421ceaSFande Kong   const char            *strreason;
297c4421ceaSFande Kong 
2989566063dSJacob Faibussowitsch   PetscCall(SNESGetConvergedReason(snes,&reason));
2999566063dSJacob Faibussowitsch   PetscCall(SNESGetConvergedReasonString(snes,&strreason));
3009566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Customized SNES converged reason view\n"));
3019566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer,1));
302c4421ceaSFande Kong   if (reason > 0) {
3039566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"Converged due to %s\n",strreason));
304c4421ceaSFande Kong   } else if (reason <= 0) {
3059566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"Did not converge due to %s\n",strreason));
306c4421ceaSFande Kong   }
3079566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISubtractTab(viewer,1));
308c4421ceaSFande Kong   return 0;
309c4421ceaSFande Kong }
310c4421ceaSFande Kong 
311c4421ceaSFande Kong PetscErrorCode MyKSPConvergedReasonView(KSP ksp,void *ctx)
312c4421ceaSFande Kong {
313c4421ceaSFande Kong   ReasonViewCtx         *monP = (ReasonViewCtx*) ctx;
314c4421ceaSFande Kong   PetscViewer           viewer = monP->viewer;
315c4421ceaSFande Kong   KSPConvergedReason    reason;
316c4421ceaSFande Kong   const char            *reasonstr;
317c4421ceaSFande Kong 
3189566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(ksp,&reason));
3199566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReasonString(ksp,&reasonstr));
3209566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer,2));
3219566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Customized KSP converged reason view\n"));
3229566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer,1));
323c4421ceaSFande Kong   if (reason > 0) {
3249566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"Converged due to %s\n",reasonstr));
325c4421ceaSFande Kong   } else if (reason <= 0) {
3269566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"Did not converge due to %s\n",reasonstr));
327c4421ceaSFande Kong   }
3289566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISubtractTab(viewer,3));
329c4421ceaSFande Kong   return 0;
330c4421ceaSFande Kong }
331c4421ceaSFande Kong 
332c4421ceaSFande Kong /*TEST
333c4421ceaSFande Kong 
334c4421ceaSFande Kong    test:
335c4421ceaSFande Kong       suffix: 1
336c4421ceaSFande Kong       nsize: 1
337bc44df21SStefano Zampini       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g"
338c4421ceaSFande Kong 
339c4421ceaSFande Kong    test:
340c4421ceaSFande Kong       suffix: 2
341c4421ceaSFande Kong       nsize: 1
342c4421ceaSFande Kong       args: -ksp_converged_reason_view_cancel
343bc44df21SStefano Zampini       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g"
344c4421ceaSFande Kong 
345c4421ceaSFande Kong    test:
346c4421ceaSFande Kong       suffix: 3
347c4421ceaSFande Kong       nsize: 1
348c4421ceaSFande Kong       args: -ksp_converged_reason_view_cancel -ksp_converged_reason
349bc44df21SStefano Zampini       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g"
350c4421ceaSFande Kong 
351c4421ceaSFande Kong    test:
352c4421ceaSFande Kong       suffix: 4
353c4421ceaSFande Kong       nsize: 1
354c4421ceaSFande Kong       args: -snes_converged_reason_view_cancel
355bc44df21SStefano Zampini       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g"
356c4421ceaSFande Kong 
357c4421ceaSFande Kong    test:
358c4421ceaSFande Kong       suffix: 5
359c4421ceaSFande Kong       nsize: 1
360c4421ceaSFande Kong       args: -snes_converged_reason_view_cancel -snes_converged_reason
361bc44df21SStefano Zampini       filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g"
362c4421ceaSFande Kong 
363c4421ceaSFande Kong TEST*/
364