Actual source code: ex19.c


  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' = '>' */

 12: /*
 13:       See src/ksp/ksp/tutorials/ex45.c
 14: */

 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*/


We thank David E. Keyes for contributing the driven cavity discretization within this example code.

This problem is modeled by the partial differential equation system

\begin{eqnarray}
- \triangle U - \nabla_y \Omega & = & 0 \\
- \triangle V + \nabla_x\Omega & = & 0 \\
- \triangle \Omega + \nabla \cdot ([U*\Omega,V*\Omega]) - GR* \nabla_x T & = & 0 \\
- \triangle T + PR* \nabla \cdot ([U*T,V*T]) & = & 0
\end{eqnarray}

in the unit square, which is uniformly discretized in each of x and y in this simple encoding.

No-slip, rigid-wall Dirichlet conditions are used for $ [U,V]$.
Dirichlet conditions are used for Omega, based on the definition of
vorticity: $ \Omega = - \nabla_y U + \nabla_x V$, where along each
constant coordinate boundary, the tangential derivative is zero.
Dirichlet conditions are used for T on the left and right walls,
and insulation homogeneous Neumann conditions are used for T on
the top and bottom walls.

A finite difference approximation with the usual 5-point stencil
is used to discretize the boundary value problem to obtain a
nonlinear system of equations. Upwinding is used for the divergence
(convective) terms and central for the gradient (source) terms.

The Jacobian can be either
* formed via finite differencing using coloring (the default), or
* applied matrix-free via the option -snes_mf
(for larger grid problems this variant may not converge
without a preconditioner due to ill-conditioning).

 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

 78: /*
 79:    User-defined routines and data structures
 80: */
 81: typedef struct {
 82:   PetscScalar u,v,omega,temp;
 83: } Field;

 85: PetscErrorCode FormFunctionLocal(DMDALocalInfo*,Field**,Field**,void*);

 87: typedef struct {
 88:   PetscReal   lidvelocity,prandtl,grashof;  /* physical parameters */
 89:   PetscBool   draw_contours;                /* flag - 1 indicates drawing contours */
 90: } AppCtx;

 92: extern PetscErrorCode FormInitialGuess(AppCtx*,DM,Vec);
 93: extern PetscErrorCode NonlinearGS(SNES,Vec,Vec,void*);

 95: int main(int argc,char **argv)
 96: {
 97:   AppCtx         user;                /* user-defined work context */
 98:   PetscInt       mx,my,its;
100:   MPI_Comm       comm;
101:   SNES           snes;
102:   DM             da;
103:   Vec            x;

105:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;

108:   comm = PETSC_COMM_WORLD;
109:   SNESCreate(comm,&snes);

111:   /*
112:       Create distributed array object to manage parallel grid and vectors
113:       for principal unknowns (x) and governing residuals (f)
114:   */
115:   DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,4,4,PETSC_DECIDE,PETSC_DECIDE,4,1,0,0,&da);
116:   DMSetFromOptions(da);
117:   DMSetUp(da);
118:   SNESSetDM(snes,(DM)da);
119:   SNESSetNGS(snes, NonlinearGS, (void*)&user);

121:   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);
122:   /*
123:      Problem parameters (velocity of lid, prandtl, and grashof numbers)
124:   */
125:   user.lidvelocity = 1.0/(mx*my);
126:   user.prandtl     = 1.0;
127:   user.grashof     = 1.0;

129:   PetscOptionsGetReal(NULL,NULL,"-lidvelocity",&user.lidvelocity,NULL);
130:   PetscOptionsGetReal(NULL,NULL,"-prandtl",&user.prandtl,NULL);
131:   PetscOptionsGetReal(NULL,NULL,"-grashof",&user.grashof,NULL);
132:   PetscOptionsHasName(NULL,NULL,"-contours",&user.draw_contours);

134:   DMDASetFieldName(da,0,"x_velocity");
135:   DMDASetFieldName(da,1,"y_velocity");
136:   DMDASetFieldName(da,2,"Omega");
137:   DMDASetFieldName(da,3,"temperature");

139:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140:      Create user context, set problem data, create vector data structures.
141:      Also, compute the initial guess.
142:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

144:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145:      Create nonlinear solver context

147:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
148:   DMSetApplicationContext(da,&user);
149:   DMDASNESSetFunctionLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))FormFunctionLocal,&user);
150:   SNESSetFromOptions(snes);
151:   PetscPrintf(comm,"lid velocity = %g, prandtl # = %g, grashof # = %g\n",(double)user.lidvelocity,(double)user.prandtl,(double)user.grashof);

153:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
154:      Solve the nonlinear system
155:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
156:   DMCreateGlobalVector(da,&x);
157:   FormInitialGuess(&user,da,x);

159:   SNESSolve(snes,NULL,x);

161:   SNESGetIterationNumber(snes,&its);
162:   PetscPrintf(comm,"Number of SNES iterations = %D\n", its);

164:   /*
165:      Visualize solution
166:   */
167:   if (user.draw_contours) {
168:     VecView(x,PETSC_VIEWER_DRAW_WORLD);
169:   }

171:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
172:      Free work space.  All PETSc objects should be destroyed when they
173:      are no longer needed.
174:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
175:   VecDestroy(&x);
176:   DMDestroy(&da);
177:   SNESDestroy(&snes);
178:   PetscFinalize();
179:   return ierr;
180: }

