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