xref: /petsc/src/snes/tutorials/ex2.c (revision 63a3b9bc7a1f24f247904ccba9383635fe6abade)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Newton method to solve u'' + u^{2} = f, sequentially.\n\
3c4762a1bSJed Brown This example employs a user-defined monitoring routine.\n\n";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown /*
6c4762a1bSJed Brown    Include "petscdraw.h" so that we can use PETSc drawing routines.
7c4762a1bSJed Brown    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
8c4762a1bSJed Brown    file automatically includes:
9c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
10c4762a1bSJed Brown      petscmat.h - matrices
11c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
12c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
13c4762a1bSJed Brown      petscksp.h   - linear solvers
14c4762a1bSJed Brown */
15c4762a1bSJed Brown 
16c4762a1bSJed Brown #include <petscsnes.h>
17c4762a1bSJed Brown 
18c4762a1bSJed Brown /*
19c4762a1bSJed Brown    User-defined routines
20c4762a1bSJed Brown */
21c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);
22c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
23c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(Vec);
24c4762a1bSJed Brown extern PetscErrorCode Monitor(SNES,PetscInt,PetscReal,void*);
25c4762a1bSJed Brown 
26c4762a1bSJed Brown /*
27c4762a1bSJed Brown    User-defined context for monitoring
28c4762a1bSJed Brown */
29c4762a1bSJed Brown typedef struct {
30c4762a1bSJed Brown   PetscViewer viewer;
31c4762a1bSJed Brown } MonitorCtx;
32c4762a1bSJed Brown 
33c4762a1bSJed Brown int main(int argc,char **argv)
34c4762a1bSJed Brown {
35c4762a1bSJed Brown   SNES           snes;                   /* SNES context */
36c4762a1bSJed Brown   Vec            x,r,F,U;             /* vectors */
37c4762a1bSJed Brown   Mat            J;                      /* Jacobian matrix */
38c4762a1bSJed Brown   MonitorCtx     monP;                   /* monitoring context */
39c4762a1bSJed Brown   PetscInt       its,n = 5,i,maxit,maxf;
40c4762a1bSJed Brown   PetscMPIInt    size;
41c4762a1bSJed Brown   PetscScalar    h,xp,v,none = -1.0;
42c4762a1bSJed Brown   PetscReal      abstol,rtol,stol,norm;
43c4762a1bSJed Brown 
449566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
459566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
46be096a46SBarry Smith   PetscCheck(size == 1,PETSC_COMM_SELF,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
48c4762a1bSJed Brown   h    = 1.0/(n-1);
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
51c4762a1bSJed Brown      Create nonlinear solver context
52c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
53c4762a1bSJed Brown 
549566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes));
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57c4762a1bSJed Brown      Create vector data structures; set function evaluation routine
58c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
59c4762a1bSJed Brown   /*
60c4762a1bSJed Brown      Note that we form 1 vector from scratch and then duplicate as needed.
61c4762a1bSJed Brown   */
629566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&x));
639566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x,PETSC_DECIDE,n));
649566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
659566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&r));
669566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&F));
679566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&U));
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   /*
70c4762a1bSJed Brown      Set function evaluation routine and vector
71c4762a1bSJed Brown   */
729566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes,r,FormFunction,(void*)F));
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine
76c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
77c4762a1bSJed Brown 
789566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&J));
799566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n));
809566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
819566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(J,3,NULL));
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   /*
84c4762a1bSJed Brown      Set Jacobian matrix data structure and default Jacobian evaluation
85c4762a1bSJed Brown      routine. User can override with:
86c4762a1bSJed Brown      -snes_fd : default finite differencing approximation of Jacobian
87c4762a1bSJed Brown      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
88c4762a1bSJed Brown                 (unless user explicitly sets preconditioner)
89c4762a1bSJed Brown      -snes_mf_operator : form preconditioning matrix as set by the user,
90c4762a1bSJed Brown                          but use matrix-free approx for Jacobian-vector
91c4762a1bSJed Brown                          products within Newton-Krylov method
92c4762a1bSJed Brown   */
93c4762a1bSJed Brown 
949566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes,J,J,FormJacobian,NULL));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
98c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /*
101c4762a1bSJed Brown      Set an optional user-defined monitoring routine
102c4762a1bSJed Brown   */
1039566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,0,0,0,400,400,&monP.viewer));
1049566063dSJacob Faibussowitsch   PetscCall(SNESMonitorSet(snes,Monitor,&monP,0));
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /*
107c4762a1bSJed Brown      Set names for some vectors to facilitate monitoring (optional)
108c4762a1bSJed Brown   */
1099566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)x,"Approximate Solution"));
1109566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)U,"Exact Solution"));
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   /*
113c4762a1bSJed Brown      Set SNES/KSP/KSP/PC runtime options, e.g.,
114c4762a1bSJed Brown          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
115c4762a1bSJed Brown   */
1169566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   /*
119c4762a1bSJed Brown      Print parameters used for convergence testing (optional) ... just
120c4762a1bSJed Brown      to demonstrate this routine; this information is also printed with
121c4762a1bSJed Brown      the option -snes_view
122c4762a1bSJed Brown   */
1239566063dSJacob Faibussowitsch   PetscCall(SNESGetTolerances(snes,&abstol,&rtol,&stol,&maxit,&maxf));
124*63a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"atol=%g, rtol=%g, stol=%g, maxit=%" PetscInt_FMT ", maxf=%" PetscInt_FMT "\n",(double)abstol,(double)rtol,(double)stol,maxit,maxf));
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127c4762a1bSJed Brown      Initialize application:
128c4762a1bSJed Brown      Store right-hand-side of PDE and exact solution
129c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   xp = 0.0;
132c4762a1bSJed Brown   for (i=0; i<n; i++) {
133c4762a1bSJed Brown     v    = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */
1349566063dSJacob Faibussowitsch     PetscCall(VecSetValues(F,1,&i,&v,INSERT_VALUES));
135c4762a1bSJed Brown     v    = xp*xp*xp;
1369566063dSJacob Faibussowitsch     PetscCall(VecSetValues(U,1,&i,&v,INSERT_VALUES));
137c4762a1bSJed Brown     xp  += h;
138c4762a1bSJed Brown   }
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
142c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
143c4762a1bSJed Brown   /*
144c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
145c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
146c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
147c4762a1bSJed Brown      this vector to zero by calling VecSet().
148c4762a1bSJed Brown   */
1499566063dSJacob Faibussowitsch   PetscCall(FormInitialGuess(x));
1509566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes,NULL,x));
1519566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes,&its));
152*63a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"number of SNES iterations = %" PetscInt_FMT "\n\n",its));
153c4762a1bSJed Brown 
154c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
155c4762a1bSJed Brown      Check solution and clean up
156c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
157c4762a1bSJed Brown 
158c4762a1bSJed Brown   /*
159c4762a1bSJed Brown      Check the error
160c4762a1bSJed Brown   */
1619566063dSJacob Faibussowitsch   PetscCall(VecAXPY(x,none,U));
1629566063dSJacob Faibussowitsch   PetscCall(VecNorm(x,NORM_2,&norm));
163*63a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g, Iterations %" PetscInt_FMT "\n",(double)norm,its));
164c4762a1bSJed Brown 
165c4762a1bSJed Brown   /*
166c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
167c4762a1bSJed Brown      are no longer needed.
168c4762a1bSJed Brown   */
1699566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));  PetscCall(VecDestroy(&r));
1709566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));  PetscCall(VecDestroy(&F));
1719566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));  PetscCall(SNESDestroy(&snes));
1729566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&monP.viewer));
1739566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
174b122ec5aSJacob Faibussowitsch   return 0;
175c4762a1bSJed Brown }
176c4762a1bSJed Brown /* ------------------------------------------------------------------- */
177c4762a1bSJed Brown /*
178c4762a1bSJed Brown    FormInitialGuess - Computes initial guess.
179c4762a1bSJed Brown 
180c4762a1bSJed Brown    Input/Output Parameter:
181c4762a1bSJed Brown .  x - the solution vector
182c4762a1bSJed Brown */
183c4762a1bSJed Brown PetscErrorCode FormInitialGuess(Vec x)
184c4762a1bSJed Brown {
185c4762a1bSJed Brown   PetscScalar    pfive = .50;
1869566063dSJacob Faibussowitsch   PetscCall(VecSet(x,pfive));
187c4762a1bSJed Brown   return 0;
188c4762a1bSJed Brown }
189c4762a1bSJed Brown /* ------------------------------------------------------------------- */
190c4762a1bSJed Brown /*
191c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
192c4762a1bSJed Brown 
193c4762a1bSJed Brown    Input Parameters:
194c4762a1bSJed Brown .  snes - the SNES context
195c4762a1bSJed Brown .  x - input vector
196c4762a1bSJed Brown .  ctx - optional user-defined context, as set by SNESSetFunction()
197c4762a1bSJed Brown 
198c4762a1bSJed Brown    Output Parameter:
199c4762a1bSJed Brown .  f - function vector
200c4762a1bSJed Brown 
201c4762a1bSJed Brown    Note:
202c4762a1bSJed Brown    The user-defined context can contain any application-specific data
203c4762a1bSJed Brown    needed for the function evaluation (such as various parameters, work
204c4762a1bSJed Brown    vectors, and grid information).  In this program the context is just
205c4762a1bSJed Brown    a vector containing the right-hand-side of the discretized PDE.
206c4762a1bSJed Brown  */
207c4762a1bSJed Brown 
208c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx)
209c4762a1bSJed Brown {
210c4762a1bSJed Brown   Vec               g = (Vec)ctx;
211c4762a1bSJed Brown   const PetscScalar *xx,*gg;
212c4762a1bSJed Brown   PetscScalar       *ff,d;
213c4762a1bSJed Brown   PetscInt          i,n;
214c4762a1bSJed Brown 
215c4762a1bSJed Brown   /*
216c4762a1bSJed Brown      Get pointers to vector data.
217c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
218c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
219c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
220c4762a1bSJed Brown          the array.
221c4762a1bSJed Brown   */
2229566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x,&xx));
2239566063dSJacob Faibussowitsch   PetscCall(VecGetArray(f,&ff));
2249566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(g,&gg));
225c4762a1bSJed Brown 
226c4762a1bSJed Brown   /*
227c4762a1bSJed Brown      Compute function
228c4762a1bSJed Brown   */
2299566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x,&n));
230c4762a1bSJed Brown   d     = (PetscReal)(n - 1); d = d*d;
231c4762a1bSJed Brown   ff[0] = xx[0];
232c4762a1bSJed Brown   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];
233c4762a1bSJed Brown   ff[n-1] = xx[n-1] - 1.0;
234c4762a1bSJed Brown 
235c4762a1bSJed Brown   /*
236c4762a1bSJed Brown      Restore vectors
237c4762a1bSJed Brown   */
2389566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x,&xx));
2399566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(f,&ff));
2409566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(g,&gg));
241c4762a1bSJed Brown   return 0;
242c4762a1bSJed Brown }
243c4762a1bSJed Brown /* ------------------------------------------------------------------- */
244c4762a1bSJed Brown /*
245c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
246c4762a1bSJed Brown 
247c4762a1bSJed Brown    Input Parameters:
248c4762a1bSJed Brown .  snes - the SNES context
249c4762a1bSJed Brown .  x - input vector
250c4762a1bSJed Brown .  dummy - optional user-defined context (not used here)
251c4762a1bSJed Brown 
252c4762a1bSJed Brown    Output Parameters:
253c4762a1bSJed Brown .  jac - Jacobian matrix
254c4762a1bSJed Brown .  B - optionally different preconditioning matrix
255c4762a1bSJed Brown 
256c4762a1bSJed Brown */
257c4762a1bSJed Brown 
258c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *dummy)
259c4762a1bSJed Brown {
260c4762a1bSJed Brown   const PetscScalar *xx;
261c4762a1bSJed Brown   PetscScalar       A[3],d;
262c4762a1bSJed Brown   PetscInt          i,n,j[3];
263c4762a1bSJed Brown 
264c4762a1bSJed Brown   /*
265c4762a1bSJed Brown      Get pointer to vector data
266c4762a1bSJed Brown   */
2679566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x,&xx));
268c4762a1bSJed Brown 
269c4762a1bSJed Brown   /*
270c4762a1bSJed Brown      Compute Jacobian entries and insert into matrix.
271c4762a1bSJed Brown       - Note that in this case we set all elements for a particular
272c4762a1bSJed Brown         row at once.
273c4762a1bSJed Brown   */
2749566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x,&n));
275c4762a1bSJed Brown   d    = (PetscReal)(n - 1); d = d*d;
276c4762a1bSJed Brown 
277c4762a1bSJed Brown   /*
278c4762a1bSJed Brown      Interior grid points
279c4762a1bSJed Brown   */
280c4762a1bSJed Brown   for (i=1; i<n-1; i++) {
281c4762a1bSJed Brown     j[0] = i - 1; j[1] = i; j[2] = i + 1;
282c4762a1bSJed Brown     A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i];
2839566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B,1,&i,3,j,A,INSERT_VALUES));
284c4762a1bSJed Brown   }
285c4762a1bSJed Brown 
286c4762a1bSJed Brown   /*
287c4762a1bSJed Brown      Boundary points
288c4762a1bSJed Brown   */
289c4762a1bSJed Brown   i = 0;   A[0] = 1.0;
290c4762a1bSJed Brown 
2919566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES));
292c4762a1bSJed Brown 
293c4762a1bSJed Brown   i = n-1; A[0] = 1.0;
294c4762a1bSJed Brown 
2959566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES));
296c4762a1bSJed Brown 
297c4762a1bSJed Brown   /*
298c4762a1bSJed Brown      Restore vector
299c4762a1bSJed Brown   */
3009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x,&xx));
301c4762a1bSJed Brown 
302c4762a1bSJed Brown   /*
303c4762a1bSJed Brown      Assemble matrix
304c4762a1bSJed Brown   */
3059566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
3069566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
307c4762a1bSJed Brown   if (jac != B) {
3089566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
3099566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
310c4762a1bSJed Brown   }
311c4762a1bSJed Brown   return 0;
312c4762a1bSJed Brown }
313c4762a1bSJed Brown /* ------------------------------------------------------------------- */
314c4762a1bSJed Brown /*
315c4762a1bSJed Brown    Monitor - User-defined monitoring routine that views the
316c4762a1bSJed Brown    current iterate with an x-window plot.
317c4762a1bSJed Brown 
318c4762a1bSJed Brown    Input Parameters:
319c4762a1bSJed Brown    snes - the SNES context
320c4762a1bSJed Brown    its - iteration number
321c4762a1bSJed Brown    norm - 2-norm function value (may be estimated)
322c4762a1bSJed Brown    ctx - optional user-defined context for private data for the
323c4762a1bSJed Brown          monitor routine, as set by SNESMonitorSet()
324c4762a1bSJed Brown 
325c4762a1bSJed Brown    Note:
326c4762a1bSJed Brown    See the manpage for PetscViewerDrawOpen() for useful runtime options,
327c4762a1bSJed Brown    such as -nox to deactivate all x-window output.
328c4762a1bSJed Brown  */
329c4762a1bSJed Brown PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal fnorm,void *ctx)
330c4762a1bSJed Brown {
331c4762a1bSJed Brown   MonitorCtx     *monP = (MonitorCtx*) ctx;
332c4762a1bSJed Brown   Vec            x;
333c4762a1bSJed Brown 
334*63a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"iter = %" PetscInt_FMT ", SNES Function norm %g\n",its,(double)fnorm));
3359566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes,&x));
3369566063dSJacob Faibussowitsch   PetscCall(VecView(x,monP->viewer));
337c4762a1bSJed Brown   return 0;
338c4762a1bSJed Brown }
339c4762a1bSJed Brown 
340c4762a1bSJed Brown /*TEST
341c4762a1bSJed Brown 
342c4762a1bSJed Brown    test:
343c4762a1bSJed Brown       args: -nox -snes_monitor_cancel -snes_monitor_short -snes_view -pc_type jacobi -ksp_gmres_cgs_refinement_type refine_always
344c4762a1bSJed Brown 
345c4762a1bSJed Brown    test:
346c4762a1bSJed Brown       suffix: 2
347c4762a1bSJed Brown       args: -nox -snes_monitor_cancel -snes_monitor_short -snes_type newtontr -snes_view
348c4762a1bSJed Brown       requires: !single
349c4762a1bSJed Brown 
350c4762a1bSJed Brown    test:
351c4762a1bSJed Brown       suffix: 3
352c4762a1bSJed Brown       args: -nox -malloc no -options_left no -snes_monitor_cancel -snes_monitor_short -snes_view -pc_type jacobi -ksp_gmres_cgs_refinement_type refine_always
353c4762a1bSJed Brown 
35441ba4c6cSHeeho Park    test:
35541ba4c6cSHeeho Park       suffix: 4
35641ba4c6cSHeeho Park       args: -nox -snes_monitor_cancel -snes_monitor_short -snes_type newtontrdc -snes_view
35741ba4c6cSHeeho Park       requires: !single
35841ba4c6cSHeeho Park 
359c4762a1bSJed Brown TEST*/
360