182: /* ------------------------------------------------------------------- */

184: /*
185:    FormInitialGuess - Forms initial approximation.

187:    Input Parameters:
188:    user - user-defined application context
189:    X - vector

191:    Output Parameter:
192:    X - vector
193: */
194: PetscErrorCode FormInitialGuess(AppCtx *user,DM da,Vec X)
195: {
196:   PetscInt       i,j,mx,xs,ys,xm,ym;
198:   PetscReal      grashof,dx;
199:   Field          **x;

202:   grashof = user->grashof;

204:   DMDAGetInfo(da,0,&mx,0,0,0,0,0,0,0,0,0,0,0);
205:   dx   = 1.0/(mx-1);

207:   /*
208:      Get local grid boundaries (for 2-dimensional DMDA):
209:        xs, ys   - starting grid indices (no ghost points)
210:        xm, ym   - widths of local grid (no ghost points)
211:   */
212:   DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);

214:   /*
215:      Get a pointer to vector data.
216:        - For default PETSc vectors, VecGetArray() returns a pointer to
217:          the data array.  Otherwise, the routine is implementation dependent.
218:        - You MUST call VecRestoreArray() when you no longer need access to
219:          the array.
220:   */
221:   DMDAVecGetArrayWrite(da,X,&x);

223:   /*
224:      Compute initial guess over the locally owned part of the grid
225:      Initial condition is motionless fluid and equilibrium temperature
226:   */
227:   for (j=ys; j<ys+ym; j++) {
228:     for (i=xs; i<xs+xm; i++) {
229:       x[j][i].u     = 0.0;
230:       x[j][i].v     = 0.0;
231:       x[j][i].omega = 0.0;
232:       x[j][i].temp  = (grashof>0)*i*dx;
233:     }
234:   }

236:   /*
237:      Restore vector
238:   */
239:   DMDAVecRestoreArrayWrite(da,X,&x);
240:   return(0);
241: }

