xref: /petsc/src/snes/tests/ex5.c (revision dd8e379b23d2ef935d8131fb74f7eb73fef09263)
1c4762a1bSJed Brown static char help[] = "Newton method to solve u'' + u^{2} = f, sequentially.\n\
2c4762a1bSJed Brown This example tests PCVPBJacobiSetBlocks().\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*
5c4762a1bSJed Brown    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
6c4762a1bSJed Brown    file automatically includes:
7c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
8c4762a1bSJed Brown      petscmat.h - matrices
9c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
10c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
11c4762a1bSJed Brown      petscksp.h   - linear solvers
12c4762a1bSJed Brown */
13c4762a1bSJed Brown 
14c4762a1bSJed Brown #include <petscsnes.h>
15c4762a1bSJed Brown 
16c4762a1bSJed Brown /*
17c4762a1bSJed Brown    User-defined routines
18c4762a1bSJed Brown */
19c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *);
20c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES, Vec, Vec, void *);
21c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(Vec);
22c4762a1bSJed Brown 
23d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
24d71ae5a4SJacob Faibussowitsch {
25c4762a1bSJed Brown   SNES        snes;       /* SNES context */
26c4762a1bSJed Brown   Vec         x, r, F, U; /* vectors */
27c4762a1bSJed Brown   Mat         J;          /* Jacobian matrix */
28c4762a1bSJed Brown   PetscInt    its, n = 5, i, maxit, maxf, lens[3] = {1, 2, 2};
29c4762a1bSJed Brown   PetscMPIInt size;
30c4762a1bSJed Brown   PetscScalar h, xp, v, none = -1.0;
31c4762a1bSJed Brown   PetscReal   abstol, rtol, stol, norm;
32c4762a1bSJed Brown   KSP         ksp;
33c4762a1bSJed Brown   PC          pc;
34c4762a1bSJed Brown 
35327415f7SBarry 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));
40c4762a1bSJed Brown   h = 1.0 / (n - 1);
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
43c4762a1bSJed Brown      Create nonlinear solver context
44c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
45c4762a1bSJed Brown 
469566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
479566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
489566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
499566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCVPBJACOBI));
50c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
51c4762a1bSJed Brown      Create vector data structures; set function evaluation routine
52c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
53c4762a1bSJed Brown   /*
54c4762a1bSJed Brown      Note that we form 1 vector from scratch and then duplicate as needed.
55c4762a1bSJed Brown   */
569566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
579566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, PETSC_DECIDE, n));
589566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
599566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
609566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &F));
619566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &U));
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   /*
64c4762a1bSJed Brown      Set function evaluation routine and vector
65c4762a1bSJed Brown   */
669566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, FormFunction, (void *)F));
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
69c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine
70c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
71c4762a1bSJed Brown 
729566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
739566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, n, n));
749566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
759566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(J, 3, NULL));
769566063dSJacob Faibussowitsch   PetscCall(MatSetVariableBlockSizes(J, 3, lens));
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   /*
79c4762a1bSJed Brown      Set Jacobian matrix data structure and default Jacobian evaluation
80c4762a1bSJed Brown      routine. User can override with:
81c4762a1bSJed Brown      -snes_fd : default finite differencing approximation of Jacobian
82c4762a1bSJed Brown      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
83c4762a1bSJed Brown                 (unless user explicitly sets preconditioner)
84c4762a1bSJed Brown      -snes_mf_operator : form preconditioning matrix as set by the user,
85c4762a1bSJed Brown                          but use matrix-free approx for Jacobian-vector
86c4762a1bSJed Brown                          products within Newton-Krylov method
87c4762a1bSJed Brown   */
88c4762a1bSJed Brown 
899566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, NULL));
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
92c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
93c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   /*
96c4762a1bSJed Brown      Set names for some vectors to facilitate monitoring (optional)
97c4762a1bSJed Brown   */
989566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)x, "Approximate Solution"));
999566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)U, "Exact Solution"));
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   /*
102c4762a1bSJed Brown      Set SNES/KSP/KSP/PC runtime options, e.g.,
103c4762a1bSJed Brown          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
104c4762a1bSJed Brown   */
1059566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   /*
108c4762a1bSJed Brown      Print parameters used for convergence testing (optional) ... just
109c4762a1bSJed Brown      to demonstrate this routine; this information is also printed with
110c4762a1bSJed Brown      the option -snes_view
111c4762a1bSJed Brown   */
1129566063dSJacob Faibussowitsch   PetscCall(SNESGetTolerances(snes, &abstol, &rtol, &stol, &maxit, &maxf));
11363a3b9bcSJacob 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));
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116c4762a1bSJed Brown      Initialize application:
117*dd8e379bSPierre Jolivet      Store right-hand side of PDE and exact solution
118c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   xp = 0.0;
121c4762a1bSJed Brown   for (i = 0; i < n; i++) {
122c4762a1bSJed Brown     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));
124c4762a1bSJed Brown     v = xp * xp * xp;
1259566063dSJacob Faibussowitsch     PetscCall(VecSetValues(U, 1, &i, &v, INSERT_VALUES));
126c4762a1bSJed Brown     xp += h;
127c4762a1bSJed Brown   }
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
131c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
132c4762a1bSJed Brown   /*
133c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
134c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
135c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
136c4762a1bSJed Brown      this vector to zero by calling VecSet().
137c4762a1bSJed Brown   */
1389566063dSJacob Faibussowitsch   PetscCall(FormInitialGuess(x));
1399566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, x));
1409566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
14163a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "number of SNES iterations = %" PetscInt_FMT "\n\n", its));
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144c4762a1bSJed Brown      Check solution and clean up
145c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   /*
148c4762a1bSJed Brown      Check the error
149c4762a1bSJed Brown   */
1509566063dSJacob Faibussowitsch   PetscCall(VecAXPY(x, none, U));
1519566063dSJacob Faibussowitsch   PetscCall(VecNorm(x, NORM_2, &norm));
15263a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of error %g, Iterations %" PetscInt_FMT "\n", (double)norm, its));
153c4762a1bSJed Brown 
154c4762a1bSJed Brown   /*
155c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
156c4762a1bSJed Brown      are no longer needed.
157c4762a1bSJed Brown   */
1589371c9d4SSatish Balay   PetscCall(VecDestroy(&x));
1599371c9d4SSatish Balay   PetscCall(VecDestroy(&r));
1609371c9d4SSatish Balay   PetscCall(VecDestroy(&U));
1619371c9d4SSatish Balay   PetscCall(VecDestroy(&F));
1629371c9d4SSatish Balay   PetscCall(MatDestroy(&J));
1639371c9d4SSatish Balay   PetscCall(SNESDestroy(&snes));
1649566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
165b122ec5aSJacob Faibussowitsch   return 0;
166c4762a1bSJed Brown }
167f6dfbefdSBarry Smith 
168c4762a1bSJed Brown /*
169c4762a1bSJed Brown    FormInitialGuess - Computes initial guess.
170c4762a1bSJed Brown 
171c4762a1bSJed Brown    Input/Output Parameter:
172c4762a1bSJed Brown .  x - the solution vector
173c4762a1bSJed Brown */
174d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(Vec x)
175d71ae5a4SJacob Faibussowitsch {
176c4762a1bSJed Brown   PetscScalar pfive = .50;
1773ba16761SJacob Faibussowitsch 
1783ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
1799566063dSJacob Faibussowitsch   PetscCall(VecSet(x, pfive));
1803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
181c4762a1bSJed Brown }
182f6dfbefdSBarry Smith 
183c4762a1bSJed Brown /*
184c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
185c4762a1bSJed Brown 
186c4762a1bSJed Brown    Input Parameters:
187c4762a1bSJed Brown .  snes - the SNES context
188c4762a1bSJed Brown .  x - input vector
189c4762a1bSJed Brown .  ctx - optional user-defined context, as set by SNESSetFunction()
190c4762a1bSJed Brown 
191c4762a1bSJed Brown    Output Parameter:
192c4762a1bSJed Brown .  f - function vector
193c4762a1bSJed Brown 
194c4762a1bSJed Brown    Note:
195c4762a1bSJed Brown    The user-defined context can contain any application-specific data
196c4762a1bSJed Brown    needed for the function evaluation (such as various parameters, work
197c4762a1bSJed Brown    vectors, and grid information).  In this program the context is just
198*dd8e379bSPierre Jolivet    a vector containing the right-hand side of the discretized PDE.
199c4762a1bSJed Brown  */
200c4762a1bSJed Brown 
201d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(SNES snes, Vec x, Vec f, void *ctx)
202d71ae5a4SJacob Faibussowitsch {
203c4762a1bSJed Brown   Vec                g = (Vec)ctx;
204c4762a1bSJed Brown   const PetscScalar *xx, *gg;
205c4762a1bSJed Brown   PetscScalar       *ff, d;
206c4762a1bSJed Brown   PetscInt           i, n;
207c4762a1bSJed Brown 
2083ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
209c4762a1bSJed Brown   /*
210c4762a1bSJed Brown      Get pointers to vector data.
211c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
212c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
213c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
214c4762a1bSJed Brown          the array.
215c4762a1bSJed Brown   */
2169566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
2179566063dSJacob Faibussowitsch   PetscCall(VecGetArray(f, &ff));
2189566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(g, &gg));
219c4762a1bSJed Brown 
220c4762a1bSJed Brown   /*
221c4762a1bSJed Brown      Compute function
222c4762a1bSJed Brown   */
2239566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x, &n));
2249371c9d4SSatish Balay   d     = (PetscReal)(n - 1);
2259371c9d4SSatish Balay   d     = d * d;
226c4762a1bSJed Brown   ff[0] = xx[0];
227c4762a1bSJed 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];
228c4762a1bSJed Brown   ff[n - 1] = xx[n - 1] - 1.0;
229c4762a1bSJed Brown 
230c4762a1bSJed Brown   /*
231c4762a1bSJed Brown      Restore vectors
232c4762a1bSJed Brown   */
2339566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
2349566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(f, &ff));
2359566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(g, &gg));
2363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
237c4762a1bSJed Brown }
238f6dfbefdSBarry Smith 
239c4762a1bSJed Brown /*
240c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
241c4762a1bSJed Brown 
242c4762a1bSJed Brown    Input Parameters:
243c4762a1bSJed Brown .  snes - the SNES context
244c4762a1bSJed Brown .  x - input vector
245c4762a1bSJed Brown .  dummy - optional user-defined context (not used here)
246c4762a1bSJed Brown 
247c4762a1bSJed Brown    Output Parameters:
248c4762a1bSJed Brown .  jac - Jacobian matrix
249c4762a1bSJed Brown .  B - optionally different preconditioning matrix
250c4762a1bSJed Brown 
251c4762a1bSJed Brown */
252c4762a1bSJed Brown 
253d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec x, Mat jac, Mat B, void *dummy)
254d71ae5a4SJacob Faibussowitsch {
255c4762a1bSJed Brown   const PetscScalar *xx;
256c4762a1bSJed Brown   PetscScalar        A[3], d;
257c4762a1bSJed Brown   PetscInt           i, n, j[3];
258c4762a1bSJed Brown 
2593ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
260c4762a1bSJed Brown   /*
261c4762a1bSJed Brown      Get pointer to vector data
262c4762a1bSJed Brown   */
2639566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
264c4762a1bSJed Brown 
265c4762a1bSJed Brown   /*
266c4762a1bSJed Brown      Compute Jacobian entries and insert into matrix.
267c4762a1bSJed Brown       - Note that in this case we set all elements for a particular
268c4762a1bSJed Brown         row at once.
269c4762a1bSJed Brown   */
2709566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x, &n));
2719371c9d4SSatish Balay   d = (PetscReal)(n - 1);
2729371c9d4SSatish Balay   d = d * d;
273c4762a1bSJed Brown 
274c4762a1bSJed Brown   /*
275c4762a1bSJed Brown      Interior grid points
276c4762a1bSJed Brown   */
277c4762a1bSJed Brown   for (i = 1; i < n - 1; i++) {
2789371c9d4SSatish Balay     j[0] = i - 1;
2799371c9d4SSatish Balay     j[1] = i;
2809371c9d4SSatish Balay     j[2] = i + 1;
2818d47944aSStefano Zampini     A[0] = d;
2829371c9d4SSatish Balay     A[1] = -2.0 * d + 2.0 * xx[i];
2838d47944aSStefano Zampini     A[2] = d;
2849566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B, 1, &i, 3, j, A, INSERT_VALUES));
285c4762a1bSJed Brown   }
286c4762a1bSJed Brown 
287c4762a1bSJed Brown   /*
288c4762a1bSJed Brown      Boundary points
289c4762a1bSJed Brown   */
2909371c9d4SSatish Balay   i    = 0;
2919371c9d4SSatish Balay   A[0] = 1.0;
292c4762a1bSJed Brown 
2939566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 1, &i, 1, &i, A, INSERT_VALUES));
294c4762a1bSJed Brown 
2959371c9d4SSatish Balay   i    = n - 1;
2969371c9d4SSatish Balay   A[0] = 1.0;
297c4762a1bSJed Brown 
2989566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 1, &i, 1, &i, A, INSERT_VALUES));
299c4762a1bSJed Brown 
300c4762a1bSJed Brown   /*
301c4762a1bSJed Brown      Restore vector
302c4762a1bSJed Brown   */
3039566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
304c4762a1bSJed Brown 
305c4762a1bSJed Brown   /*
306c4762a1bSJed Brown      Assemble matrix
307c4762a1bSJed Brown   */
3089566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
3099566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
310c4762a1bSJed Brown   if (jac != B) {
3119566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
3129566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
313c4762a1bSJed Brown   }
3143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
315c4762a1bSJed Brown }
316c4762a1bSJed Brown 
317c4762a1bSJed Brown /*TEST
318c4762a1bSJed Brown 
319f1be3500SJunchao Zhang    testset:
320c4762a1bSJed Brown      args: -snes_monitor_short -snes_view -ksp_monitor
321f1be3500SJunchao Zhang      output_file: output/ex5_1.out
322f1be3500SJunchao Zhang      filter: grep -v "type: seqaij"
323f1be3500SJunchao Zhang 
324f1be3500SJunchao Zhang      test:
325f1be3500SJunchao Zhang       suffix: 1
326f1be3500SJunchao Zhang 
327f1be3500SJunchao Zhang      test:
328f1be3500SJunchao Zhang       suffix: cuda
329f1be3500SJunchao Zhang       requires: cuda
330f1be3500SJunchao Zhang       args: -mat_type aijcusparse -vec_type cuda
331f1be3500SJunchao Zhang 
332f1be3500SJunchao Zhang      test:
333f1be3500SJunchao Zhang       suffix: kok
334f1be3500SJunchao Zhang       requires: kokkos_kernels
335f1be3500SJunchao Zhang       args: -mat_type aijkokkos -vec_type kokkos
336c4762a1bSJed Brown 
337c4762a1bSJed Brown    # this is just a test for SNESKSPTRASPOSEONLY and KSPSolveTranspose to behave properly
338c4762a1bSJed Brown    # the solution is wrong on purpose
339c4762a1bSJed Brown    test:
340c4762a1bSJed Brown       requires: !single !complex
341c4762a1bSJed Brown       suffix: transpose_only
342c4762a1bSJed Brown       args: -snes_monitor_short -snes_view -ksp_monitor -snes_type ksptransposeonly -pc_type ilu -snes_test_jacobian -snes_test_jacobian_view -ksp_view_rhs -ksp_view_solution -ksp_view_mat_explicit -ksp_view_preconditioned_operator_explicit
343c4762a1bSJed Brown 
344c4762a1bSJed Brown TEST*/
345