Actual source code: ex19.c

petsc-3.10.5 2019-03-28
Report Typos and Errors

  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/examples/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).

 60: /*
 61:    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
 62:    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
 63:    file automatically includes:
 64:      petscsys.h       - base PETSc routines   petscvec.h - vectors
 65:      petscmat.h - matrices
 66:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 67:      petscviewer.h - viewers               petscpc.h  - preconditioners
 68:      petscksp.h   - linear solvers
 69: */
 70: #if defined(PETSC_APPLE_FRAMEWORK)
 71: #import <PETSc/petscsnes.h>
 72: #import <PETSc/petscdmda.h>
 73: #else
 74:  #include <petscsnes.h>
 75:  #include <petscdm.h>
 76:  #include <petscdmda.h>
 77: #endif

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

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

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

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

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

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

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

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

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

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

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

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

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

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


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

161:   SNESSolve(snes,NULL,x);

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

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

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

184: /* ------------------------------------------------------------------- */

186: /*
187:    FormInitialGuess - Forms initial approximation.

189:    Input Parameters:
190:    user - user-defined application context
191:    X - vector

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

204:   grashof = user->grashof;

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

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

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

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

238:   /*
239:      Restore vector
240:   */
241:   DMDAVecRestoreArray(da,X,&x);
242:   return(0);
243: }

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

255:   grashof = user->grashof;
256:   prandtl = user->prandtl;
257:   lid     = user->lidvelocity;

259:   /*
260:      Define mesh intervals ratios for uniform grid.

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


266:   */
267:   dhx   = (PetscReal)(info->mx-1);  dhy = (PetscReal)(info->my-1);
268:   hx    = 1.0/dhx;                   hy = 1.0/dhy;
269:   hxdhy = hx*dhy;                 hydhx = hy*dhx;

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

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

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

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

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

325:   /* Compute over the interior points */
326:   for (j=yints; j<yinte; j++) {
327:     for (i=xints; i<xinte; i++) {

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

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

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

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

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

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

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:   DMDAVecGetArray(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:         for (l = 0; l < max_its && !ptconverged; l++) {
512:           if (B) {
513:             bjiu     = b[j][i].u;
514:             bjiv     = b[j][i].v;
515:             bjiomega = b[j][i].omega;
516:             bjitemp  = b[j][i].temp;
517:           } else {
518:             bjiu     = 0.0;
519:             bjiv     = 0.0;
520:             bjiomega = 0.0;
521:             bjitemp  = 0.0;
522:           }

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

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

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

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

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

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


650: /*TEST

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

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

663:    test:
664:       suffix: 11
665:       nsize: 4
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 2 -da_refine 4 -ksp_type fgmres

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

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

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

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

694:    test:
695:       suffix: 17
696:       args: -snes_monitor_short -ksp_pc_side right
697:       requires: !single

699:    test:
700:       suffix: 18
701:       args: -ksp_monitor_snes_lg -ksp_pc_side right
702:       requires: x !single

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

710:    test:
711:       suffix: 2_bcols1
712:       nsize: 4
713:       args: -da_refine 3 -snes_converged_reason -pc_type mg -mat_fd_type ds -mat_fd_coloring_bcols 1> ex19_1.tmp 2>&1
714:       output_file: output/ex19_2.out
715:       requires: !single

717:    test:
718:       suffix: 3
719:       nsize: 4
720:       requires: mumps
721:       args: -da_refine 3 -snes_monitor_short -pc_type redundant -dm_mat_type mpiaij -redundant_ksp_type preonly -redundant_pc_factor_mat_solver_type mumps -pc_redundant_number 2

723:    test:
724:       suffix: 4
725:       nsize: 12
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 5
728:       output_file: output/ex19_3.out

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

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

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

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

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

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

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

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

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

783:    test:
784:       suffix: cgs
785:       args: -da_refine 1 -ksp_monitor_short -ksp_type cgs
786:       requires: !single

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

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

798:    test:
799:       suffix: composite_fieldsplit_bjacobi_2
800:       nsize: 4
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

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

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

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

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

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

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

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

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

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

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

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

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

871:    test:
872:       suffix: fieldsplit_2
873:       args: -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_block_size 4 -pc_fieldsplit_type additive -pc_fieldsplit_0_fields 0,1,2 -pc_fieldsplit_1_fields 3 -snes_monitor_short -ksp_monitor_short
874:       requires: !single

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

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

886:    test:
887:       suffix: fieldsplit_hypre
888:       nsize: 2
889:       requires: hypre mumps
890:       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

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

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

905:    test:
906:       suffix: hypre
907:       nsize: 2
908:       requires: hypre
909:       args: -da_refine 3 -snes_monitor_short -pc_type hypre

911:    test:
912:       suffix: ibcgs
913:       nsize: 2
914:       args: -ksp_type ibcgs -ksp_monitor_short -da_refine 2 -snes_view
915:       requires: !complex !single

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1014:    test:
1015:       suffix: superlu_equil_sell
1016:       requires: superlu
1017:       args: -da_grid_x 20 -da_grid_y 20 -{snes,ksp}_monitor_short -pc_type lu -pc_factor_mat_solver_type superlu -mat_superlu_equil -dm_mat_type sell -pc_factor_mat_ordering_type natural
1018:       output_file: output/ex19_superlu_equil.out

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

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

1030:    test:
1031:       suffix: umfpack
1032:       requires: suitesparse
1033:       args: -da_refine 2 -pc_type lu -pc_factor_mat_solver_type umfpack -snes_view -snes_monitor_short -ksp_monitor_short

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

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

1047:    test:
1048:       suffix: tut_3
1049:       nsize: 4
1050:       requires: hypre !single
1051:       args: -da_refine 5 -snes_monitor -ksp_monitor -snes_view -pc_type hypre

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

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

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

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

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

1091: TEST*/