243: PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,Field **x,Field **f,void *ptr)
244: {
245:   AppCtx         *user = (AppCtx*)ptr;
247:   PetscInt       xints,xinte,yints,yinte,i,j;
248:   PetscReal      hx,hy,dhx,dhy,hxdhy,hydhx;
249:   PetscReal      grashof,prandtl,lid;
250:   PetscScalar    u,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;

253:   grashof = user->grashof;
254:   prandtl = user->prandtl;
255:   lid     = user->lidvelocity;

257:   /*
258:      Define mesh intervals ratios for uniform grid.

260:      Note: FD formulae below are normalized by multiplying through by
261:      local volume element (i.e. hx*hy) to obtain coefficients O(1) in two dimensions.

263:   */
264:   dhx   = (PetscReal)(info->mx-1);  dhy = (PetscReal)(info->my-1);
265:   hx    = 1.0/dhx;                   hy = 1.0/dhy;
266:   hxdhy = hx*dhy;                 hydhx = hy*dhx;

268:   xints = info->xs; xinte = info->xs+info->xm; yints = info->ys; yinte = info->ys+info->ym;

270:   /* Test whether we are on the bottom edge of the global array */
271:   if (yints == 0) {
272:     j     = 0;
273:     yints = yints + 1;
274:     /* bottom edge */
275:     for (i=info->xs; i<info->xs+info->xm; i++) {
276:       f[j][i].u     = x[j][i].u;
277:       f[j][i].v     = x[j][i].v;
278:       f[j][i].omega = x[j][i].omega + (x[j+1][i].u - x[j][i].u)*dhy;
279:       f[j][i].temp  = x[j][i].temp-x[j+1][i].temp;
280:     }
281:   }

283:   /* Test whether we are on the top edge of the global array */
284:   if (yinte == info->my) {
285:     j     = info->my - 1;
286:     yinte = yinte - 1;
287:     /* top edge */
288:     for (i=info->xs; i<info->xs+info->xm; i++) {
289:       f[j][i].u     = x[j][i].u - lid;
290:       f[j][i].v     = x[j][i].v;
291:       f[j][i].omega = x[j][i].omega + (x[j][i].u - x[j-1][i].u)*dhy;
292:       f[j][i].temp  = x[j][i].temp-x[j-1][i].temp;
293:     }
294:   }

296:   /* Test whether we are on the left edge of the global array */
297:   if (xints == 0) {
298:     i     = 0;
299:     xints = xints + 1;
300:     /* left edge */
301:     for (j=info->ys; j<info->ys+info->ym; j++) {
302:       f[j][i].u     = x[j][i].u;
303:       f[j][i].v     = x[j][i].v;
304:       f[j][i].omega = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx;
305:       f[j][i].temp  = x[j][i].temp;
306:     }
307:   }

309:   /* Test whether we are on the right edge of the global array */
310:   if (xinte == info->mx) {
311:     i     = info->mx - 1;
312:     xinte = xinte - 1;
313:     /* right edge */
314:     for (j=info->ys; j<info->ys+info->ym; j++) {
315:       f[j][i].u     = x[j][i].u;
316:       f[j][i].v     = x[j][i].v;
317:       f[j][i].omega = x[j][i].omega - (x[j][i].v - x[j][i-1].v)*dhx;
318:       f[j][i].temp  = x[j][i].temp - (PetscReal)(grashof>0);
319:     }
320:   }

322:   /* Compute over the interior points */
323:   for (j=yints; j<yinte; j++) {
324:     for (i=xints; i<xinte; i++) {

326:       /*
327:        convective coefficients for upwinding
328:       */
329:       vx  = x[j][i].u; avx = PetscAbsScalar(vx);
330:       vxp = .5*(vx+avx); vxm = .5*(vx-avx);
331:       vy  = x[j][i].v; avy = PetscAbsScalar(vy);
332:       vyp = .5*(vy+avy); vym = .5*(vy-avy);

334:       /* U velocity */
335:       u         = x[j][i].u;
336:       uxx       = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
337:       uyy       = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
338:       f[j][i].u = uxx + uyy - .5*(x[j+1][i].omega-x[j-1][i].omega)*hx;

340:       /* V velocity */
341:       u         = x[j][i].v;
342:       uxx       = (2.0*u - x[j][i-1].v - x[j][i+1].v)*hydhx;
343:       uyy       = (2.0*u - x[j-1][i].v - x[j+1][i].v)*hxdhy;
344:       f[j][i].v = uxx + uyy + .5*(x[j][i+1].omega-x[j][i-1].omega)*hy;

346:       /* Omega */
347:       u             = x[j][i].omega;
348:       uxx           = (2.0*u - x[j][i-1].omega - x[j][i+1].omega)*hydhx;
349:       uyy           = (2.0*u - x[j-1][i].omega - x[j+1][i].omega)*hxdhy;
350:       f[j][i].omega = uxx + uyy + (vxp*(u - x[j][i-1].omega) + vxm*(x[j][i+1].omega - u))*hy +
351:                       (vyp*(u - x[j-1][i].omega) + vym*(x[j+1][i].omega - u))*hx -
352:                       .5*grashof*(x[j][i+1].temp - x[j][i-1].temp)*hy;

354:       /* Temperature */
355:       u            = x[j][i].temp;
356:       uxx          = (2.0*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx;
357:       uyy          = (2.0*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy;
358:       f[j][i].temp =  uxx + uyy  + prandtl*((vxp*(u - x[j][i-1].temp) + vxm*(x[j][i+1].temp - u))*hy +
359:                                             (vyp*(u - x[j-1][i].temp) + vym*(x[j+1][i].temp - u))*hx);
360:     }
361:   }

363:   /*
364:      Flop count (multiply-adds are counted as 2 operations)
365:   */
366:   PetscLogFlops(84.0*info->ym*info->xm);
367:   return(0);
368: }

370: /*
371:     Performs sweeps of point block nonlinear Gauss-Seidel on all the local grid points
372: */
373: PetscErrorCode NonlinearGS(SNES snes, Vec X, Vec B, void *ctx)
374: {
375:   DMDALocalInfo  info;
376:   Field          **x,**b;
378:   Vec            localX, localB;
379:   DM             da;
380:   PetscInt       xints,xinte,yints,yinte,i,j,k,l;
381:   PetscInt       max_its,tot_its;
382:   PetscInt       sweeps;
383:   PetscReal      rtol,atol,stol;
384:   PetscReal      hx,hy,dhx,dhy,hxdhy,hydhx;
385:   PetscReal      grashof,prandtl,lid;
386:   PetscScalar    u,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;
387:   PetscScalar    fu, fv, fomega, ftemp;
388:   PetscScalar    dfudu;
389:   PetscScalar    dfvdv;
390:   PetscScalar    dfodu, dfodv, dfodo;
391:   PetscScalar    dftdu, dftdv, dftdt;
392:   PetscScalar    yu=0, yv=0, yo=0, yt=0;
393:   PetscScalar    bjiu, bjiv, bjiomega, bjitemp;
394:   PetscBool      ptconverged;
395:   PetscReal      pfnorm,pfnorm0,pynorm,pxnorm;
396:   AppCtx         *user = (AppCtx*)ctx;

399:   grashof = user->grashof;
400:   prandtl = user->prandtl;
401:   lid     = user->lidvelocity;
402:   tot_its = 0;
403:   SNESNGSGetTolerances(snes,&rtol,&atol,&stol,&max_its);
404:   SNESNGSGetSweeps(snes,&sweeps);
405:   SNESGetDM(snes,(DM*)&da);
406:   DMGetLocalVector(da,&localX);
407:   if (B) {
408:     DMGetLocalVector(da,&localB);
409:   }
410:   /*
411:      Scatter ghost points to local vector, using the 2-step process
412:         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
413:   */
414:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
415:   DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
416:   if (B) {
417:     DMGlobalToLocalBegin(da,B,INSERT_VALUES,localB);
418:     DMGlobalToLocalEnd(da,B,INSERT_VALUES,localB);
419:   }
420:   DMDAGetLocalInfo(da,&info);
421:   DMDAVecGetArrayWrite(da,localX,&x);
422:   if (B) {
423:     DMDAVecGetArrayRead(da,localB,&b);
424:   }
425:   /* looks like a combination of the formfunction / formjacobian routines */
426:   dhx   = (PetscReal)(info.mx-1);dhy   = (PetscReal)(info.my-1);
427:   hx    = 1.0/dhx;               hy    = 1.0/dhy;
428:   hxdhy = hx*dhy;                hydhx = hy*dhx;

430:   xints = info.xs; xinte = info.xs+info.xm; yints = info.ys; yinte = info.ys+info.ym;

432:   /* Set the boundary conditions on the momentum equations */
433:   /* Test whether we are on the bottom edge of the global array */
434:   if (yints == 0) {
435:     j     = 0;
436:     /* bottom edge */
437:     for (i=info.xs; i<info.xs+info.xm; i++) {

439:       if (B) {
440:         bjiu = b[j][i].u;
441:         bjiv = b[j][i].v;
442:       } else {
443:         bjiu = 0.0;
444:         bjiv = 0.0;
445:       }
446:       x[j][i].u = 0.0 + bjiu;
447:       x[j][i].v = 0.0 + bjiv;
448:     }
449:   }

451:   /* Test whether we are on the top edge of the global array */
452:   if (yinte == info.my) {
453:     j     = info.my - 1;
454:     /* top edge */
455:     for (i=info.xs; i<info.xs+info.xm; i++) {
456:       if (B) {
457:         bjiu = b[j][i].u;
458:         bjiv = b[j][i].v;
459:       } else {
460:         bjiu = 0.0;
461:         bjiv = 0.0;
462:       }
463:       x[j][i].u = lid + bjiu;
464:       x[j][i].v = bjiv;
465:     }
466:   }

468:   /* Test whether we are on the left edge of the global array */
469:   if (xints == 0) {
470:     i     = 0;
471:     /* left edge */
472:     for (j=info.ys; j<info.ys+info.ym; j++) {
473:       if (B) {
474:         bjiu = b[j][i].u;
475:         bjiv = b[j][i].v;
476:       } else {
477:         bjiu = 0.0;
478:         bjiv = 0.0;
479:       }
480:       x[j][i].u = 0.0 + bjiu;
481:       x[j][i].v = 0.0 + bjiv;
482:     }
483:   }

485:   /* Test whether we are on the right edge of the global array */
486:   if (xinte == info.mx) {
487:     i     = info.mx - 1;
488:     /* right edge */
489:     for (j=info.ys; j<info.ys+info.ym; j++) {
490:       if (B) {
491:         bjiu = b[j][i].u;
492:         bjiv = b[j][i].v;
493:       } else {
494:         bjiu = 0.0;
495:         bjiv = 0.0;
496:       }
497:       x[j][i].u = 0.0 + bjiu;
498:       x[j][i].v = 0.0 + bjiv;
499:     }
500:   }

502:   for (k=0; k < sweeps; k++) {
503:     for (j=info.ys; j<info.ys + info.ym; j++) {
504:       for (i=info.xs; i<info.xs + info.xm; i++) {
505:         ptconverged = PETSC_FALSE;
506:         pfnorm0     = 0.0;
507:         fu          = 0.0;
508:         fv          = 0.0;
509:         fomega      = 0.0;
510:         ftemp       = 0.0;
511:         /*  Run Newton's method on a single grid point */
512:         for (l = 0; l < max_its && !ptconverged; l++) {
513:           if (B) {
514:             bjiu     = b[j][i].u;
515:             bjiv     = b[j][i].v;
516:             bjiomega = b[j][i].omega;
517:             bjitemp  = b[j][i].temp;
518:           } else {
519:             bjiu     = 0.0;
520:             bjiv     = 0.0;
521:             bjiomega = 0.0;
522:             bjitemp  = 0.0;
523:           }

525:           if (i != 0 && i != info.mx - 1 && j != 0 && j != info.my-1) {
526:             /* U velocity */
527:             u     = x[j][i].u;
528:             uxx   = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
529:             uyy   = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
530:             fu    = uxx + uyy - .5*(x[j+1][i].omega-x[j-1][i].omega)*hx - bjiu;
531:             dfudu = 2.0*(hydhx + hxdhy);
532:             /* V velocity */
533:             u     = x[j][i].v;
534:             uxx   = (2.0*u - x[j][i-1].v - x[j][i+1].v)*hydhx;
535:             uyy   = (2.0*u - x[j-1][i].v - x[j+1][i].v)*hxdhy;
536:             fv    = uxx + uyy + .5*(x[j][i+1].omega-x[j][i-1].omega)*hy - bjiv;
537:             dfvdv = 2.0*(hydhx + hxdhy);
538:             /*
539:              convective coefficients for upwinding
540:              */
541:             vx  = x[j][i].u; avx = PetscAbsScalar(vx);
542:             vxp = .5*(vx+avx); vxm = .5*(vx-avx);
543:             vy  = x[j][i].v; avy = PetscAbsScalar(vy);
544:             vyp = .5*(vy+avy); vym = .5*(vy-avy);
545:             /* Omega */
546:             u      = x[j][i].omega;
547:             uxx    = (2.0*u - x[j][i-1].omega - x[j][i+1].omega)*hydhx;
548:             uyy    = (2.0*u - x[j-1][i].omega - x[j+1][i].omega)*hxdhy;
549:             fomega = uxx + uyy +  (vxp*(u - x[j][i-1].omega) + vxm*(x[j][i+1].omega - u))*hy +
550:                      (vyp*(u - x[j-1][i].omega) + vym*(x[j+1][i].omega - u))*hx -
551:                      .5*grashof*(x[j][i+1].temp - x[j][i-1].temp)*hy - bjiomega;
552:             /* convective coefficient derivatives */
553:             dfodo = 2.0*(hydhx + hxdhy) + ((vxp - vxm)*hy + (vyp - vym)*hx);
554:             if (PetscRealPart(vx) > 0.0) dfodu = (u - x[j][i-1].omega)*hy;
555:             else dfodu = (x[j][i+1].omega - u)*hy;

557:             if (PetscRealPart(vy) > 0.0) dfodv = (u - x[j-1][i].omega)*hx;
558:             else dfodv = (x[j+1][i].omega - u)*hx;

560:             /* Temperature */
561:             u     = x[j][i].temp;
562:             uxx   = (2.0*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx;
563:             uyy   = (2.0*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy;
564:             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;
565:             dftdt = 2.0*(hydhx + hxdhy) + prandtl*((vxp - vxm)*hy + (vyp - vym)*hx);
566:             if (PetscRealPart(vx) > 0.0) dftdu = prandtl*(u - x[j][i-1].temp)*hy;
567:             else dftdu = prandtl*(x[j][i+1].temp - u)*hy;

569:             if (PetscRealPart(vy) > 0.0) dftdv = prandtl*(u - x[j-1][i].temp)*hx;
570:             else dftdv = prandtl*(x[j+1][i].temp - u)*hx;

572:             /* invert the system:
573:              [ dfu / du     0        0        0    ][yu] = [fu]
574:              [     0    dfv / dv     0        0    ][yv]   [fv]
575:              [ dfo / du dfo / dv dfo / do     0    ][yo]   [fo]
576:              [ dft / du dft / dv     0    dft / dt ][yt]   [ft]
577:              by simple back-substitution
578:            */
579:             yu = fu / dfudu;
580:             yv = fv / dfvdv;
581:             yo = (fomega - (dfodu*yu + dfodv*yv)) / dfodo;
582:             yt = (ftemp - (dftdu*yu + dftdv*yv)) / dftdt;

584:             x[j][i].u     = x[j][i].u - yu;
585:             x[j][i].v     = x[j][i].v - yv;
586:             x[j][i].temp  = x[j][i].temp - yt;
587:             x[j][i].omega = x[j][i].omega - yo;
588:           }
589:           if (i == 0) {
590:             fomega        = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx - bjiomega;
591:             ftemp         = x[j][i].temp - bjitemp;
592:             yo            = fomega;
593:             yt            = ftemp;
594:             x[j][i].omega = x[j][i].omega - fomega;
595:             x[j][i].temp  = x[j][i].temp - ftemp;
596:           }
597:           if (i == info.mx - 1) {
598:             fomega        = x[j][i].omega - (x[j][i].v - x[j][i-1].v)*dhx - bjiomega;
599:             ftemp         = x[j][i].temp - (PetscReal)(grashof>0) - bjitemp;
600:             yo            = fomega;
601:             yt            = ftemp;
602:             x[j][i].omega = x[j][i].omega - fomega;
603:             x[j][i].temp  = x[j][i].temp - ftemp;
604:           }
605:           if (j == 0) {
606:             fomega        = x[j][i].omega + (x[j+1][i].u - x[j][i].u)*dhy - bjiomega;
607:             ftemp         = x[j][i].temp-x[j+1][i].temp - bjitemp;
608:             yo            = fomega;
609:             yt            = ftemp;
610:             x[j][i].omega = x[j][i].omega - fomega;
611:             x[j][i].temp  = x[j][i].temp - ftemp;
612:           }
613:           if (j == info.my - 1) {
614:             fomega        = x[j][i].omega + (x[j][i].u - x[j-1][i].u)*dhy - bjiomega;
615:             ftemp         = x[j][i].temp-x[j-1][i].temp - bjitemp;
616:             yo            = fomega;
617:             yt            = ftemp;
618:             x[j][i].omega = x[j][i].omega - fomega;
619:             x[j][i].temp  = x[j][i].temp - ftemp;
620:           }
621:           tot_its++;
622:           pfnorm = PetscRealPart(fu*fu + fv*fv + fomega*fomega + ftemp*ftemp);
623:           pfnorm = PetscSqrtReal(pfnorm);
624:           pynorm = PetscRealPart(yu*yu + yv*yv + yo*yo + yt*yt);
625:           pynorm = PetscSqrtReal(pynorm);
626:           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);
627:           pxnorm = PetscSqrtReal(pxnorm);
628:           if (l == 0) pfnorm0 = pfnorm;
629:           if (rtol*pfnorm0 >pfnorm || atol > pfnorm || pxnorm*stol > pynorm) ptconverged = PETSC_TRUE;
630:         }
631:       }
632:     }
633:   }
634:   DMDAVecRestoreArrayWrite(da,localX,&x);
635:   if (B) {
636:     DMDAVecRestoreArrayRead(da,localB,&b);
637:   }
638:   DMLocalToGlobalBegin(da,localX,INSERT_VALUES,X);
639:   DMLocalToGlobalEnd(da,localX,INSERT_VALUES,X);
640:   PetscLogFlops(tot_its*(84.0 + 41.0 + 26.0));
641:   DMRestoreLocalVector(da,&localX);
642:   if (B) {
643:     DMRestoreLocalVector(da,&localB);
644:   }
645:   return(0);
646: }

648: /*TEST

650:    test:
651:       nsize: 2
652:       args: -da_refine 3 -snes_monitor_short -pc_type mg -ksp_type fgmres -pc_mg_type full
653:       requires: !single

655:    test:
656:       suffix: 10
657:       nsize: 3
658:       args: -snes_monitor_short -ksp_monitor_short -pc_type fieldsplit -pc_fieldsplit_type symmetric_multiplicative -snes_view -da_refine 1 -ksp_type fgmres
659:       requires: !single

661:    test:
662:       suffix: 11
663:       nsize: 4
664:       requires: pastix
665:       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

667:    test:
668:       suffix: 12
669:       nsize: 12
670:       requires: pastix
671:       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

673:    test:
674:       suffix: 13
675:       nsize: 3
676:       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
677:       requires: !single

679:    test:
680:       suffix: 14
681:       nsize: 4
682:       args: -snes_monitor_short -pc_type mg -dm_mat_type baij -mg_coarse_pc_type bjacobi -da_refine 3 -ksp_type fgmres
683:       requires: !single

685:    test:
686:       suffix: 14_ds
687:       nsize: 4
688:       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
689:       output_file: output/ex19_2.out
690:       requires: !single

692:    test:
693:       suffix: 17
694:       args: -snes_monitor_short -ksp_pc_side right
695:       requires: !single

697:    test:
698:       suffix: 18
699:       args: -snes_monitor_ksp draw::draw_lg -ksp_pc_side right
700:       requires: x !single

702:    test:
703:       suffix: 2
704:       nsize: 4
705:       args: -da_refine 3 -snes_converged_reason -pc_type mg -mat_fd_type ds
706:       requires: !single

708:    test:
709:       suffix: 2_bcols1
710:       nsize: 4
711:       args: -da_refine 3 -snes_converged_reason -pc_type mg -mat_fd_type ds -mat_fd_coloring_bcols
712:       output_file: output/ex19_2.out
713:       requires: !single

715:    test:
716:       suffix: 3
717:       nsize: 4
718:       requires: mumps
719:       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

721:    test:
722:       suffix: 4
723:       nsize: 12
724:       requires: mumps
725:       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
726:       output_file: output/ex19_3.out

728:    test:
729:       suffix: 6
730:       args: -snes_monitor_short -ksp_monitor_short -pc_type fieldsplit -snes_view -ksp_type fgmres -da_refine 1
731:       requires: !single

733:    test:
734:       suffix: 7
735:       nsize: 3
736:       args: -snes_monitor_short -ksp_monitor_short -pc_type fieldsplit -snes_view -da_refine 1 -ksp_type fgmres

738:       requires: !single
739:    test:
740:       suffix: 8
741:       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
742:       requires: !single

744:    test:
745:       suffix: 9
746:       nsize: 3
747:       args: -snes_monitor_short -ksp_monitor_short -pc_type fieldsplit -pc_fieldsplit_type multiplicative -snes_view -da_refine 1 -ksp_type fgmres
748:       requires: !single

750:    test:
751:       suffix: aspin
752:       nsize: 4
753:       args: -da_refine 3 -da_overlap 2 -snes_monitor_short -snes_type aspin -grashof 4e4 -lidvelocity 100 -ksp_monitor_short
754:       requires: !single

756:    test:
757:       suffix: bcgsl
758:       nsize: 2
759:       args: -ksp_type bcgsl -ksp_monitor_short -da_refine 2 -ksp_bcgsl_ell 3 -snes_view
760:       requires: !single

762:    test:
763:       suffix: bcols1
764:       nsize: 2
765:       args: -da_refine 3 -snes_monitor_short -pc_type mg -ksp_type fgmres -pc_mg_type full -mat_fd_coloring_bcols 1
766:       output_file: output/ex19_1.out
767:       requires: !single

769:    test:
770:       suffix: bjacobi
771:       nsize: 4
772:       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
773:       requires: !single

775:    test:
776:       suffix: cgne
777:       args: -da_refine 2 -pc_type lu -ksp_type cgne -ksp_monitor_short -ksp_converged_reason -ksp_view -ksp_norm_type unpreconditioned
778:       filter: grep -v HERMITIAN
779:       requires: !single

781:    test:
782:       suffix: cgs
783:       args: -da_refine 1 -ksp_monitor_short -ksp_type cgs
784:       requires: !single

786:    test:
787:       suffix: composite_fieldsplit
788:       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
789:       requires: !single

791:    test:
792:       suffix: composite_fieldsplit_bjacobi
793:       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
794:       requires: !single

796:    test:
797:       suffix: composite_fieldsplit_bjacobi_2
798:       nsize: 4
799:       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
800:       requires: !single

802:    test:
803:       suffix: composite_gs_newton
804:       nsize: 2
805:       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
806:       requires: !single

808:    test:
809:       suffix: cuda
810:       requires: cuda !single
811:       args: -dm_vec_type cuda -dm_mat_type aijcusparse -pc_type none -ksp_type fgmres -snes_monitor_short -snes_rtol 1.e-5

813:    test:
814:       suffix: draw
815:       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
816:       requires: x !single

818:    test:
819:       suffix: drawports
820:       args: -snes_monitor_solution draw::draw_ports -da_refine 1
821:       output_file: output/ex19_draw.out
822:       requires: x !single

824:    test:
825:       suffix: fas
826:       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
827:       requires: !single

829:    test:
830:       suffix: fas_full
831:       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
832:       requires: !single

834:    test:
835:       suffix: fdcoloring_ds
836:       args: -da_refine 3 -snes_converged_reason -pc_type mg -mat_fd_type ds
837:       output_file: output/ex19_2.out
838:       requires: !single

840:    test:
841:       suffix: fdcoloring_ds_baij
842:       args: -da_refine 3 -snes_converged_reason -pc_type mg -mat_fd_type ds -dm_mat_type baij
843:       output_file: output/ex19_2.out
844:       requires: !single

846:    test:
847:       suffix: fdcoloring_ds_bcols1
848:       args: -da_refine 3 -snes_converged_reason -pc_type mg -mat_fd_type ds -mat_fd_coloring_bcols 1
849:       output_file: output/ex19_2.out
850:       requires: !single

852:    test:
853:       suffix: fdcoloring_wp
854:       args: -da_refine 3 -snes_monitor_short -pc_type mg
855:       requires: !single

857:    test:
858:       suffix: fdcoloring_wp_baij
859:       args: -da_refine 3 -snes_monitor_short -pc_type mg -dm_mat_type baij
860:       output_file: output/ex19_fdcoloring_wp.out
861:       requires: !single

863:    test:
864:       suffix: fdcoloring_wp_bcols1
865:       args: -da_refine 3 -snes_monitor_short -pc_type mg -mat_fd_coloring_bcols 1
866:       output_file: output/ex19_fdcoloring_wp.out
867:       requires: !single

869:    test:
870:       suffix: fieldsplit_2
871:       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
872:       requires: !single

874:    test:
875:       suffix: fieldsplit_3
876:       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
877:       requires: !single

879:    test:
880:       suffix: fieldsplit_4
881:       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
882:       requires: !single

884:    # HYPRE PtAP broken with complex numbers
885:    test:
886:       suffix: fieldsplit_hypre
887:       nsize: 2
888:       requires: hypre mumps !complex !defined(PETSC_HAVE_HYPRE_DEVICE)
889:       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

891:    test:
892:       suffix: fieldsplit_mumps
893:       nsize: 2
894:       requires: mumps
895:       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
896:       output_file: output/ex19_fieldsplit_5.out

898:    test:
899:       suffix: greedy_coloring
900:       nsize: 2
901:       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
902:       requires: !single

904:    # HYPRE PtAP broken with complex numbers
905:    test:
906:       suffix: hypre
907:       nsize: 2
908:       requires: hypre !complex !defined(PETSC_HAVE_HYPRE_DEVICE)
909:       args: -da_refine 3 -snes_monitor_short -pc_type hypre -ksp_norm_type unpreconditioned

911:    # ibcgs is broken when using device vectors
912:    test:
913:       suffix: ibcgs
914:       nsize: 2
915:       args: -ksp_type ibcgs -ksp_monitor_short -da_refine 2 -snes_view
916:       requires: !complex !single

918:    test:
919:       suffix: kaczmarz
920:       nsize: 2
921:       args: -pc_type kaczmarz -ksp_monitor_short -snes_monitor_short -snes_view
922:       requires: !single

924:    test:
925:       suffix: klu
926:       requires: suitesparse
927:       args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type klu
928:       output_file: output/ex19_superlu.out

930:    test:
931:       suffix: klu_2
932:       requires: suitesparse
933:       args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type klu -pc_factor_mat_ordering_type nd
934:       output_file: output/ex19_superlu.out

936:    test:
937:       suffix: klu_3
938:       requires: suitesparse
939:       args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type klu -mat_klu_use_btf 0
940:       output_file: output/ex19_superlu.out

942:    test:
943:       suffix: ml
944:       nsize: 2
945:       requires: ml
946:       args: -da_refine 3 -snes_monitor_short -pc_type ml

948:    test:
949:       suffix: ngmres_fas
950:       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
951:       requires: !single

953:    test:
954:       suffix: ngmres_fas_gssecant
955:       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
956:       requires: !single

958:    test:
959:       suffix: ngmres_fas_ms
960:       nsize: 2
961:       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
962:       requires: !single

964:    test:
965:       suffix: ngmres_nasm
966:       nsize: 4
967:       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
968:       requires: !single

970:    test:
971:       suffix: ngs
972:       args: -snes_type ngs -snes_view -snes_monitor -snes_rtol 1e-4
973:       requires: !single

975:    test:
976:       suffix: ngs_fd
977:       args: -snes_type ngs -snes_ngs_secant -snes_view -snes_monitor -snes_rtol 1e-4
978:       requires: !single

980:    test:
981:       suffix: parms
982:       nsize: 2
983:       requires: parms
984:       args: -pc_type parms -ksp_monitor_short -snes_view

986:    test:
987:       suffix: superlu
988:       requires: superlu
989:       args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type superlu

991:    test:
992:       suffix: superlu_sell
993:       requires: superlu
994:       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
995:       output_file: output/ex19_superlu.out

997:    test:
998:       suffix: superlu_dist
999:       requires: superlu_dist
1000:       args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type superlu_dist
1001:       output_file: output/ex19_superlu.out

1003:    test:
1004:       suffix: superlu_dist_2
1005:       nsize: 2
1006:       requires: superlu_dist
1007:       args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type superlu_dist
1008:       output_file: output/ex19_superlu.out

1010:    test:
1011:       suffix: superlu_equil
1012:       requires: superlu
1013:       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

1015:    test:
1016:       suffix: superlu_equil_sell
1017:       requires: superlu
1018:       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
1019:       output_file: output/ex19_superlu_equil.out

1021:    test:
1022:       suffix: tcqmr
1023:       args: -da_refine 1 -ksp_monitor_short -ksp_type tcqmr
1024:       requires: !single

1026:    test:
1027:       suffix: tfqmr
1028:       args: -da_refine 1 -ksp_monitor_short -ksp_type tfqmr
1029:       requires: !single

1031:    test:
1032:       suffix: umfpack
1033:       requires: suitesparse
1034:       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

1036:    test:
1037:       suffix: tut_1
1038:       nsize: 4
1039:       requires: !single
1040:       args: -da_refine 5 -snes_monitor -ksp_monitor -snes_view

1042:    test:
1043:       suffix: tut_2
1044:       nsize: 4
1045:       requires: !single
1046:       args: -da_refine 5 -snes_monitor -ksp_monitor -snes_view -pc_type mg

1048:    # HYPRE PtAP broken with complex numbers
1049:    test:
1050:       suffix: tut_3
1051:       nsize: 4
1052:       requires: hypre !single !complex !defined(PETSC_HAVE_HYPRE_DEVICE)
1053:       args: -da_refine 5 -snes_monitor -ksp_monitor -snes_view -pc_type hypre

1055:    test:
1056:       suffix: tut_8
1057:       nsize: 4
1058:       requires: ml !single
1059:       args: -da_refine 5 -snes_monitor -ksp_monitor -snes_view -pc_type ml

1061:    test:
1062:       suffix: tut_4
1063:       nsize: 1
1064:       requires: !single
1065:       args: -da_refine 5 -log_view
1066:       filter: head -n 2
1067:       filter_output: head -n 2

1069:    test:
1070:       suffix: tut_5
1071:       nsize: 1
1072:       requires: !single
1073:       args: -da_refine 5 -log_view -pc_type mg
1074:       filter: head -n 2
1075:       filter_output: head -n 2

1077:    test:
1078:       suffix: tut_6
1079:       nsize: 4
1080:       requires: !single
1081:       args: -da_refine 5 -log_view
1082:       filter: head -n 2
1083:       filter_output: head -n 2

1085:    test:
1086:       suffix: tut_7
1087:       nsize: 4
1088:       requires: !single
1089:       args: -da_refine 5 -log_view -pc_type mg
1090:       filter: head -n 2
1091:       filter_output: head -n 2

1093:    test:
1094:       suffix: cuda_1
1095:       nsize: 1
1096:       requires: cuda
1097:       args: -snes_monitor -dm_mat_type seqaijcusparse -dm_vec_type seqcuda -pc_type gamg -ksp_monitor -mg_levels_ksp_max_it 3

1099:    test:
1100:       suffix: cuda_2
1101:       nsize: 3
1102:       requires: cuda !single
1103:       args: -snes_monitor -dm_mat_type mpiaijcusparse -dm_vec_type mpicuda -pc_type gamg -ksp_monitor  -mg_levels_ksp_max_it 3

1105:    test:
1106:       suffix: seqbaijmkl
1107:       nsize: 1
1108:       requires: defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE)
1109:       args: -dm_mat_type baij -snes_monitor -ksp_monitor -snes_view

1111:    test:
1112:       suffix: mpibaijmkl
1113:       nsize: 2
1114:       requires:  defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE)
1115:       args: -dm_mat_type baij -snes_monitor -ksp_monitor -snes_view

1117:    test:
1118:      suffix: cpardiso
1119:      nsize: 4
1120:      requires: mkl_cpardiso
1121:      args: -pc_type lu -pc_factor_mat_solver_type mkl_cpardiso -ksp_monitor

1123:    test:
1124:      suffix: logviewmemory
1125:      requires: defined(PETSC_USE_LOG) !defined(PETSCTEST_VALGRIND)
1126:      args: -log_view -log_view_memory -da_refine 4
1127:      filter: grep MatFDColorSetUp | wc -w | xargs  -I % sh -c "expr % \> 21"

1129:    test:
1130:      suffix: fs
1131:      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

1133:    test:
1134:      suffix: asm_matconvert
1135:      args: -mat_type aij -pc_type asm -pc_asm_sub_mat_type dense -snes_view

1137:    test:
1138:       suffix: euclid
1139:       nsize: 2
1140:       requires: hypre !single !complex !defined(PETSC_HAVE_HYPRE_MIXEDINT) !defined(PETSC_HAVE_HYPRE_DEVICE)
1141:       args: -da_refine 2 -ksp_monitor -snes_monitor -snes_view -pc_type hypre -pc_hypre_type euclid

1143:    test:
1144:       suffix: euclid_bj
1145:       nsize: 2
1146:       requires: hypre !single !complex !defined(PETSC_HAVE_HYPRE_MIXEDINT) !defined(PETSC_HAVE_HYPRE_DEVICE)
1147:       args: -da_refine 2 -ksp_monitor -snes_monitor -snes_view -pc_type hypre -pc_hypre_type euclid -pc_hypre_euclid_bj

1149:    test:
1150:       suffix: euclid_droptolerance
1151:       nsize: 1
1152:       requires: hypre !single !complex !defined(PETSC_HAVE_HYPRE_MIXEDINT) !defined(PETSC_HAVE_HYPRE_DEVICE)
1153:       args: -da_refine 2 -ksp_monitor -snes_monitor -snes_view -pc_type hypre -pc_hypre_type euclid -pc_hypre_euclid_droptolerance .1

1155: TEST*/