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