xref: /petsc/src/snes/tutorials/ex19.c (revision 0faed27c4f59cfdcd00951238dac7483cf8e56f1)
1 
2 static char help[] = "Nonlinear driven cavity with multigrid in 2d.\n \
3   \n\
4 The 2D driven cavity problem is solved in a velocity-vorticity formulation.\n\
5 The flow can be driven with the lid or with bouyancy or both:\n\
6   -lidvelocity &ltlid&gt, where &ltlid&gt = dimensionless velocity of lid\n\
7   -grashof &ltgr&gt, where &ltgr&gt = dimensionless temperature gradent\n\
8   -prandtl &ltpr&gt, where &ltpr&gt = dimensionless thermal/momentum diffusity ratio\n\
9  -contours : draw contour plots of solution\n\n";
10 /* in HTML, '&lt' = '<' and '&gt' = '>' */
11 
12 /*
13       See src/ksp/ksp/tutorials/ex45.c
14 */
15 
16 /*F-----------------------------------------------------------------------
17 
18     We thank David E. Keyes for contributing the driven cavity discretization within this example code.
19 
20     This problem is modeled by the partial differential equation system
21 
22 \begin{eqnarray}
23         - \triangle U - \nabla_y \Omega & = & 0  \\
24         - \triangle V + \nabla_x\Omega & = & 0  \\
25         - \triangle \Omega + \nabla \cdot ([U*\Omega,V*\Omega]) - GR* \nabla_x T & = & 0  \\
26         - \triangle T + PR* \nabla \cdot ([U*T,V*T]) & = & 0
27 \end{eqnarray}
28 
29     in the unit square, which is uniformly discretized in each of x and y in this simple encoding.
30 
31     No-slip, rigid-wall Dirichlet conditions are used for $ [U,V]$.
32     Dirichlet conditions are used for Omega, based on the definition of
33     vorticity: $ \Omega = - \nabla_y U + \nabla_x V$, where along each
34     constant coordinate boundary, the tangential derivative is zero.
35     Dirichlet conditions are used for T on the left and right walls,
36     and insulation homogeneous Neumann conditions are used for T on
37     the top and bottom walls.
38 
39     A finite difference approximation with the usual 5-point stencil
40     is used to discretize the boundary value problem to obtain a
41     nonlinear system of equations.  Upwinding is used for the divergence
42     (convective) terms and central for the gradient (source) terms.
43 
44     The Jacobian can be either
45       * formed via finite differencing using coloring (the default), or
46       * applied matrix-free via the option -snes_mf
47         (for larger grid problems this variant may not converge
48         without a preconditioner due to ill-conditioning).
49 
50   ------------------------------------------------------------------------F*/
51 
52 /*
53    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
54    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
55    file automatically includes:
56      petscsys.h       - base PETSc routines   petscvec.h - vectors
57      petscmat.h - matrices
58      petscis.h     - index sets            petscksp.h - Krylov subspace methods
59      petscviewer.h - viewers               petscpc.h  - preconditioners
60      petscksp.h   - linear solvers
61 */
62 #if defined(PETSC_APPLE_FRAMEWORK)
63 #import <PETSc/petscsnes.h>
64 #import <PETSc/petscdmda.h>
65 #else
66 #include <petscsnes.h>
67 #include <petscdm.h>
68 #include <petscdmda.h>
69 #endif
70 
71 /*
72    User-defined routines and data structures
73 */
74 typedef struct {
75   PetscScalar u,v,omega,temp;
76 } Field;
77 
78 PetscErrorCode FormFunctionLocal(DMDALocalInfo*,Field**,Field**,void*);
79 
80 typedef struct {
81   PetscReal   lidvelocity,prandtl,grashof;  /* physical parameters */
82   PetscBool   draw_contours;                /* flag - 1 indicates drawing contours */
83 } AppCtx;
84 
85 extern PetscErrorCode FormInitialGuess(AppCtx*,DM,Vec);
86 extern PetscErrorCode NonlinearGS(SNES,Vec,Vec,void*);
87 
88 int main(int argc,char **argv)
89 {
90   AppCtx         user;                /* user-defined work context */
91   PetscInt       mx,my,its;
92   MPI_Comm       comm;
93   SNES           snes;
94   DM             da;
95   Vec            x;
96 
97   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
98   comm = PETSC_COMM_WORLD;
99   PetscCall(SNESCreate(comm,&snes));
100 
101   /*
102       Create distributed array object to manage parallel grid and vectors
103       for principal unknowns (x) and governing residuals (f)
104   */
105   PetscCall(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,4,4,PETSC_DECIDE,PETSC_DECIDE,4,1,0,0,&da));
106   PetscCall(DMSetFromOptions(da));
107   PetscCall(DMSetUp(da));
108   PetscCall(SNESSetDM(snes,(DM)da));
109   PetscCall(SNESSetNGS(snes, NonlinearGS, (void*)&user));
110 
111   PetscCall(DMDAGetInfo(da,0,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE));
112   /*
113      Problem parameters (velocity of lid, prandtl, and grashof numbers)
114   */
115   user.lidvelocity = 1.0/(mx*my);
116   user.prandtl     = 1.0;
117   user.grashof     = 1.0;
118 
119   PetscCall(PetscOptionsGetReal(NULL,NULL,"-lidvelocity",&user.lidvelocity,NULL));
120   PetscCall(PetscOptionsGetReal(NULL,NULL,"-prandtl",&user.prandtl,NULL));
121   PetscCall(PetscOptionsGetReal(NULL,NULL,"-grashof",&user.grashof,NULL));
122   PetscCall(PetscOptionsHasName(NULL,NULL,"-contours",&user.draw_contours));
123 
124   PetscCall(DMDASetFieldName(da,0,"x_velocity"));
125   PetscCall(DMDASetFieldName(da,1,"y_velocity"));
126   PetscCall(DMDASetFieldName(da,2,"Omega"));
127   PetscCall(DMDASetFieldName(da,3,"temperature"));
128 
129   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130      Create user context, set problem data, create vector data structures.
131      Also, compute the initial guess.
132      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133 
134   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135      Create nonlinear solver context
136 
137      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138   PetscCall(DMSetApplicationContext(da,&user));
139   PetscCall(DMDASNESSetFunctionLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))FormFunctionLocal,&user));
140   PetscCall(SNESSetFromOptions(snes));
141   PetscCall(PetscPrintf(comm,"lid velocity = %g, prandtl # = %g, grashof # = %g\n",(double)user.lidvelocity,(double)user.prandtl,(double)user.grashof));
142 
143   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144      Solve the nonlinear system
145      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
146   PetscCall(DMCreateGlobalVector(da,&x));
147   PetscCall(FormInitialGuess(&user,da,x));
148 
149   PetscCall(SNESSolve(snes,NULL,x));
150 
151   PetscCall(SNESGetIterationNumber(snes,&its));
152   PetscCall(PetscPrintf(comm,"Number of SNES iterations = %" PetscInt_FMT "\n", its));
153 
154   /*
155      Visualize solution
156   */
157   if (user.draw_contours) PetscCall(VecView(x,PETSC_VIEWER_DRAW_WORLD));
158 
159   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
160      Free work space.  All PETSc objects should be destroyed when they
161      are no longer needed.
162      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
163   PetscCall(VecDestroy(&x));
164   PetscCall(DMDestroy(&da));
165   PetscCall(SNESDestroy(&snes));
166   PetscCall(PetscFinalize());
167   return 0;
168 }
169 
170 /* ------------------------------------------------------------------- */
171 
172 /*
173    FormInitialGuess - Forms initial approximation.
174 
175    Input Parameters:
176    user - user-defined application context
177    X - vector
178 
179    Output Parameter:
180    X - vector
181 */
182 PetscErrorCode FormInitialGuess(AppCtx *user,DM da,Vec X)
183 {
184   PetscInt       i,j,mx,xs,ys,xm,ym;
185   PetscReal      grashof,dx;
186   Field          **x;
187 
188   PetscFunctionBeginUser;
189   grashof = user->grashof;
190 
191   PetscCall(DMDAGetInfo(da,0,&mx,0,0,0,0,0,0,0,0,0,0,0));
192   dx   = 1.0/(mx-1);
193 
194   /*
195      Get local grid boundaries (for 2-dimensional DMDA):
196        xs, ys   - starting grid indices (no ghost points)
197        xm, ym   - widths of local grid (no ghost points)
198   */
199   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
200 
201   /*
202      Get a pointer to vector data.
203        - For default PETSc vectors, VecGetArray() returns a pointer to
204          the data array.  Otherwise, the routine is implementation dependent.
205        - You MUST call VecRestoreArray() when you no longer need access to
206          the array.
207   */
208   PetscCall(DMDAVecGetArrayWrite(da,X,&x));
209 
210   /*
211      Compute initial guess over the locally owned part of the grid
212      Initial condition is motionless fluid and equilibrium temperature
213   */
214   for (j=ys; j<ys+ym; j++) {
215     for (i=xs; i<xs+xm; i++) {
216       x[j][i].u     = 0.0;
217       x[j][i].v     = 0.0;
218       x[j][i].omega = 0.0;
219       x[j][i].temp  = (grashof>0)*i*dx;
220     }
221   }
222 
223   /*
224      Restore vector
225   */
226   PetscCall(DMDAVecRestoreArrayWrite(da,X,&x));
227   PetscFunctionReturn(0);
228 }
229 
230 PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,Field **x,Field **f,void *ptr)
231 {
232   AppCtx         *user = (AppCtx*)ptr;
233   PetscInt       xints,xinte,yints,yinte,i,j;
234   PetscReal      hx,hy,dhx,dhy,hxdhy,hydhx;
235   PetscReal      grashof,prandtl,lid;
236   PetscScalar    u,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;
237 
238   PetscFunctionBeginUser;
239   grashof = user->grashof;
240   prandtl = user->prandtl;
241   lid     = user->lidvelocity;
242 
243   /*
244      Define mesh intervals ratios for uniform grid.
245 
246      Note: FD formulae below are normalized by multiplying through by
247      local volume element (i.e. hx*hy) to obtain coefficients O(1) in two dimensions.
248 
249   */
250   dhx   = (PetscReal)(info->mx-1);  dhy = (PetscReal)(info->my-1);
251   hx    = 1.0/dhx;                   hy = 1.0/dhy;
252   hxdhy = hx*dhy;                 hydhx = hy*dhx;
253 
254   xints = info->xs; xinte = info->xs+info->xm; yints = info->ys; yinte = info->ys+info->ym;
255 
256   /* Test whether we are on the bottom edge of the global array */
257   if (yints == 0) {
258     j     = 0;
259     yints = yints + 1;
260     /* bottom edge */
261     for (i=info->xs; i<info->xs+info->xm; i++) {
262       f[j][i].u     = x[j][i].u;
263       f[j][i].v     = x[j][i].v;
264       f[j][i].omega = x[j][i].omega + (x[j+1][i].u - x[j][i].u)*dhy;
265       f[j][i].temp  = x[j][i].temp-x[j+1][i].temp;
266     }
267   }
268 
269   /* Test whether we are on the top edge of the global array */
270   if (yinte == info->my) {
271     j     = info->my - 1;
272     yinte = yinte - 1;
273     /* top edge */
274     for (i=info->xs; i<info->xs+info->xm; i++) {
275       f[j][i].u     = x[j][i].u - lid;
276       f[j][i].v     = x[j][i].v;
277       f[j][i].omega = x[j][i].omega + (x[j][i].u - x[j-1][i].u)*dhy;
278       f[j][i].temp  = x[j][i].temp-x[j-1][i].temp;
279     }
280   }
281 
282   /* Test whether we are on the left edge of the global array */
283   if (xints == 0) {
284     i     = 0;
285     xints = xints + 1;
286     /* left edge */
287     for (j=info->ys; j<info->ys+info->ym; j++) {
288       f[j][i].u     = x[j][i].u;
289       f[j][i].v     = x[j][i].v;
290       f[j][i].omega = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx;
291       f[j][i].temp  = x[j][i].temp;
292     }
293   }
294 
295   /* Test whether we are on the right edge of the global array */
296   if (xinte == info->mx) {
297     i     = info->mx - 1;
298     xinte = xinte - 1;
299     /* right edge */
300     for (j=info->ys; j<info->ys+info->ym; j++) {
301       f[j][i].u     = x[j][i].u;
302       f[j][i].v     = x[j][i].v;
303       f[j][i].omega = x[j][i].omega - (x[j][i].v - x[j][i-1].v)*dhx;
304       f[j][i].temp  = x[j][i].temp - (PetscReal)(grashof>0);
305     }
306   }
307 
308   /* Compute over the interior points */
309   for (j=yints; j<yinte; j++) {
310     for (i=xints; i<xinte; i++) {
311 
312       /*
313        convective coefficients for upwinding
314       */
315       vx  = x[j][i].u; avx = PetscAbsScalar(vx);
316       vxp = .5*(vx+avx); vxm = .5*(vx-avx);
317       vy  = x[j][i].v; avy = PetscAbsScalar(vy);
318       vyp = .5*(vy+avy); vym = .5*(vy-avy);
319 
320       /* U velocity */
321       u         = x[j][i].u;
322       uxx       = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
323       uyy       = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
324       f[j][i].u = uxx + uyy - .5*(x[j+1][i].omega-x[j-1][i].omega)*hx;
325 
326       /* V velocity */
327       u         = x[j][i].v;
328       uxx       = (2.0*u - x[j][i-1].v - x[j][i+1].v)*hydhx;
329       uyy       = (2.0*u - x[j-1][i].v - x[j+1][i].v)*hxdhy;
330       f[j][i].v = uxx + uyy + .5*(x[j][i+1].omega-x[j][i-1].omega)*hy;
331 
332       /* Omega */
333       u             = x[j][i].omega;
334       uxx           = (2.0*u - x[j][i-1].omega - x[j][i+1].omega)*hydhx;
335       uyy           = (2.0*u - x[j-1][i].omega - x[j+1][i].omega)*hxdhy;
336       f[j][i].omega = uxx + uyy + (vxp*(u - x[j][i-1].omega) + vxm*(x[j][i+1].omega - u))*hy +
337                       (vyp*(u - x[j-1][i].omega) + vym*(x[j+1][i].omega - u))*hx -
338                       .5*grashof*(x[j][i+1].temp - x[j][i-1].temp)*hy;
339 
340       /* Temperature */
341       u            = x[j][i].temp;
342       uxx          = (2.0*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx;
343       uyy          = (2.0*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy;
344       f[j][i].temp =  uxx + uyy  + prandtl*((vxp*(u - x[j][i-1].temp) + vxm*(x[j][i+1].temp - u))*hy +
345                                             (vyp*(u - x[j-1][i].temp) + vym*(x[j+1][i].temp - u))*hx);
346     }
347   }
348 
349   /*
350      Flop count (multiply-adds are counted as 2 operations)
351   */
352   PetscCall(PetscLogFlops(84.0*info->ym*info->xm));
353   PetscFunctionReturn(0);
354 }
355 
356 /*
357     Performs sweeps of point block nonlinear Gauss-Seidel on all the local grid points
358 */
359 PetscErrorCode NonlinearGS(SNES snes, Vec X, Vec B, void *ctx)
360 {
361   DMDALocalInfo  info;
362   Field          **x,**b;
363   Vec            localX, localB;
364   DM             da;
365   PetscInt       xints,xinte,yints,yinte,i,j,k,l;
366   PetscInt       max_its,tot_its;
367   PetscInt       sweeps;
368   PetscReal      rtol,atol,stol;
369   PetscReal      hx,hy,dhx,dhy,hxdhy,hydhx;
370   PetscReal      grashof,prandtl,lid;
371   PetscScalar    u,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;
372   PetscScalar    fu, fv, fomega, ftemp;
373   PetscScalar    dfudu;
374   PetscScalar    dfvdv;
375   PetscScalar    dfodu, dfodv, dfodo;
376   PetscScalar    dftdu, dftdv, dftdt;
377   PetscScalar    yu=0, yv=0, yo=0, yt=0;
378   PetscScalar    bjiu, bjiv, bjiomega, bjitemp;
379   PetscBool      ptconverged;
380   PetscReal      pfnorm,pfnorm0,pynorm,pxnorm;
381   AppCtx         *user = (AppCtx*)ctx;
382 
383   PetscFunctionBeginUser;
384   grashof = user->grashof;
385   prandtl = user->prandtl;
386   lid     = user->lidvelocity;
387   tot_its = 0;
388   PetscCall(SNESNGSGetTolerances(snes,&rtol,&atol,&stol,&max_its));
389   PetscCall(SNESNGSGetSweeps(snes,&sweeps));
390   PetscCall(SNESGetDM(snes,(DM*)&da));
391   PetscCall(DMGetLocalVector(da,&localX));
392   if (B) {
393     PetscCall(DMGetLocalVector(da,&localB));
394   }
395   /*
396      Scatter ghost points to local vector, using the 2-step process
397         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
398   */
399   PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
400   PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
401   if (B) {
402     PetscCall(DMGlobalToLocalBegin(da,B,INSERT_VALUES,localB));
403     PetscCall(DMGlobalToLocalEnd(da,B,INSERT_VALUES,localB));
404   }
405   PetscCall(DMDAGetLocalInfo(da,&info));
406   PetscCall(DMDAVecGetArrayWrite(da,localX,&x));
407   if (B) {
408     PetscCall(DMDAVecGetArrayRead(da,localB,&b));
409   }
410   /* looks like a combination of the formfunction / formjacobian routines */
411   dhx   = (PetscReal)(info.mx-1);dhy   = (PetscReal)(info.my-1);
412   hx    = 1.0/dhx;               hy    = 1.0/dhy;
413   hxdhy = hx*dhy;                hydhx = hy*dhx;
414 
415   xints = info.xs; xinte = info.xs+info.xm; yints = info.ys; yinte = info.ys+info.ym;
416 
417   /* Set the boundary conditions on the momentum equations */
418   /* Test whether we are on the bottom edge of the global array */
419   if (yints == 0) {
420     j     = 0;
421     /* bottom edge */
422     for (i=info.xs; i<info.xs+info.xm; i++) {
423 
424       if (B) {
425         bjiu = b[j][i].u;
426         bjiv = b[j][i].v;
427       } else {
428         bjiu = 0.0;
429         bjiv = 0.0;
430       }
431       x[j][i].u = 0.0 + bjiu;
432       x[j][i].v = 0.0 + bjiv;
433     }
434   }
435 
436   /* Test whether we are on the top edge of the global array */
437   if (yinte == info.my) {
438     j     = info.my - 1;
439     /* top edge */
440     for (i=info.xs; i<info.xs+info.xm; i++) {
441       if (B) {
442         bjiu = b[j][i].u;
443         bjiv = b[j][i].v;
444       } else {
445         bjiu = 0.0;
446         bjiv = 0.0;
447       }
448       x[j][i].u = lid + bjiu;
449       x[j][i].v = bjiv;
450     }
451   }
452 
453   /* Test whether we are on the left edge of the global array */
454   if (xints == 0) {
455     i     = 0;
456     /* left edge */
457     for (j=info.ys; j<info.ys+info.ym; j++) {
458       if (B) {
459         bjiu = b[j][i].u;
460         bjiv = b[j][i].v;
461       } else {
462         bjiu = 0.0;
463         bjiv = 0.0;
464       }
465       x[j][i].u = 0.0 + bjiu;
466       x[j][i].v = 0.0 + bjiv;
467     }
468   }
469 
470   /* Test whether we are on the right edge of the global array */
471   if (xinte == info.mx) {
472     i     = info.mx - 1;
473     /* right edge */
474     for (j=info.ys; j<info.ys+info.ym; j++) {
475       if (B) {
476         bjiu = b[j][i].u;
477         bjiv = b[j][i].v;
478       } else {
479         bjiu = 0.0;
480         bjiv = 0.0;
481       }
482       x[j][i].u = 0.0 + bjiu;
483       x[j][i].v = 0.0 + bjiv;
484     }
485   }
486 
487   for (k=0; k < sweeps; k++) {
488     for (j=info.ys; j<info.ys + info.ym; j++) {
489       for (i=info.xs; i<info.xs + info.xm; i++) {
490         ptconverged = PETSC_FALSE;
491         pfnorm0     = 0.0;
492         fu          = 0.0;
493         fv          = 0.0;
494         fomega      = 0.0;
495         ftemp       = 0.0;
496         /*  Run Newton's method on a single grid point */
497         for (l = 0; l < max_its && !ptconverged; l++) {
498           if (B) {
499             bjiu     = b[j][i].u;
500             bjiv     = b[j][i].v;
501             bjiomega = b[j][i].omega;
502             bjitemp  = b[j][i].temp;
503           } else {
504             bjiu     = 0.0;
505             bjiv     = 0.0;
506             bjiomega = 0.0;
507             bjitemp  = 0.0;
508           }
509 
510           if (i != 0 && i != info.mx - 1 && j != 0 && j != info.my-1) {
511             /* U velocity */
512             u     = x[j][i].u;
513             uxx   = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
514             uyy   = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
515             fu    = uxx + uyy - .5*(x[j+1][i].omega-x[j-1][i].omega)*hx - bjiu;
516             dfudu = 2.0*(hydhx + hxdhy);
517             /* V velocity */
518             u     = x[j][i].v;
519             uxx   = (2.0*u - x[j][i-1].v - x[j][i+1].v)*hydhx;
520             uyy   = (2.0*u - x[j-1][i].v - x[j+1][i].v)*hxdhy;
521             fv    = uxx + uyy + .5*(x[j][i+1].omega-x[j][i-1].omega)*hy - bjiv;
522             dfvdv = 2.0*(hydhx + hxdhy);
523             /*
524              convective coefficients for upwinding
525              */
526             vx  = x[j][i].u; avx = PetscAbsScalar(vx);
527             vxp = .5*(vx+avx); vxm = .5*(vx-avx);
528             vy  = x[j][i].v; avy = PetscAbsScalar(vy);
529             vyp = .5*(vy+avy); vym = .5*(vy-avy);
530             /* Omega */
531             u      = x[j][i].omega;
532             uxx    = (2.0*u - x[j][i-1].omega - x[j][i+1].omega)*hydhx;
533             uyy    = (2.0*u - x[j-1][i].omega - x[j+1][i].omega)*hxdhy;
534             fomega = uxx + uyy +  (vxp*(u - x[j][i-1].omega) + vxm*(x[j][i+1].omega - u))*hy +
535                      (vyp*(u - x[j-1][i].omega) + vym*(x[j+1][i].omega - u))*hx -
536                      .5*grashof*(x[j][i+1].temp - x[j][i-1].temp)*hy - bjiomega;
537             /* convective coefficient derivatives */
538             dfodo = 2.0*(hydhx + hxdhy) + ((vxp - vxm)*hy + (vyp - vym)*hx);
539             if (PetscRealPart(vx) > 0.0) dfodu = (u - x[j][i-1].omega)*hy;
540             else dfodu = (x[j][i+1].omega - u)*hy;
541 
542             if (PetscRealPart(vy) > 0.0) dfodv = (u - x[j-1][i].omega)*hx;
543             else dfodv = (x[j+1][i].omega - u)*hx;
544 
545             /* Temperature */
546             u     = x[j][i].temp;
547             uxx   = (2.0*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx;
548             uyy   = (2.0*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy;
549             ftemp =  uxx + uyy  + prandtl*((vxp*(u - x[j][i-1].temp) + vxm*(x[j][i+1].temp - u))*hy + (vyp*(u - x[j-1][i].temp) + vym*(x[j+1][i].temp - u))*hx) - bjitemp;
550             dftdt = 2.0*(hydhx + hxdhy) + prandtl*((vxp - vxm)*hy + (vyp - vym)*hx);
551             if (PetscRealPart(vx) > 0.0) dftdu = prandtl*(u - x[j][i-1].temp)*hy;
552             else dftdu = prandtl*(x[j][i+1].temp - u)*hy;
553 
554             if (PetscRealPart(vy) > 0.0) dftdv = prandtl*(u - x[j-1][i].temp)*hx;
555             else dftdv = prandtl*(x[j+1][i].temp - u)*hx;
556 
557             /* invert the system:
558              [ dfu / du     0        0        0    ][yu] = [fu]
559              [     0    dfv / dv     0        0    ][yv]   [fv]
560              [ dfo / du dfo / dv dfo / do     0    ][yo]   [fo]
561              [ dft / du dft / dv     0    dft / dt ][yt]   [ft]
562              by simple back-substitution
563            */
564             yu = fu / dfudu;
565             yv = fv / dfvdv;
566             yo = (fomega - (dfodu*yu + dfodv*yv)) / dfodo;
567             yt = (ftemp - (dftdu*yu + dftdv*yv)) / dftdt;
568 
569             x[j][i].u     = x[j][i].u - yu;
570             x[j][i].v     = x[j][i].v - yv;
571             x[j][i].temp  = x[j][i].temp - yt;
572             x[j][i].omega = x[j][i].omega - yo;
573           }
574           if (i == 0) {
575             fomega        = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx - bjiomega;
576             ftemp         = x[j][i].temp - bjitemp;
577             yo            = fomega;
578             yt            = ftemp;
579             x[j][i].omega = x[j][i].omega - fomega;
580             x[j][i].temp  = x[j][i].temp - ftemp;
581           }
582           if (i == info.mx - 1) {
583             fomega        = x[j][i].omega - (x[j][i].v - x[j][i-1].v)*dhx - bjiomega;
584             ftemp         = x[j][i].temp - (PetscReal)(grashof>0) - bjitemp;
585             yo            = fomega;
586             yt            = ftemp;
587             x[j][i].omega = x[j][i].omega - fomega;
588             x[j][i].temp  = x[j][i].temp - ftemp;
589           }
590           if (j == 0) {
591             fomega        = x[j][i].omega + (x[j+1][i].u - x[j][i].u)*dhy - bjiomega;
592             ftemp         = x[j][i].temp-x[j+1][i].temp - bjitemp;
593             yo            = fomega;
594             yt            = ftemp;
595             x[j][i].omega = x[j][i].omega - fomega;
596             x[j][i].temp  = x[j][i].temp - ftemp;
597           }
598           if (j == info.my - 1) {
599             fomega        = x[j][i].omega + (x[j][i].u - x[j-1][i].u)*dhy - bjiomega;
600             ftemp         = x[j][i].temp-x[j-1][i].temp - bjitemp;
601             yo            = fomega;
602             yt            = ftemp;
603             x[j][i].omega = x[j][i].omega - fomega;
604             x[j][i].temp  = x[j][i].temp - ftemp;
605           }
606           tot_its++;
607           pfnorm = PetscRealPart(fu*fu + fv*fv + fomega*fomega + ftemp*ftemp);
608           pfnorm = PetscSqrtReal(pfnorm);
609           pynorm = PetscRealPart(yu*yu + yv*yv + yo*yo + yt*yt);
610           pynorm = PetscSqrtReal(pynorm);
611           pxnorm = PetscRealPart(x[j][i].u*x[j][i].u + x[j][i].v*x[j][i].v + x[j][i].omega*x[j][i].omega + x[j][i].temp*x[j][i].temp);
612           pxnorm = PetscSqrtReal(pxnorm);
613           if (l == 0) pfnorm0 = pfnorm;
614           if (rtol*pfnorm0 >pfnorm || atol > pfnorm || pxnorm*stol > pynorm) ptconverged = PETSC_TRUE;
615         }
616       }
617     }
618   }
619   PetscCall(DMDAVecRestoreArrayWrite(da,localX,&x));
620   if (B) {
621     PetscCall(DMDAVecRestoreArrayRead(da,localB,&b));
622   }
623   PetscCall(DMLocalToGlobalBegin(da,localX,INSERT_VALUES,X));
624   PetscCall(DMLocalToGlobalEnd(da,localX,INSERT_VALUES,X));
625   PetscCall(PetscLogFlops(tot_its*(84.0 + 41.0 + 26.0)));
626   PetscCall(DMRestoreLocalVector(da,&localX));
627   if (B) {
628     PetscCall(DMRestoreLocalVector(da,&localB));
629   }
630   PetscFunctionReturn(0);
631 }
632 
633 /*TEST
634 
635    test:
636       nsize: 2
637       args: -da_refine 3 -snes_monitor_short -pc_type mg -ksp_type fgmres -pc_mg_type full
638       requires: !single
639 
640    test:
641       suffix: 10
642       nsize: 3
643       args: -snes_monitor_short -ksp_monitor_short -pc_type fieldsplit -pc_fieldsplit_type symmetric_multiplicative -snes_view -da_refine 1 -ksp_type fgmres
644       requires: !single
645 
646    test:
647       suffix: 11
648       nsize: 4
649       requires: pastix
650       args: -snes_monitor_short -pc_type redundant -dm_mat_type mpiaij -redundant_pc_factor_mat_solver_type pastix -pc_redundant_number 2 -da_refine 4 -ksp_type fgmres
651 
652    test:
653       suffix: 12
654       nsize: 12
655       requires: pastix
656       args: -snes_monitor_short -pc_type redundant -dm_mat_type mpiaij -redundant_pc_factor_mat_solver_type pastix -pc_redundant_number 5 -da_refine 4 -ksp_type fgmres
657 
658    test:
659       suffix: 13
660       nsize: 3
661       args: -snes_monitor_short -ksp_monitor_short -pc_type fieldsplit -pc_fieldsplit_type multiplicative -snes_view -da_refine 1 -ksp_type fgmres -snes_mf_operator
662       requires: !single
663 
664    test:
665       suffix: 14
666       nsize: 4
667       args: -snes_monitor_short -pc_type mg -dm_mat_type baij -mg_coarse_pc_type bjacobi -da_refine 3 -ksp_type fgmres
668       requires: !single
669 
670    test:
671       suffix: 14_ds
672       nsize: 4
673       args: -snes_converged_reason -pc_type mg -dm_mat_type baij -mg_coarse_pc_type bjacobi -da_refine 3 -ksp_type fgmres -mat_fd_type ds
674       output_file: output/ex19_2.out
675       requires: !single
676 
677    test:
678       suffix: 17
679       args: -snes_monitor_short -ksp_pc_side right
680       requires: !single
681 
682    test:
683       suffix: 18
684       args: -snes_monitor_ksp draw::draw_lg -ksp_pc_side right
685       requires: x !single
686 
687    test:
688       suffix: 19
689       nsize: 2
690       args: -da_refine 3 -snes_monitor_short -pc_type mg -ksp_type fgmres -pc_mg_type full -snes_type newtontrdc
691       requires: !single
692 
693    test:
694       suffix: 20
695       nsize: 2
696       args: -da_refine 3 -snes_monitor_short -pc_type mg -ksp_type fgmres -pc_mg_type full -snes_type newtontrdc -snes_trdc_use_cauchy false
697       requires: !single
698 
699    test:
700       suffix: 2
701       nsize: 4
702       args: -da_refine 3 -snes_converged_reason -pc_type mg -mat_fd_type ds
703       requires: !single
704 
705    test:
706       suffix: 2_bcols1
707       nsize: 4
708       args: -da_refine 3 -snes_converged_reason -pc_type mg -mat_fd_type ds -mat_fd_coloring_bcols
709       output_file: output/ex19_2.out
710       requires: !single
711 
712    test:
713       suffix: 3
714       nsize: 4
715       requires: mumps
716       args: -da_refine 3 -snes_monitor_short -pc_type redundant -dm_mat_type mpiaij -redundant_ksp_type preonly -redundant_pc_factor_mat_solver_type mumps -pc_redundant_number 2
717 
718    test:
719       suffix: 4
720       nsize: 12
721       requires: mumps
722       args: -da_refine 3 -snes_monitor_short -pc_type redundant -dm_mat_type mpiaij -redundant_ksp_type preonly -redundant_pc_factor_mat_solver_type mumps -pc_redundant_number 5
723       output_file: output/ex19_3.out
724 
725    test:
726       suffix: 6
727       args: -snes_monitor_short -ksp_monitor_short -pc_type fieldsplit -snes_view -ksp_type fgmres -da_refine 1
728       requires: !single
729 
730    test:
731       suffix: 7
732       nsize: 3
733       args: -snes_monitor_short -ksp_monitor_short -pc_type fieldsplit -snes_view -da_refine 1 -ksp_type fgmres
734 
735       requires: !single
736    test:
737       suffix: 8
738       args: -snes_monitor_short -ksp_monitor_short -pc_type fieldsplit -pc_fieldsplit_block_size 2 -pc_fieldsplit_0_fields 0,1 -pc_fieldsplit_1_fields 0,1 -pc_fieldsplit_type multiplicative -snes_view -fieldsplit_pc_type lu -da_refine 1 -ksp_type fgmres
739       requires: !single
740 
741    test:
742       suffix: 9
743       nsize: 3
744       args: -snes_monitor_short -ksp_monitor_short -pc_type fieldsplit -pc_fieldsplit_type multiplicative -snes_view -da_refine 1 -ksp_type fgmres
745       requires: !single
746 
747    test:
748       suffix: aspin
749       nsize: 4
750       args: -da_refine 3 -da_overlap 2 -snes_monitor_short -snes_type aspin -grashof 4e4 -lidvelocity 100 -ksp_monitor_short
751       requires: !single
752 
753    test:
754       suffix: bcgsl
755       nsize: 2
756       args: -ksp_type bcgsl -ksp_monitor_short -da_refine 2 -ksp_bcgsl_ell 3 -snes_view
757       requires: !single
758 
759    test:
760       suffix: bcols1
761       nsize: 2
762       args: -da_refine 3 -snes_monitor_short -pc_type mg -ksp_type fgmres -pc_mg_type full -mat_fd_coloring_bcols 1
763       output_file: output/ex19_1.out
764       requires: !single
765 
766    test:
767       suffix: bjacobi
768       nsize: 4
769       args: -da_refine 4 -ksp_type fgmres -pc_type bjacobi -pc_bjacobi_blocks 2 -sub_ksp_type gmres -sub_ksp_max_it 2 -sub_pc_type bjacobi -sub_sub_ksp_type preonly -sub_sub_pc_type ilu -snes_monitor_short
770       requires: !single
771 
772    test:
773       suffix: cgne
774       args: -da_refine 2 -pc_type lu -ksp_type cgne -ksp_monitor_short -ksp_converged_reason -ksp_view -ksp_norm_type unpreconditioned
775       filter: grep -v HERMITIAN
776       requires: !single
777 
778    test:
779       suffix: cgs
780       args: -da_refine 1 -ksp_monitor_short -ksp_type cgs
781       requires: !single
782 
783    test:
784       suffix: composite_fieldsplit
785       args: -ksp_type fgmres -pc_type composite -pc_composite_type MULTIPLICATIVE -pc_composite_pcs fieldsplit,none -sub_0_pc_fieldsplit_block_size 4 -sub_0_pc_fieldsplit_type additive -sub_0_pc_fieldsplit_0_fields 0,1,2 -sub_0_pc_fieldsplit_1_fields 3 -snes_monitor_short -ksp_monitor_short
786       requires: !single
787 
788    test:
789       suffix: composite_fieldsplit_bjacobi
790       args: -ksp_type fgmres -pc_type composite -pc_composite_type MULTIPLICATIVE -pc_composite_pcs fieldsplit,bjacobi -sub_0_pc_fieldsplit_block_size 4 -sub_0_pc_fieldsplit_type additive -sub_0_pc_fieldsplit_0_fields 0,1,2 -sub_0_pc_fieldsplit_1_fields 3 -sub_1_pc_bjacobi_blocks 16 -sub_1_sub_pc_type lu -snes_monitor_short -ksp_monitor_short
791       requires: !single
792 
793    test:
794       suffix: composite_fieldsplit_bjacobi_2
795       nsize: 4
796       args: -ksp_type fgmres -pc_type composite -pc_composite_type MULTIPLICATIVE -pc_composite_pcs fieldsplit,bjacobi -sub_0_pc_fieldsplit_block_size 4 -sub_0_pc_fieldsplit_type additive -sub_0_pc_fieldsplit_0_fields 0,1,2 -sub_0_pc_fieldsplit_1_fields 3 -sub_1_pc_bjacobi_blocks 16 -sub_1_sub_pc_type lu -snes_monitor_short -ksp_monitor_short
797       requires: !single
798 
799    test:
800       suffix: composite_gs_newton
801       nsize: 2
802       args: -da_refine 3 -grashof 4e4 -lidvelocity 100 -snes_monitor_short -snes_type composite -snes_composite_type additiveoptimal -snes_composite_sneses ngs,newtonls -sub_0_snes_max_it 20 -sub_1_pc_type mg
803       requires: !single
804 
805    test:
806       suffix: cuda
807       requires: cuda !single
808       args: -dm_vec_type cuda -dm_mat_type aijcusparse -pc_type none -ksp_type fgmres -snes_monitor_short -snes_rtol 1.e-5
809 
810    test:
811       suffix: draw
812       args: -pc_type fieldsplit -snes_view draw -fieldsplit_x_velocity_pc_type mg -fieldsplit_x_velocity_pc_mg_galerkin pmat -fieldsplit_x_velocity_pc_mg_levels 2 -da_refine 1 -fieldsplit_x_velocity_mg_coarse_pc_type svd
813       requires: x !single
814 
815    test:
816       suffix: drawports
817       args: -snes_monitor_solution draw::draw_ports -da_refine 1
818       output_file: output/ex19_draw.out
819       requires: x !single
820 
821    test:
822       suffix: fas
823       args: -da_refine 4 -snes_monitor_short -snes_type fas -fas_levels_snes_type ngs -fas_levels_snes_ngs_sweeps 3 -fas_levels_snes_ngs_atol 0.0 -fas_levels_snes_ngs_stol 0.0 -grashof 4e4 -snes_fas_smoothup 6 -snes_fas_smoothdown 6 -lidvelocity 100
824       requires: !single
825 
826    test:
827       suffix: fas_full
828       args: -da_refine 4 -snes_monitor_short -snes_type fas -snes_fas_type full -snes_fas_full_downsweep -fas_levels_snes_type ngs -fas_levels_snes_ngs_sweeps 3 -fas_levels_snes_ngs_atol 0.0 -fas_levels_snes_ngs_stol 0.0 -grashof 4e4 -snes_fas_smoothup 6 -snes_fas_smoothdown 6 -lidvelocity 100
829       requires: !single
830 
831    test:
832       suffix: fdcoloring_ds
833       args: -da_refine 3 -snes_converged_reason -pc_type mg -mat_fd_type ds
834       output_file: output/ex19_2.out
835       requires: !single
836 
837    test:
838       suffix: fdcoloring_ds_baij
839       args: -da_refine 3 -snes_converged_reason -pc_type mg -mat_fd_type ds -dm_mat_type baij
840       output_file: output/ex19_2.out
841       requires: !single
842 
843    test:
844       suffix: fdcoloring_ds_bcols1
845       args: -da_refine 3 -snes_converged_reason -pc_type mg -mat_fd_type ds -mat_fd_coloring_bcols 1
846       output_file: output/ex19_2.out
847       requires: !single
848 
849    test:
850       suffix: fdcoloring_wp
851       args: -da_refine 3 -snes_monitor_short -pc_type mg
852       requires: !single
853 
854    test:
855       suffix: fdcoloring_wp_baij
856       args: -da_refine 3 -snes_monitor_short -pc_type mg -dm_mat_type baij
857       output_file: output/ex19_fdcoloring_wp.out
858       requires: !single
859 
860    test:
861       suffix: fdcoloring_wp_bcols1
862       args: -da_refine 3 -snes_monitor_short -pc_type mg -mat_fd_coloring_bcols 1
863       output_file: output/ex19_fdcoloring_wp.out
864       requires: !single
865 
866    test:
867       suffix: fieldsplit_2
868       args: -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_block_size 4 -pc_fieldsplit_type additive -pc_fieldsplit_0_fields 0,1,2 -pc_fieldsplit_1_fields 3 -snes_monitor_short -ksp_monitor_short
869       requires: !single
870 
871    test:
872       suffix: fieldsplit_3
873       args: -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_block_size 4 -pc_fieldsplit_type additive -pc_fieldsplit_0_fields 0,1,2 -pc_fieldsplit_1_fields 3 -fieldsplit_0_pc_type lu -fieldsplit_1_pc_type lu -snes_monitor_short -ksp_monitor_short
874       requires: !single
875 
876    test:
877       suffix: fieldsplit_4
878       args: -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_block_size 4 -pc_fieldsplit_type SCHUR -pc_fieldsplit_0_fields 0,1,2 -pc_fieldsplit_1_fields 3 -fieldsplit_0_pc_type lu -fieldsplit_1_pc_type lu -snes_monitor_short -ksp_monitor_short
879       requires: !single
880 
881    # HYPRE PtAP broken with complex numbers
882    test:
883       suffix: fieldsplit_hypre
884       nsize: 2
885       requires: hypre mumps !complex !defined(PETSC_HAVE_HYPRE_DEVICE)
886       args: -pc_type fieldsplit -pc_fieldsplit_block_size 4 -pc_fieldsplit_type SCHUR -pc_fieldsplit_0_fields 0,1,2 -pc_fieldsplit_1_fields 3 -fieldsplit_0_pc_type lu -fieldsplit_0_pc_factor_mat_solver_type mumps -fieldsplit_1_pc_type hypre -fieldsplit_1_pc_hypre_type boomeramg -snes_monitor_short -ksp_monitor_short
887 
888    test:
889       suffix: fieldsplit_mumps
890       nsize: 2
891       requires: mumps
892       args: -pc_type fieldsplit -pc_fieldsplit_block_size 4 -pc_fieldsplit_type SCHUR -pc_fieldsplit_0_fields 0,1,2 -pc_fieldsplit_1_fields 3 -fieldsplit_0_pc_type lu -fieldsplit_1_pc_type lu -snes_monitor_short -ksp_monitor_short -fieldsplit_0_pc_factor_mat_solver_type mumps -fieldsplit_1_pc_factor_mat_solver_type mumps
893       output_file: output/ex19_fieldsplit_5.out
894 
895    test:
896       suffix: greedy_coloring
897       nsize: 2
898       args: -da_refine 3 -snes_monitor_short -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -mat_coloring_weight_type lf -mat_coloring_view> ex19_greedy_coloring.tmp 2>&1
899       requires: !single
900 
901    # HYPRE PtAP broken with complex numbers
902    test:
903       suffix: hypre
904       nsize: 2
905       requires: hypre !complex !defined(PETSC_HAVE_HYPRE_DEVICE)
906       args: -da_refine 3 -snes_monitor_short -pc_type hypre -ksp_norm_type unpreconditioned
907 
908    # ibcgs is broken when using device vectors
909    test:
910       suffix: ibcgs
911       nsize: 2
912       args: -ksp_type ibcgs -ksp_monitor_short -da_refine 2 -snes_view
913       requires: !complex !single
914 
915    test:
916       suffix: kaczmarz
917       nsize: 2
918       args: -pc_type kaczmarz -ksp_monitor_short -snes_monitor_short -snes_view
919       requires: !single
920 
921    test:
922       suffix: klu
923       requires: suitesparse
924       args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type klu
925       output_file: output/ex19_superlu.out
926 
927    test:
928       suffix: klu_2
929       requires: suitesparse
930       args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type klu -pc_factor_mat_ordering_type nd
931       output_file: output/ex19_superlu.out
932 
933    test:
934       suffix: klu_3
935       requires: suitesparse
936       args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type klu -mat_klu_use_btf 0
937       output_file: output/ex19_superlu.out
938 
939    test:
940       suffix: ml
941       nsize: 2
942       requires: ml
943       args: -da_refine 3 -snes_monitor_short -pc_type ml
944 
945    test:
946       suffix: ngmres_fas
947       args: -da_refine 4 -snes_monitor_short -snes_type ngmres -npc_fas_levels_snes_type ngs -npc_fas_levels_snes_ngs_sweeps 3 -npc_fas_levels_snes_ngs_atol 0.0 -npc_fas_levels_snes_ngs_stol 0.0 -npc_snes_type fas -npc_fas_levels_snes_type ngs -npc_snes_max_it 1 -npc_snes_fas_smoothup 6 -npc_snes_fas_smoothdown 6 -lidvelocity 100 -grashof 4e4
948       requires: !single
949 
950    test:
951       suffix: ngmres_fas_gssecant
952       args: -da_refine 3 -snes_monitor_short -snes_type ngmres -npc_snes_type fas -npc_fas_levels_snes_type ngs -npc_fas_levels_snes_max_it 6 -npc_fas_levels_snes_ngs_secant -npc_fas_levels_snes_ngs_max_it 1 -npc_fas_coarse_snes_max_it 1 -lidvelocity 100 -grashof 4e4
953       requires: !single
954 
955    test:
956       suffix: ngmres_fas_ms
957       nsize: 2
958       args: -snes_grid_sequence 2 -lidvelocity 200 -grashof 1e4 -snes_monitor_short -snes_view -snes_converged_reason -snes_type ngmres -npc_snes_type fas -npc_fas_coarse_snes_type newtonls -npc_fas_coarse_ksp_type preonly -npc_snes_max_it 1
959       requires: !single
960 
961    test:
962       suffix: ngmres_nasm
963       nsize: 4
964       args: -da_refine 4 -da_overlap 2 -snes_monitor_short -snes_type ngmres -snes_max_it 10 -npc_snes_type nasm -npc_snes_nasm_type basic -grashof 4e4 -lidvelocity 100
965       requires: !single
966 
967    test:
968       suffix: ngs
969       args: -snes_type ngs -snes_view -snes_monitor -snes_rtol 1e-4
970       requires: !single
971 
972    test:
973       suffix: ngs_fd
974       args: -snes_type ngs -snes_ngs_secant -snes_view -snes_monitor -snes_rtol 1e-4
975       requires: !single
976 
977    test:
978       suffix: parms
979       nsize: 2
980       requires: parms
981       args: -pc_type parms -ksp_monitor_short -snes_view
982 
983    test:
984       suffix: superlu
985       requires: superlu
986       args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type superlu
987 
988    test:
989       suffix: superlu_sell
990       requires: superlu
991       args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type superlu -dm_mat_type sell -pc_factor_mat_ordering_type natural
992       output_file: output/ex19_superlu.out
993 
994    test:
995       suffix: superlu_dist
996       requires: superlu_dist
997       args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type superlu_dist
998       output_file: output/ex19_superlu.out
999 
1000    test:
1001       suffix: superlu_dist_2
1002       nsize: 2
1003       requires: superlu_dist
1004       args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type superlu_dist
1005       output_file: output/ex19_superlu.out
1006 
1007    test:
1008       suffix: superlu_dist_3d
1009       nsize: 4
1010       requires: superlu_dist !defined(PETSCTEST_VALGRIND)
1011       filter: grep -v iam | grep -v openMP
1012       args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type superlu_dist -mat_superlu_dist_3d -mat_superlu_dist_d 2 -snes_view -snes_monitor -ksp_monitor
1013 
1014    test:
1015       suffix: superlu_equil
1016       requires: superlu
1017       args: -da_grid_x 20 -da_grid_y 20 -{snes,ksp}_monitor_short -pc_type lu -pc_factor_mat_solver_type superlu -mat_superlu_equil
1018 
1019    test:
1020       suffix: superlu_equil_sell
1021       requires: superlu
1022       args: -da_grid_x 20 -da_grid_y 20 -{snes,ksp}_monitor_short -pc_type lu -pc_factor_mat_solver_type superlu -mat_superlu_equil -dm_mat_type sell -pc_factor_mat_ordering_type natural
1023       output_file: output/ex19_superlu_equil.out
1024 
1025    test:
1026       suffix: tcqmr
1027       args: -da_refine 1 -ksp_monitor_short -ksp_type tcqmr
1028       requires: !single
1029 
1030    test:
1031       suffix: tfqmr
1032       args: -da_refine 1 -ksp_monitor_short -ksp_type tfqmr
1033       requires: !single
1034 
1035    test:
1036       suffix: umfpack
1037       requires: suitesparse
1038       args: -da_refine 2 -pc_type lu -pc_factor_mat_solver_type umfpack -snes_view -snes_monitor_short -ksp_monitor_short -pc_factor_mat_ordering_type external
1039 
1040    test:
1041       suffix: tut_1
1042       nsize: 4
1043       requires: !single
1044       args: -da_refine 5 -snes_monitor -ksp_monitor -snes_view
1045 
1046    test:
1047       suffix: tut_2
1048       nsize: 4
1049       requires: !single
1050       args: -da_refine 5 -snes_monitor -ksp_monitor -snes_view -pc_type mg
1051 
1052    # HYPRE PtAP broken with complex numbers
1053    test:
1054       suffix: tut_3
1055       nsize: 4
1056       requires: hypre !single !complex !defined(PETSC_HAVE_HYPRE_DEVICE)
1057       args: -da_refine 5 -snes_monitor -ksp_monitor -snes_view -pc_type hypre
1058 
1059    test:
1060       suffix: tut_8
1061       nsize: 4
1062       requires: ml !single
1063       args: -da_refine 5 -snes_monitor -ksp_monitor -snes_view -pc_type ml
1064 
1065    test:
1066       suffix: tut_4
1067       nsize: 1
1068       requires: !single
1069       args: -da_refine 5 -log_view
1070       filter: head -n 2
1071       filter_output: head -n 2
1072 
1073    test:
1074       suffix: tut_5
1075       nsize: 1
1076       requires: !single
1077       args: -da_refine 5 -log_view -pc_type mg
1078       filter: head -n 2
1079       filter_output: head -n 2
1080 
1081    test:
1082       suffix: tut_6
1083       nsize: 4
1084       requires: !single
1085       args: -da_refine 5 -log_view
1086       filter: head -n 2
1087       filter_output: head -n 2
1088 
1089    test:
1090       suffix: tut_7
1091       nsize: 4
1092       requires: !single
1093       args: -da_refine 5 -log_view -pc_type mg
1094       filter: head -n 2
1095       filter_output: head -n 2
1096 
1097    test:
1098       suffix: cuda_1
1099       nsize: 1
1100       requires: cuda
1101       args: -snes_monitor -dm_mat_type seqaijcusparse -dm_vec_type seqcuda -pc_type gamg -pc_gamg_esteig_ksp_max_it 10 -ksp_monitor -mg_levels_ksp_max_it 3
1102 
1103    test:
1104       suffix: cuda_2
1105       nsize: 3
1106       requires: cuda !single
1107       args: -snes_monitor -dm_mat_type mpiaijcusparse -dm_vec_type mpicuda -pc_type gamg -pc_gamg_esteig_ksp_max_it 10 -ksp_monitor  -mg_levels_ksp_max_it 3
1108 
1109    test:
1110       suffix: cuda_dm_bind_below
1111       nsize: 2
1112       requires: cuda
1113       args: -dm_mat_type aijcusparse -dm_vec_type cuda -da_refine 3 -pc_type mg -mg_levels_ksp_type chebyshev -mg_levels_pc_type jacobi -log_view -pc_mg_log -dm_bind_below 10000
1114       filter: awk "/Level/ {print \$NF}"
1115 
1116    test:
1117       suffix: viennacl_dm_bind_below
1118       nsize: 2
1119       requires: viennacl
1120       args: -dm_mat_type aijviennacl -dm_vec_type viennacl -da_refine 3 -pc_type mg -mg_levels_ksp_type chebyshev -mg_levels_pc_type jacobi -log_view -pc_mg_log -dm_bind_below 10000
1121       filter: awk "/Level/ {print \$NF}"
1122 
1123    test:
1124       suffix: seqbaijmkl
1125       nsize: 1
1126       requires: defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE)
1127       args: -dm_mat_type baij -snes_monitor -ksp_monitor -snes_view
1128 
1129    test:
1130       suffix: mpibaijmkl
1131       nsize: 2
1132       requires:  defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE)
1133       args: -dm_mat_type baij -snes_monitor -ksp_monitor -snes_view
1134 
1135    test:
1136      suffix: cpardiso
1137      nsize: 4
1138      requires: mkl_cpardiso
1139      args: -pc_type lu -pc_factor_mat_solver_type mkl_cpardiso -ksp_monitor
1140 
1141    test:
1142      suffix: logviewmemory
1143      requires: defined(PETSC_USE_LOG) !defined(PETSCTEST_VALGRIND)
1144      args: -log_view -log_view_memory -da_refine 4
1145      filter: grep MatFDColorSetUp | wc -w | xargs  -I % sh -c "expr % \> 21"
1146 
1147    test:
1148      suffix: fs
1149      args: -pc_type fieldsplit -da_refine 3  -all_ksp_monitor -fieldsplit_y_velocity_pc_type lu  -fieldsplit_temperature_pc_type lu -fieldsplit_x_velocity_pc_type lu  -snes_view
1150 
1151    test:
1152      suffix: asm_matconvert
1153      args: -mat_type aij -pc_type asm -pc_asm_sub_mat_type dense -snes_view
1154 
1155    test:
1156       suffix: euclid
1157       nsize: 2
1158       requires: hypre !single !complex !defined(PETSC_HAVE_HYPRE_MIXEDINT) !defined(PETSC_HAVE_HYPRE_DEVICE)
1159       args: -da_refine 2 -ksp_monitor -snes_monitor -snes_view -pc_type hypre -pc_hypre_type euclid
1160 
1161    test:
1162       suffix: euclid_bj
1163       nsize: 2
1164       requires: hypre !single !complex !defined(PETSC_HAVE_HYPRE_MIXEDINT) !defined(PETSC_HAVE_HYPRE_DEVICE)
1165       args: -da_refine 2 -ksp_monitor -snes_monitor -snes_view -pc_type hypre -pc_hypre_type euclid -pc_hypre_euclid_bj
1166 
1167    test:
1168       suffix: euclid_droptolerance
1169       nsize: 1
1170       requires: hypre !single !complex !defined(PETSC_HAVE_HYPRE_MIXEDINT) !defined(PETSC_HAVE_HYPRE_DEVICE)
1171       args: -da_refine 2 -ksp_monitor -snes_monitor -snes_view -pc_type hypre -pc_hypre_type euclid -pc_hypre_euclid_droptolerance .1
1172 
1173    test:
1174       suffix: failure_size
1175       nsize: 1
1176       requires: !defined(PETSC_USE_64BIT_INDICES) !defined(PETSCTEST_VALGRIND)
1177       args: -da_refine 100 -petsc_ci_portable_error_output -error_output_stdout
1178       filter: egrep -v "(options_left|memory block|leaked context|is not freed before MPI_Finalize)"
1179 
1180 TEST*/
1181