xref: /petsc/src/snes/tutorials/ex1.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Newton's method for a two-variable system, sequential.\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*T
5c4762a1bSJed Brown    Concepts: SNES^basic example
6c4762a1bSJed Brown T*/
7c4762a1bSJed Brown 
8c4762a1bSJed Brown /*
9c4762a1bSJed Brown    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
10c4762a1bSJed Brown    file automatically includes:
11c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
12c4762a1bSJed Brown      petscmat.h - matrices
13c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
14c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
15c4762a1bSJed Brown      petscksp.h   - linear solvers
16c4762a1bSJed Brown */
17c4762a1bSJed Brown /*F
18c4762a1bSJed Brown This examples solves either
19c4762a1bSJed Brown \begin{equation}
20c4762a1bSJed Brown   F\genfrac{(}{)}{0pt}{}{x_0}{x_1} = \genfrac{(}{)}{0pt}{}{x^2_0 + x_0 x_1 - 3}{x_0 x_1 + x^2_1 - 6}
21c4762a1bSJed Brown \end{equation}
22c4762a1bSJed Brown or if the {\tt -hard} options is given
23c4762a1bSJed Brown \begin{equation}
24c4762a1bSJed Brown   F\genfrac{(}{)}{0pt}{}{x_0}{x_1} = \genfrac{(}{)}{0pt}{}{\sin(3 x_0) + x_0}{x_1}
25c4762a1bSJed Brown \end{equation}
26c4762a1bSJed Brown F*/
27c4762a1bSJed Brown #include <petscsnes.h>
28c4762a1bSJed Brown 
29c4762a1bSJed Brown /*
30c4762a1bSJed Brown    User-defined routines
31c4762a1bSJed Brown */
32c4762a1bSJed Brown extern PetscErrorCode FormJacobian1(SNES,Vec,Mat,Mat,void*);
33c4762a1bSJed Brown extern PetscErrorCode FormFunction1(SNES,Vec,Vec,void*);
34c4762a1bSJed Brown extern PetscErrorCode FormJacobian2(SNES,Vec,Mat,Mat,void*);
35c4762a1bSJed Brown extern PetscErrorCode FormFunction2(SNES,Vec,Vec,void*);
36c4762a1bSJed Brown 
37c4762a1bSJed Brown int main(int argc,char **argv)
38c4762a1bSJed Brown {
39c4762a1bSJed Brown   SNES           snes;         /* nonlinear solver context */
40c4762a1bSJed Brown   KSP            ksp;          /* linear solver context */
41c4762a1bSJed Brown   PC             pc;           /* preconditioner context */
42c4762a1bSJed Brown   Vec            x,r;          /* solution, residual vectors */
43c4762a1bSJed Brown   Mat            J;            /* Jacobian matrix */
44c4762a1bSJed Brown   PetscErrorCode ierr;
45c4762a1bSJed Brown   PetscMPIInt    size;
46c4762a1bSJed Brown   PetscScalar    pfive = .5,*xx;
47c4762a1bSJed Brown   PetscBool      flg;
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
50*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
512c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size > 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Example is only for sequential runs");
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54c4762a1bSJed Brown      Create nonlinear solver context
55c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
56*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes));
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
59c4762a1bSJed Brown      Create matrix and vector data structures; set corresponding routines
60c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61c4762a1bSJed Brown   /*
62c4762a1bSJed Brown      Create vectors for solution and nonlinear function
63c4762a1bSJed Brown   */
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x));
65*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(x,PETSC_DECIDE,2));
66*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(x));
67*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&r));
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   /*
70c4762a1bSJed Brown      Create Jacobian matrix data structure
71c4762a1bSJed Brown   */
72*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&J));
73*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,2,2));
74*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(J));
75*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(J));
76c4762a1bSJed Brown 
77*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-hard",&flg));
78c4762a1bSJed Brown   if (!flg) {
79c4762a1bSJed Brown     /*
80c4762a1bSJed Brown      Set function evaluation routine and vector.
81c4762a1bSJed Brown     */
82*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESSetFunction(snes,r,FormFunction1,NULL));
83c4762a1bSJed Brown 
84c4762a1bSJed Brown     /*
85c4762a1bSJed Brown      Set Jacobian matrix data structure and Jacobian evaluation routine
86c4762a1bSJed Brown     */
87*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian1,NULL));
88c4762a1bSJed Brown   } else {
89*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESSetFunction(snes,r,FormFunction2,NULL));
90*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian2,NULL));
91c4762a1bSJed Brown   }
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
95c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96c4762a1bSJed Brown   /*
97c4762a1bSJed Brown      Set linear solver defaults for this problem. By extracting the
98c4762a1bSJed Brown      KSP and PC contexts from the SNES context, we can then
99c4762a1bSJed Brown      directly call any KSP and PC routines to set various options.
100c4762a1bSJed Brown   */
101*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetKSP(snes,&ksp));
102*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetPC(ksp,&pc));
103*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetType(pc,PCNONE));
104*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetTolerances(ksp,1.e-4,PETSC_DEFAULT,PETSC_DEFAULT,20));
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /*
107c4762a1bSJed Brown      Set SNES/KSP/KSP/PC runtime options, e.g.,
108c4762a1bSJed Brown          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
109c4762a1bSJed Brown      These options will override those specified above as long as
110c4762a1bSJed Brown      SNESSetFromOptions() is called _after_ any other customization
111c4762a1bSJed Brown      routines.
112c4762a1bSJed Brown   */
113*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(snes));
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
117c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
118c4762a1bSJed Brown   if (!flg) {
119*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(x,pfive));
120c4762a1bSJed Brown   } else {
121*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(x,&xx));
122c4762a1bSJed Brown     xx[0] = 2.0; xx[1] = 3.0;
123*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(x,&xx));
124c4762a1bSJed Brown   }
125c4762a1bSJed Brown   /*
126c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
127c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
128c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
129c4762a1bSJed Brown      this vector to zero by calling VecSet().
130c4762a1bSJed Brown   */
131c4762a1bSJed Brown 
132*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(snes,NULL,x));
133c4762a1bSJed Brown   if (flg) {
134c4762a1bSJed Brown     Vec f;
135*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
136*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESGetFunction(snes,&f,0,0));
137*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(r,PETSC_VIEWER_STDOUT_WORLD));
138c4762a1bSJed Brown   }
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
142c4762a1bSJed Brown      are no longer needed.
143c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144c4762a1bSJed Brown 
145*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x)); CHKERRQ(VecDestroy(&r));
146*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&J)); CHKERRQ(SNESDestroy(&snes));
147c4762a1bSJed Brown   ierr = PetscFinalize();
148c4762a1bSJed Brown   return ierr;
149c4762a1bSJed Brown }
150c4762a1bSJed Brown /* ------------------------------------------------------------------- */
151c4762a1bSJed Brown /*
152c4762a1bSJed Brown    FormFunction1 - Evaluates nonlinear function, F(x).
153c4762a1bSJed Brown 
154c4762a1bSJed Brown    Input Parameters:
155c4762a1bSJed Brown .  snes - the SNES context
156c4762a1bSJed Brown .  x    - input vector
157c4762a1bSJed Brown .  ctx  - optional user-defined context
158c4762a1bSJed Brown 
159c4762a1bSJed Brown    Output Parameter:
160c4762a1bSJed Brown .  f - function vector
161c4762a1bSJed Brown  */
162c4762a1bSJed Brown PetscErrorCode FormFunction1(SNES snes,Vec x,Vec f,void *ctx)
163c4762a1bSJed Brown {
164c4762a1bSJed Brown   const PetscScalar *xx;
165c4762a1bSJed Brown   PetscScalar       *ff;
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   /*
168c4762a1bSJed Brown    Get pointers to vector data.
169c4762a1bSJed Brown       - For default PETSc vectors, VecGetArray() returns a pointer to
170c4762a1bSJed Brown         the data array.  Otherwise, the routine is implementation dependent.
171c4762a1bSJed Brown       - You MUST call VecRestoreArray() when you no longer need access to
172c4762a1bSJed Brown         the array.
173c4762a1bSJed Brown    */
174*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(x,&xx));
175*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(f,&ff));
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   /* Compute function */
178c4762a1bSJed Brown   ff[0] = xx[0]*xx[0] + xx[0]*xx[1] - 3.0;
179c4762a1bSJed Brown   ff[1] = xx[0]*xx[1] + xx[1]*xx[1] - 6.0;
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   /* Restore vectors */
182*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(x,&xx));
183*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(f,&ff));
184c4762a1bSJed Brown   return 0;
185c4762a1bSJed Brown }
186c4762a1bSJed Brown /* ------------------------------------------------------------------- */
187c4762a1bSJed Brown /*
188c4762a1bSJed Brown    FormJacobian1 - Evaluates Jacobian matrix.
189c4762a1bSJed Brown 
190c4762a1bSJed Brown    Input Parameters:
191c4762a1bSJed Brown .  snes - the SNES context
192c4762a1bSJed Brown .  x - input vector
193c4762a1bSJed Brown .  dummy - optional user-defined context (not used here)
194c4762a1bSJed Brown 
195c4762a1bSJed Brown    Output Parameters:
196c4762a1bSJed Brown .  jac - Jacobian matrix
197c4762a1bSJed Brown .  B - optionally different preconditioning matrix
198c4762a1bSJed Brown .  flag - flag indicating matrix structure
199c4762a1bSJed Brown */
200c4762a1bSJed Brown PetscErrorCode FormJacobian1(SNES snes,Vec x,Mat jac,Mat B,void *dummy)
201c4762a1bSJed Brown {
202c4762a1bSJed Brown   const PetscScalar *xx;
203c4762a1bSJed Brown   PetscScalar       A[4];
204c4762a1bSJed Brown   PetscInt          idx[2] = {0,1};
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   /*
207c4762a1bSJed Brown      Get pointer to vector data
208c4762a1bSJed Brown   */
209*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(x,&xx));
210c4762a1bSJed Brown 
211c4762a1bSJed Brown   /*
212c4762a1bSJed Brown      Compute Jacobian entries and insert into matrix.
213c4762a1bSJed Brown       - Since this is such a small problem, we set all entries for
214c4762a1bSJed Brown         the matrix at once.
215c4762a1bSJed Brown   */
216c4762a1bSJed Brown   A[0]  = 2.0*xx[0] + xx[1]; A[1] = xx[0];
217c4762a1bSJed Brown   A[2]  = xx[1]; A[3] = xx[0] + 2.0*xx[1];
218*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(B,2,idx,2,idx,A,INSERT_VALUES));
219c4762a1bSJed Brown 
220c4762a1bSJed Brown   /*
221c4762a1bSJed Brown      Restore vector
222c4762a1bSJed Brown   */
223*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(x,&xx));
224c4762a1bSJed Brown 
225c4762a1bSJed Brown   /*
226c4762a1bSJed Brown      Assemble matrix
227c4762a1bSJed Brown   */
228*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
229*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
230c4762a1bSJed Brown   if (jac != B) {
231*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
232*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
233c4762a1bSJed Brown   }
234c4762a1bSJed Brown   return 0;
235c4762a1bSJed Brown }
236c4762a1bSJed Brown 
237c4762a1bSJed Brown /* ------------------------------------------------------------------- */
238c4762a1bSJed Brown PetscErrorCode FormFunction2(SNES snes,Vec x,Vec f,void *dummy)
239c4762a1bSJed Brown {
240c4762a1bSJed Brown   const PetscScalar *xx;
241c4762a1bSJed Brown   PetscScalar       *ff;
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   /*
244c4762a1bSJed Brown      Get pointers to vector data.
245c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
246c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
247c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
248c4762a1bSJed Brown          the array.
249c4762a1bSJed Brown   */
250*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(x,&xx));
251*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(f,&ff));
252c4762a1bSJed Brown 
253c4762a1bSJed Brown   /*
254c4762a1bSJed Brown      Compute function
255c4762a1bSJed Brown   */
256c4762a1bSJed Brown   ff[0] = PetscSinScalar(3.0*xx[0]) + xx[0];
257c4762a1bSJed Brown   ff[1] = xx[1];
258c4762a1bSJed Brown 
259c4762a1bSJed Brown   /*
260c4762a1bSJed Brown      Restore vectors
261c4762a1bSJed Brown   */
262*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(x,&xx));
263*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(f,&ff));
264c4762a1bSJed Brown   return 0;
265c4762a1bSJed Brown }
266c4762a1bSJed Brown /* ------------------------------------------------------------------- */
267c4762a1bSJed Brown PetscErrorCode FormJacobian2(SNES snes,Vec x,Mat jac,Mat B,void *dummy)
268c4762a1bSJed Brown {
269c4762a1bSJed Brown   const PetscScalar *xx;
270c4762a1bSJed Brown   PetscScalar       A[4];
271c4762a1bSJed Brown   PetscInt          idx[2] = {0,1};
272c4762a1bSJed Brown 
273c4762a1bSJed Brown   /*
274c4762a1bSJed Brown      Get pointer to vector data
275c4762a1bSJed Brown   */
276*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(x,&xx));
277c4762a1bSJed Brown 
278c4762a1bSJed Brown   /*
279c4762a1bSJed Brown      Compute Jacobian entries and insert into matrix.
280c4762a1bSJed Brown       - Since this is such a small problem, we set all entries for
281c4762a1bSJed Brown         the matrix at once.
282c4762a1bSJed Brown   */
283c4762a1bSJed Brown   A[0]  = 3.0*PetscCosScalar(3.0*xx[0]) + 1.0; A[1] = 0.0;
284c4762a1bSJed Brown   A[2]  = 0.0;                     A[3] = 1.0;
285*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(B,2,idx,2,idx,A,INSERT_VALUES));
286c4762a1bSJed Brown 
287c4762a1bSJed Brown   /*
288c4762a1bSJed Brown      Restore vector
289c4762a1bSJed Brown   */
290*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(x,&xx));
291c4762a1bSJed Brown 
292c4762a1bSJed Brown   /*
293c4762a1bSJed Brown      Assemble matrix
294c4762a1bSJed Brown   */
295*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
296*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
297c4762a1bSJed Brown   if (jac != B) {
298*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
299*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
300c4762a1bSJed Brown   }
301c4762a1bSJed Brown   return 0;
302c4762a1bSJed Brown }
303c4762a1bSJed Brown 
304c4762a1bSJed Brown /*TEST
305c4762a1bSJed Brown 
306c4762a1bSJed Brown    test:
307c4762a1bSJed Brown       args: -ksp_gmres_cgs_refinement_type refine_always -snes_monitor_short
308c4762a1bSJed Brown       requires: !single
309c4762a1bSJed Brown 
310c4762a1bSJed Brown    test:
311c4762a1bSJed Brown       suffix: 2
312c4762a1bSJed Brown       requires: !single
313c4762a1bSJed Brown       args:  -snes_monitor_short
314c4762a1bSJed Brown       output_file: output/ex1_1.out
315c4762a1bSJed Brown 
316c4762a1bSJed Brown    test:
317c4762a1bSJed Brown       suffix: 3
318c4762a1bSJed Brown       args: -ksp_view_solution ascii:ex1_2_sol.tmp:ascii_matlab  -snes_monitor_short
319c4762a1bSJed Brown       requires: !single
320c4762a1bSJed Brown       output_file: output/ex1_1.out
321c4762a1bSJed Brown 
322c4762a1bSJed Brown    test:
323c4762a1bSJed Brown       suffix: 4
324c4762a1bSJed Brown       args: -ksp_view_solution ascii:ex1_2_sol.tmp::append  -snes_monitor_short
325c4762a1bSJed Brown       requires: !single
326c4762a1bSJed Brown       output_file: output/ex1_1.out
327c4762a1bSJed Brown 
328c4762a1bSJed Brown    test:
329c4762a1bSJed Brown       suffix: 5
330c4762a1bSJed Brown       args: -ksp_view_solution ascii:ex1_2_sol.tmp:ascii_matlab:append  -snes_monitor_short
331c4762a1bSJed Brown       requires: !single
332c4762a1bSJed Brown       output_file: output/ex1_1.out
333c4762a1bSJed Brown 
334c4762a1bSJed Brown    test:
335c4762a1bSJed Brown       suffix: 6
336c4762a1bSJed Brown       args: -ksp_view_solution ascii:ex1_2_sol.tmp:default:append  -snes_monitor_short
337c4762a1bSJed Brown       requires: !single
338c4762a1bSJed Brown       output_file: output/ex1_1.out
339c4762a1bSJed Brown 
340c4762a1bSJed Brown    test:
341c4762a1bSJed Brown       suffix: X
342c4762a1bSJed Brown       args: -ksp_monitor_short -ksp_type gmres -ksp_gmres_krylov_monitor -snes_monitor_short -snes_rtol 1.e-4
343c4762a1bSJed Brown       requires: !single x
344c4762a1bSJed Brown 
345c4762a1bSJed Brown TEST*/
346