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