xref: /petsc/src/snes/tests/ex5.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Newton method to solve u'' + u^{2} = f, sequentially.\n\
3*c4762a1bSJed Brown This example tests PCVPBJacobiSetBlocks().\n\n";
4*c4762a1bSJed Brown 
5*c4762a1bSJed Brown /*T
6*c4762a1bSJed Brown    Concepts: SNES^basic uniprocessor example
7*c4762a1bSJed Brown    Processors: 1
8*c4762a1bSJed Brown T*/
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 
20*c4762a1bSJed Brown #include <petscsnes.h>
21*c4762a1bSJed Brown 
22*c4762a1bSJed Brown /*
23*c4762a1bSJed Brown    User-defined routines
24*c4762a1bSJed Brown */
25*c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);
26*c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
27*c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(Vec);
28*c4762a1bSJed Brown 
29*c4762a1bSJed Brown int main(int argc,char **argv)
30*c4762a1bSJed Brown {
31*c4762a1bSJed Brown   SNES           snes;                   /* SNES context */
32*c4762a1bSJed Brown   Vec            x,r,F,U;                /* vectors */
33*c4762a1bSJed Brown   Mat            J;                      /* Jacobian matrix */
34*c4762a1bSJed Brown   PetscErrorCode ierr;
35*c4762a1bSJed Brown   PetscInt       its,n = 5,i,maxit,maxf,lens[3] = {1,2,2};
36*c4762a1bSJed Brown   PetscMPIInt    size;
37*c4762a1bSJed Brown   PetscScalar    h,xp,v,none = -1.0;
38*c4762a1bSJed Brown   PetscReal      abstol,rtol,stol,norm;
39*c4762a1bSJed Brown   KSP            ksp;
40*c4762a1bSJed Brown   PC             pc;
41*c4762a1bSJed Brown 
42*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
43*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
44*c4762a1bSJed Brown   if (size != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"This is a uniprocessor example only!");
45*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr);
46*c4762a1bSJed Brown   h    = 1.0/(n-1);
47*c4762a1bSJed Brown 
48*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
49*c4762a1bSJed Brown      Create nonlinear solver context
50*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
51*c4762a1bSJed Brown 
52*c4762a1bSJed Brown   ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr);
53*c4762a1bSJed Brown   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
54*c4762a1bSJed Brown   ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
55*c4762a1bSJed Brown   ierr = PCSetType(pc,PCVPBJACOBI);CHKERRQ(ierr);
56*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57*c4762a1bSJed Brown      Create vector data structures; set function evaluation routine
58*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
59*c4762a1bSJed Brown   /*
60*c4762a1bSJed Brown      Note that we form 1 vector from scratch and then duplicate as needed.
61*c4762a1bSJed Brown   */
62*c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr);
63*c4762a1bSJed Brown   ierr = VecSetSizes(x,PETSC_DECIDE,n);CHKERRQ(ierr);
64*c4762a1bSJed Brown   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
65*c4762a1bSJed Brown   ierr = VecDuplicate(x,&r);CHKERRQ(ierr);
66*c4762a1bSJed Brown   ierr = VecDuplicate(x,&F);CHKERRQ(ierr);
67*c4762a1bSJed Brown   ierr = VecDuplicate(x,&U);CHKERRQ(ierr);
68*c4762a1bSJed Brown 
69*c4762a1bSJed Brown   /*
70*c4762a1bSJed Brown      Set function evaluation routine and vector
71*c4762a1bSJed Brown   */
72*c4762a1bSJed Brown   ierr = SNESSetFunction(snes,r,FormFunction,(void*)F);CHKERRQ(ierr);
73*c4762a1bSJed Brown 
74*c4762a1bSJed Brown 
75*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
76*c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine
77*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
78*c4762a1bSJed Brown 
79*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&J);CHKERRQ(ierr);
80*c4762a1bSJed Brown   ierr = MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n);CHKERRQ(ierr);
81*c4762a1bSJed Brown   ierr = MatSetFromOptions(J);CHKERRQ(ierr);
82*c4762a1bSJed Brown   ierr = MatSeqAIJSetPreallocation(J,3,NULL);CHKERRQ(ierr);
83*c4762a1bSJed Brown   ierr = MatSetVariableBlockSizes(J,3,lens);CHKERRQ(ierr);
84*c4762a1bSJed Brown 
85*c4762a1bSJed Brown   /*
86*c4762a1bSJed Brown      Set Jacobian matrix data structure and default Jacobian evaluation
87*c4762a1bSJed Brown      routine. User can override with:
88*c4762a1bSJed Brown      -snes_fd : default finite differencing approximation of Jacobian
89*c4762a1bSJed Brown      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
90*c4762a1bSJed Brown                 (unless user explicitly sets preconditioner)
91*c4762a1bSJed Brown      -snes_mf_operator : form preconditioning matrix as set by the user,
92*c4762a1bSJed Brown                          but use matrix-free approx for Jacobian-vector
93*c4762a1bSJed Brown                          products within Newton-Krylov method
94*c4762a1bSJed Brown   */
95*c4762a1bSJed Brown 
96*c4762a1bSJed Brown   ierr = SNESSetJacobian(snes,J,J,FormJacobian,NULL);CHKERRQ(ierr);
97*c4762a1bSJed Brown 
98*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
99*c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
100*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
101*c4762a1bSJed Brown 
102*c4762a1bSJed Brown   /*
103*c4762a1bSJed Brown      Set names for some vectors to facilitate monitoring (optional)
104*c4762a1bSJed Brown   */
105*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject)x,"Approximate Solution");CHKERRQ(ierr);
106*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject)U,"Exact Solution");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   */
112*c4762a1bSJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
113*c4762a1bSJed Brown 
114*c4762a1bSJed Brown   /*
115*c4762a1bSJed Brown      Print parameters used for convergence testing (optional) ... just
116*c4762a1bSJed Brown      to demonstrate this routine; this information is also printed with
117*c4762a1bSJed Brown      the option -snes_view
118*c4762a1bSJed Brown   */
119*c4762a1bSJed Brown   ierr = SNESGetTolerances(snes,&abstol,&rtol,&stol,&maxit,&maxf);CHKERRQ(ierr);
120*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"atol=%g, rtol=%g, stol=%g, maxit=%D, maxf=%D\n",(double)abstol,(double)rtol,(double)stol,maxit,maxf);CHKERRQ(ierr);
121*c4762a1bSJed Brown 
122*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123*c4762a1bSJed Brown      Initialize application:
124*c4762a1bSJed Brown      Store right-hand-side of PDE and exact solution
125*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
126*c4762a1bSJed Brown 
127*c4762a1bSJed Brown   xp = 0.0;
128*c4762a1bSJed Brown   for (i=0; i<n; i++) {
129*c4762a1bSJed Brown     v    = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */
130*c4762a1bSJed Brown     ierr = VecSetValues(F,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr);
131*c4762a1bSJed Brown     v    = xp*xp*xp;
132*c4762a1bSJed Brown     ierr = VecSetValues(U,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr);
133*c4762a1bSJed Brown     xp  += h;
134*c4762a1bSJed Brown   }
135*c4762a1bSJed Brown 
136*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137*c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
138*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139*c4762a1bSJed Brown   /*
140*c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
141*c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
142*c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
143*c4762a1bSJed Brown      this vector to zero by calling VecSet().
144*c4762a1bSJed Brown   */
145*c4762a1bSJed Brown   ierr = FormInitialGuess(x);CHKERRQ(ierr);
146*c4762a1bSJed Brown   ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr);
147*c4762a1bSJed Brown   ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
148*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"number of SNES iterations = %D\n\n",its);CHKERRQ(ierr);
149*c4762a1bSJed Brown 
150*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151*c4762a1bSJed Brown      Check solution and clean up
152*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153*c4762a1bSJed Brown 
154*c4762a1bSJed Brown   /*
155*c4762a1bSJed Brown      Check the error
156*c4762a1bSJed Brown   */
157*c4762a1bSJed Brown   ierr = VecAXPY(x,none,U);CHKERRQ(ierr);
158*c4762a1bSJed Brown   ierr = VecNorm(x,NORM_2,&norm);CHKERRQ(ierr);
159*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g, Iterations %D\n",(double)norm,its);CHKERRQ(ierr);
160*c4762a1bSJed Brown 
161*c4762a1bSJed Brown 
162*c4762a1bSJed Brown   /*
163*c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
164*c4762a1bSJed Brown      are no longer needed.
165*c4762a1bSJed Brown   */
166*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);  ierr = VecDestroy(&r);CHKERRQ(ierr);
167*c4762a1bSJed Brown   ierr = VecDestroy(&U);CHKERRQ(ierr);  ierr = VecDestroy(&F);CHKERRQ(ierr);
168*c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);  ierr = SNESDestroy(&snes);CHKERRQ(ierr);
169*c4762a1bSJed Brown   ierr = PetscFinalize();
170*c4762a1bSJed Brown   return ierr;
171*c4762a1bSJed Brown }
172*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
173*c4762a1bSJed Brown /*
174*c4762a1bSJed Brown    FormInitialGuess - Computes initial guess.
175*c4762a1bSJed Brown 
176*c4762a1bSJed Brown    Input/Output Parameter:
177*c4762a1bSJed Brown .  x - the solution vector
178*c4762a1bSJed Brown */
179*c4762a1bSJed Brown PetscErrorCode FormInitialGuess(Vec x)
180*c4762a1bSJed Brown {
181*c4762a1bSJed Brown   PetscErrorCode ierr;
182*c4762a1bSJed Brown   PetscScalar    pfive = .50;
183*c4762a1bSJed Brown   ierr = VecSet(x,pfive);CHKERRQ(ierr);
184*c4762a1bSJed Brown   return 0;
185*c4762a1bSJed Brown }
186*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
187*c4762a1bSJed Brown /*
188*c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
189*c4762a1bSJed Brown 
190*c4762a1bSJed Brown    Input Parameters:
191*c4762a1bSJed Brown .  snes - the SNES context
192*c4762a1bSJed Brown .  x - input vector
193*c4762a1bSJed Brown .  ctx - optional user-defined context, as set by SNESSetFunction()
194*c4762a1bSJed Brown 
195*c4762a1bSJed Brown    Output Parameter:
196*c4762a1bSJed Brown .  f - function vector
197*c4762a1bSJed Brown 
198*c4762a1bSJed Brown    Note:
199*c4762a1bSJed Brown    The user-defined context can contain any application-specific data
200*c4762a1bSJed Brown    needed for the function evaluation (such as various parameters, work
201*c4762a1bSJed Brown    vectors, and grid information).  In this program the context is just
202*c4762a1bSJed Brown    a vector containing the right-hand-side of the discretized PDE.
203*c4762a1bSJed Brown  */
204*c4762a1bSJed Brown 
205*c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx)
206*c4762a1bSJed Brown {
207*c4762a1bSJed Brown   Vec               g = (Vec)ctx;
208*c4762a1bSJed Brown   const PetscScalar *xx,*gg;
209*c4762a1bSJed Brown   PetscScalar       *ff,d;
210*c4762a1bSJed Brown   PetscErrorCode    ierr;
211*c4762a1bSJed Brown   PetscInt          i,n;
212*c4762a1bSJed Brown 
213*c4762a1bSJed Brown   /*
214*c4762a1bSJed Brown      Get pointers to vector data.
215*c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
216*c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
217*c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
218*c4762a1bSJed Brown          the array.
219*c4762a1bSJed Brown   */
220*c4762a1bSJed Brown   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
221*c4762a1bSJed Brown   ierr = VecGetArray(f,&ff);CHKERRQ(ierr);
222*c4762a1bSJed Brown   ierr = VecGetArrayRead(g,&gg);CHKERRQ(ierr);
223*c4762a1bSJed Brown 
224*c4762a1bSJed Brown   /*
225*c4762a1bSJed Brown      Compute function
226*c4762a1bSJed Brown   */
227*c4762a1bSJed Brown   ierr  = VecGetSize(x,&n);CHKERRQ(ierr);
228*c4762a1bSJed Brown   d     = (PetscReal)(n - 1); d = d*d;
229*c4762a1bSJed Brown   ff[0] = xx[0];
230*c4762a1bSJed 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];
231*c4762a1bSJed Brown   ff[n-1] = xx[n-1] - 1.0;
232*c4762a1bSJed Brown 
233*c4762a1bSJed Brown   /*
234*c4762a1bSJed Brown      Restore vectors
235*c4762a1bSJed Brown   */
236*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
237*c4762a1bSJed Brown   ierr = VecRestoreArray(f,&ff);CHKERRQ(ierr);
238*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(g,&gg);CHKERRQ(ierr);
239*c4762a1bSJed Brown   return 0;
240*c4762a1bSJed Brown }
241*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
242*c4762a1bSJed Brown /*
243*c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
244*c4762a1bSJed Brown 
245*c4762a1bSJed Brown    Input Parameters:
246*c4762a1bSJed Brown .  snes - the SNES context
247*c4762a1bSJed Brown .  x - input vector
248*c4762a1bSJed Brown .  dummy - optional user-defined context (not used here)
249*c4762a1bSJed Brown 
250*c4762a1bSJed Brown    Output Parameters:
251*c4762a1bSJed Brown .  jac - Jacobian matrix
252*c4762a1bSJed Brown .  B - optionally different preconditioning matrix
253*c4762a1bSJed Brown 
254*c4762a1bSJed Brown */
255*c4762a1bSJed Brown 
256*c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *dummy)
257*c4762a1bSJed Brown {
258*c4762a1bSJed Brown   const PetscScalar *xx;
259*c4762a1bSJed Brown   PetscScalar       A[3],d;
260*c4762a1bSJed Brown   PetscErrorCode    ierr;
261*c4762a1bSJed Brown   PetscInt          i,n,j[3];
262*c4762a1bSJed Brown 
263*c4762a1bSJed Brown   /*
264*c4762a1bSJed Brown      Get pointer to vector data
265*c4762a1bSJed Brown   */
266*c4762a1bSJed Brown   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
267*c4762a1bSJed Brown 
268*c4762a1bSJed Brown   /*
269*c4762a1bSJed Brown      Compute Jacobian entries and insert into matrix.
270*c4762a1bSJed Brown       - Note that in this case we set all elements for a particular
271*c4762a1bSJed Brown         row at once.
272*c4762a1bSJed Brown   */
273*c4762a1bSJed Brown   ierr = VecGetSize(x,&n);CHKERRQ(ierr);
274*c4762a1bSJed Brown   d    = (PetscReal)(n - 1); d = d*d;
275*c4762a1bSJed Brown 
276*c4762a1bSJed Brown   /*
277*c4762a1bSJed Brown      Interior grid points
278*c4762a1bSJed Brown   */
279*c4762a1bSJed Brown   for (i=1; i<n-1; i++) {
280*c4762a1bSJed Brown     j[0] = i - 1; j[1] = i; j[2] = i + 1;
281*c4762a1bSJed Brown     A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i];
282*c4762a1bSJed Brown     ierr = MatSetValues(B,1,&i,3,j,A,INSERT_VALUES);CHKERRQ(ierr);
283*c4762a1bSJed Brown   }
284*c4762a1bSJed Brown 
285*c4762a1bSJed Brown   /*
286*c4762a1bSJed Brown      Boundary points
287*c4762a1bSJed Brown   */
288*c4762a1bSJed Brown   i = 0;   A[0] = 1.0;
289*c4762a1bSJed Brown 
290*c4762a1bSJed Brown   ierr = MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES);CHKERRQ(ierr);
291*c4762a1bSJed Brown 
292*c4762a1bSJed Brown   i = n-1; A[0] = 1.0;
293*c4762a1bSJed Brown 
294*c4762a1bSJed Brown   ierr = MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES);CHKERRQ(ierr);
295*c4762a1bSJed Brown 
296*c4762a1bSJed Brown   /*
297*c4762a1bSJed Brown      Restore vector
298*c4762a1bSJed Brown   */
299*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
300*c4762a1bSJed Brown 
301*c4762a1bSJed Brown   /*
302*c4762a1bSJed Brown      Assemble matrix
303*c4762a1bSJed Brown   */
304*c4762a1bSJed Brown   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
305*c4762a1bSJed Brown   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
306*c4762a1bSJed Brown   if (jac != B) {
307*c4762a1bSJed Brown     ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
308*c4762a1bSJed Brown     ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
309*c4762a1bSJed Brown   }
310*c4762a1bSJed Brown   return 0;
311*c4762a1bSJed Brown }
312*c4762a1bSJed Brown 
313*c4762a1bSJed Brown /*TEST
314*c4762a1bSJed Brown 
315*c4762a1bSJed Brown    test:
316*c4762a1bSJed Brown       args: -snes_monitor_short -snes_view -ksp_monitor
317*c4762a1bSJed Brown 
318*c4762a1bSJed Brown    # this is just a test for SNESKSPTRASPOSEONLY and KSPSolveTranspose to behave properly
319*c4762a1bSJed Brown    # the solution is wrong on purpose
320*c4762a1bSJed Brown    test:
321*c4762a1bSJed Brown       requires: !single !complex
322*c4762a1bSJed Brown       suffix: transpose_only
323*c4762a1bSJed 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
324*c4762a1bSJed Brown 
325*c4762a1bSJed Brown TEST*/
326