Actual source code: ex19.c

petsc-3.3-p7 2013-05-11
  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 <lid>, where <lid> = dimensionless velocity of lid\n\
  7:   -grashof <gr>, where <gr> = dimensionless temperature gradent\n\
  8:   -prandtl <pr>, where <pr> = dimensionless thermal/momentum diffusity ratio\n\
  9:   -contours : draw contour plots of solution\n\n";
 10: /*
 11:       See src/ksp/ksp/examples/tutorials/ex45.c
 12: */

 14: /*T
 15:    Concepts: SNES^solving a system of nonlinear equations (parallel multicomponent example);
 16:    Concepts: DMDA^using distributed arrays;
 17:    Concepts: multicomponent
 18:    Processors: n
 19: T*/
 20: /*F-----------------------------------------------------------------------

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

 54:   ------------------------------------------------------------------------F*/

 56: /* 
 57:    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
 58:    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
 59:    file automatically includes:
 60:      petscsys.h       - base PETSc routines   petscvec.h - vectors
 61:      petscmat.h - matrices
 62:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 63:      petscviewer.h - viewers               petscpc.h  - preconditioners
 64:      petscksp.h   - linear solvers 
 65: */
 66: #include <petscsnes.h>
 67: #include <petscdmda.h>

 69: /* 
 70:    User-defined routines and data structures
 71: */
 72: typedef struct {
 73:   PetscScalar u,v,omega,temp;
 74: } Field;

 76: PetscErrorCode FormFunctionLocal(DMDALocalInfo*,Field**,Field**,void*);
 77: PetscErrorCode FormFunction(SNES,Vec,Vec,void*);

 79: typedef struct {
 80:    PassiveReal  lidvelocity,prandtl,grashof;  /* physical parameters */
 81:    PetscBool    draw_contours;                /* flag - 1 indicates drawing contours */
 82: } AppCtx;

 84: PetscErrorCode FormInitialGuess(AppCtx*,DM,Vec);
 85: extern PetscErrorCode NonlinearGS(SNES,Vec,Vec,void*);

 89: int main(int argc,char **argv)
 90: {
 91:   AppCtx         user;                /* user-defined work context */
 92:   PetscInt       mx,my,its;
 94:   MPI_Comm       comm;
 95:   SNES           snes;
 96:   DM             da;
 97:   Vec            x;

 99:   PetscInitialize(&argc,&argv,(char *)0,help);if (ierr) return(1);
100:   comm = PETSC_COMM_WORLD;

102:   SNESCreate(comm,&snes);
103: 
104:   /*
105:       Create distributed array object to manage parallel grid and vectors
106:       for principal unknowns (x) and governing residuals (f)
107:   */
108:   DMDACreate2d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_STENCIL_STAR,-4,-4,PETSC_DECIDE,PETSC_DECIDE,4,1,0,0,&da);
109:   SNESSetDM(snes,(DM)da);
110:   SNESSetGS(snes, NonlinearGS, (void *)&user);

112:   DMDAGetInfo(da,0,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
113:                    PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
114:   /* 
115:      Problem parameters (velocity of lid, prandtl, and grashof numbers)
116:   */
117:   user.lidvelocity = 1.0/(mx*my);
118:   user.prandtl     = 1.0;
119:   user.grashof     = 1.0;
120:   PetscOptionsGetReal(PETSC_NULL,"-lidvelocity",&user.lidvelocity,PETSC_NULL);
121:   PetscOptionsGetReal(PETSC_NULL,"-prandtl",&user.prandtl,PETSC_NULL);
122:   PetscOptionsGetReal(PETSC_NULL,"-grashof",&user.grashof,PETSC_NULL);
123:   PetscOptionsHasName(PETSC_NULL,"-contours",&user.draw_contours);

125:   DMDASetFieldName(da,0,"x_velocity");
126:   DMDASetFieldName(da,1,"y_velocity");
127:   DMDASetFieldName(da,2,"Omega");
128:   DMDASetFieldName(da,3,"temperature");

130:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131:      Create user context, set problem data, create vector data structures.
132:      Also, compute the initial guess.
133:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

135:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136:      Create nonlinear solver context

138:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139:   DMSetApplicationContext(da,&user);
140:   DMDASetLocalFunction(da,(DMDALocalFunction1)FormFunctionLocal);
141:   SNESSetFromOptions(snes);
142:   PetscPrintf(comm,"lid velocity = %G, prandtl # = %G, grashof # = %G\n",user.lidvelocity,user.prandtl,user.grashof);


145:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146:      Solve the nonlinear system
147:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
148:   DMCreateGlobalVector(da,&x);
149:   FormInitialGuess(&user,da,x);

151:   SNESSolve(snes,PETSC_NULL,x);

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

156:   /*
157:      Visualize solution
158:   */
159:   if (user.draw_contours) {
160:     VecView(x,PETSC_VIEWER_DRAW_WORLD);
161:   }

163:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
164:      Free work space.  All PETSc objects should be destroyed when they
165:      are no longer needed.
166:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
167:   VecDestroy(&x);
168:   DMDestroy(&da);
169:   SNESDestroy(&snes);
170:   PetscFinalize();
171:   return 0;
172: }

174: /* ------------------------------------------------------------------- */


179: /* 
180:    FormInitialGuess - Forms initial approximation.

182:    Input Parameters:
183:    user - user-defined application context
184:    X - vector

186:    Output Parameter:
187:    X - vector
188:  */
189: PetscErrorCode FormInitialGuess(AppCtx *user,DM da,Vec X)
190: {
191:   PetscInt       i,j,mx,xs,ys,xm,ym;
193:   PetscReal      grashof,dx;
194:   Field          **x;

196:   grashof = user->grashof;

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

201:   /*
202:      Get local grid boundaries (for 2-dimensional DMDA):
203:        xs, ys   - starting grid indices (no ghost points)
204:        xm, ym   - widths of local grid (no ghost points)
205:   */
206:   DMDAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);

208:   /*
209:      Get a pointer to vector data.
210:        - For default PETSc vectors, VecGetArray() returns a pointer to
211:          the data array.  Otherwise, the routine is implementation dependent.
212:        - You MUST call VecRestoreArray() when you no longer need access to
213:          the array.
214:   */
215:   DMDAVecGetArray(da,X,&x);

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

230:   /*
231:      Restore vector
232:   */
233:   DMDAVecRestoreArray(da,X,&x);
234:   return 0;
235: }
236: 
239: PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,Field **x,Field **f,void *ptr)
240:  {
241:   AppCtx         *user = (AppCtx*)ptr;
243:   PetscInt       xints,xinte,yints,yinte,i,j;
244:   PetscReal      hx,hy,dhx,dhy,hxdhy,hydhx;
245:   PetscReal      grashof,prandtl,lid;
246:   PetscScalar    u,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;

249:   grashof = user->grashof;
250:   prandtl = user->prandtl;
251:   lid     = user->lidvelocity;

253:   /* 
254:      Define mesh intervals ratios for uniform grid.

256:      Note: FD formulae below are normalized by multiplying through by
257:      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;

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

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:   }

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:   }

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:   }

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:   }

319:   /* Compute over the interior points */
320:   for (j=yints; j<yinte; j++) {
321:     for (i=xints; i<xinte; i++) {

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);

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;

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;

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 +
348:                         (vxp*(u - x[j][i-1].omega) +
349:                           vxm*(x[j][i+1].omega - u)) * hy +
350:                         (vyp*(u - x[j-1][i].omega) +
351:                           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 * (
359:                         (vxp*(u - x[j][i-1].temp) +
360:                           vxm*(x[j][i+1].temp - u)) * hy +
361:                         (vyp*(u - x[j-1][i].temp) +
362:                                  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: }


376: PetscErrorCode FormFunction(SNES snes, Vec X, Vec F, void *user)
377: {
378:   DMDALocalInfo  info;
379:   Field          **u,**fu;
381:   Vec            localX;
382:   DM             da;

385:   SNESGetDM(snes,(DM*)&da);
386:   DMGetLocalVector(da,&localX);
387:   /*
388:      Scatter ghost points to local vector, using the 2-step process
389:         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
390:   */
391:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
392:   DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
393:   DMDAGetLocalInfo(da,&info);
394:   DMDAVecGetArray(da,localX,&u);
395:   DMDAVecGetArray(da,F,&fu);
396:   FormFunctionLocal(&info,u,fu,user);
397:   DMDAVecRestoreArray(da,localX,&u);
398:   DMDAVecRestoreArray(da,F,&fu);
399:   DMRestoreLocalVector(da,&localX);
400:   return(0);
401: }

405: PetscErrorCode NonlinearGS(SNES snes, Vec X, Vec B, void *ctx)
406: {
407:   DMDALocalInfo  info;
408:   Field          **x,**b;
410:   Vec            localX, localB;
411:   DM             da;
412:   PetscInt       xints,xinte,yints,yinte,i,j,k,l;
413:   PetscInt       n_pointwise = 50;
414:   PetscInt       n_sweeps = 3;
415:   PetscReal      hx,hy,dhx,dhy,hxdhy,hydhx;
416:   PetscReal      grashof,prandtl,lid;
417:   PetscScalar    u,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;
418:   PetscScalar    fu, fv, fomega, ftemp;
419:   PetscScalar    dfudu;
420:   PetscScalar    dfvdv;
421:   PetscScalar    dfodu, dfodv, dfodo;
422:   PetscScalar    dftdu, dftdv, dftdt;
423:   PetscScalar    yu, yv, yo, yt;
424:   PetscScalar    bjiu, bjiv, bjiomega, bjitemp;
425:   PetscBool      ptconverged;
426:   PetscScalar    pfnorm, pfnorm0;
427:   AppCtx         *user = (AppCtx*)ctx;

430:   grashof = user->grashof;
431:   prandtl = user->prandtl;
432:   lid     = user->lidvelocity;
433:   SNESGetDM(snes,(DM*)&da);
434:   DMGetLocalVector(da,&localX);
435:   if (B) {
436:     DMGetLocalVector(da,&localB);
437:   }
438:   /*
439:      Scatter ghost points to local vector, using the 2-step process
440:         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
441:   */
442:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
443:   DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
444:   if (B) {
445:     DMGlobalToLocalBegin(da,B,INSERT_VALUES,localB);
446:     DMGlobalToLocalEnd(da,B,INSERT_VALUES,localB);
447:   }
448:   DMDAGetLocalInfo(da,&info);
449:   DMDAVecGetArray(da,localX,&x);
450:   if (B) {
451:     DMDAVecGetArray(da,localB,&b);
452:   }
453:   /* looks like a combination of the formfunction / formjacobian routines */
454:   dhx = (PetscReal)(info.mx-1);  dhy = (PetscReal)(info.my-1);
455:   hx = 1.0/dhx;                   hy = 1.0/dhy;
456:   hxdhy = hx*dhy;                 hydhx = hy*dhx;

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

460:   /* Set the boundary conditions on the momentum equations */
461:   /* Test whether we are on the bottom edge of the global array */
462:   if (yints == 0) {
463:     j = 0;
464:     yints = yints + 1;
465:     /* bottom edge */
466:     for (i=info.xs; i<info.xs+info.xm; i++) {

468:       if (B) {
469:         bjiu = b[j][i].u;
470:         bjiv = b[j][i].v;
471:       } else {
472:         bjiu = 0.0;
473:         bjiv = 0.0;
474:       }
475:       x[j][i].u     = 0.0 + bjiu;
476:       x[j][i].v     = 0.0 + bjiv;
477:     }
478:   }

480:   /* Test whether we are on the top edge of the global array */
481:   if (yinte == info.my) {
482:     j = info.my - 1;
483:     yinte = yinte - 1;
484:     /* top edge */
485:     for (i=info.xs; i<info.xs+info.xm; i++) {
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     = lid + bjiu;
494:       x[j][i].v     = bjiv;
495:     }
496:   }

498:   /* Test whether we are on the left edge of the global array */
499:   if (xints == 0) {
500:     i = 0;
501:     xints = xints + 1;
502:     /* left edge */
503:     for (j=info.ys; j<info.ys+info.ym; j++) {
504:       if (B) {
505:         bjiu = b[j][i].u;
506:         bjiv = b[j][i].v;
507:       } else {
508:         bjiu = 0.0;
509:         bjiv = 0.0;
510:       }
511:       x[j][i].u     = 0.0 + bjiu;
512:       x[j][i].v     = 0.0 + bjiv;
513:     }
514:   }

516:   /* Test whether we are on the right edge of the global array */
517:   if (xinte == info.mx) {
518:     i = info.mx - 1;
519:     xinte = xinte - 1;
520:     /* right edge */
521:     for (j=info.ys; j<info.ys+info.ym; j++) {
522:       if (B) {
523:         bjiu = b[j][i].u;
524:         bjiv = b[j][i].v;
525:       } else {
526:         bjiu = 0.0;
527:         bjiv = 0.0;
528:       }
529:       x[j][i].u     = 0.0 + bjiu;
530:       x[j][i].v     = 0.0 + bjiv;
531:     }
532:   }

534:   for (k=0; k < n_sweeps; k++) {
535:     for (j=info.ys; j<info.ys + info.ym; j++) {
536:       for (i=info.xs; i<info.xs + info.xm; i++) {
537:         ptconverged = PETSC_FALSE;
538:         pfnorm0 = 0.0;
539:         pfnorm = 0.0;
540:         fu = 0.0;
541:         fv = 0.0;
542:         fomega = 0.0;
543:         ftemp = 0.0;
544:         for (l = 0; l < n_pointwise && !ptconverged; l++) {
545:           if (B) {
546:             bjiu = b[j][i].u;
547:             bjiv = b[j][i].v;
548:             bjiomega = b[j][i].omega;
549:             bjitemp = b[j][i].temp;
550:           } else {
551:             bjiu = 0.0;
552:             bjiv = 0.0;
553:             bjiomega = 0.0;
554:             bjitemp = 0.0;
555:           }

557:           if (i != 0 && i != info.mx - 1 && j != 0 && j != info.my-1) {
558:             /* U velocity */
559:             u          = x[j][i].u;
560:             uxx        = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
561:             uyy        = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
562:             fu    = uxx + uyy - .5*(x[j+1][i].omega-x[j-1][i].omega)*hx - bjiu;
563:             dfudu = 2.0*(hydhx + hxdhy);
564:             /* V velocity */
565:             u          = x[j][i].v;
566:             uxx        = (2.0*u - x[j][i-1].v - x[j][i+1].v)*hydhx;
567:             uyy        = (2.0*u - x[j-1][i].v - x[j+1][i].v)*hxdhy;
568:             fv    = uxx + uyy + .5*(x[j][i+1].omega-x[j][i-1].omega)*hy - bjiv;
569:             dfvdv = 2.0*(hydhx + hxdhy);
570:             /*
571:              convective coefficients for upwinding
572:              */
573:             vx = x[j][i].u; avx = PetscAbsScalar(vx);
574:             vxp = .5*(vx+avx); vxm = .5*(vx-avx);
575:             vy = x[j][i].v; avy = PetscAbsScalar(vy);
576:             vyp = .5*(vy+avy); vym = .5*(vy-avy);
577:             /* Omega */
578:             u          = x[j][i].omega;
579:             uxx        = (2.0*u - x[j][i-1].omega - x[j][i+1].omega)*hydhx;
580:             uyy        = (2.0*u - x[j-1][i].omega - x[j+1][i].omega)*hxdhy;
581:             fomega = uxx + uyy +
582:               (vxp*(u - x[j][i-1].omega) +
583:                vxm*(x[j][i+1].omega - u)) * hy +
584:               (vyp*(u - x[j-1][i].omega) +
585:                vym*(x[j+1][i].omega - u)) * hx -
586:               .5 * grashof * (x[j][i+1].temp - x[j][i-1].temp) * hy - bjiomega;
587:             /* convective coefficient derivatives */
588:             dfodo = 2.0*(hydhx + hxdhy) + ((vxp - vxm)*hy + (vyp - vym)*hx);
589:             if (PetscRealPart(vx) > 0.0) {
590:               dfodu = (u - x[j][i-1].omega)*hy;
591:             } else {
592:               dfodu = (x[j][i+1].omega - u)*hy;
593:             }
594:             if (PetscRealPart(vy) > 0.0) {
595:               dfodv = (u - x[j-1][i].omega)*hx;
596:             } else {
597:               dfodv = (x[j+1][i].omega - u)*hx;
598:             }
599:             /* Temperature */
600:             u             = x[j][i].temp;
601:             uxx           = (2.0*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx;
602:             uyy           = (2.0*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy;
603:             ftemp =  uxx + uyy  + prandtl * (
604:               (vxp*(u - x[j][i-1].temp) +
605:                vxm*(x[j][i+1].temp - u)) * hy +
606:               (vyp*(u - x[j-1][i].temp) +
607:                vym*(x[j+1][i].temp - u)) * hx) - bjitemp;
608:             dftdt = 2.0*(hydhx + hxdhy) + prandtl*((vxp - vxm)*hy + (vyp - vym)*hx);
609:             if (PetscRealPart(vx) > 0.0) {
610:               dftdu = prandtl*(u - x[j][i-1].temp)*hy;
611:             } else {
612:               dftdu = prandtl*(x[j][i+1].temp - u)*hy;
613:             }
614:             if (PetscRealPart(vy) > 0.0) {
615:               dftdv = prandtl*(u - x[j-1][i].temp)*hx;
616:             } else {
617:               dftdv = prandtl*(x[j+1][i].temp - u)*hx;
618:             }
619:             /* invert the system:
620:              [ dfu / du     0        0        0    ][yu] = [fu]
621:              [     0    dfv / dv     0        0    ][yv]   [fv]
622:              [ dfo / du dfo / dv dfo / do     0    ][yo]   [fo]
623:              [ dft / du dft / dv     0    dft / dt ][yt]   [ft]
624:              by simple back-substitution
625:            */
626:             yu = fu / dfudu;
627:             yv = fv / dfvdv;
628:             yo = fomega / dfodo;
629:             yt = ftemp / dftdt;
630:             yo = (fomega - (dfodu*yu + dfodv*yv)) / dfodo;
631:             yt = (ftemp - (dftdu*yu + dftdv*yv)) / dftdt;

633:             x[j][i].u = x[j][i].u - yu;
634:             x[j][i].v = x[j][i].v - yv;
635:             x[j][i].temp  = x[j][i].temp - yt;
636:             x[j][i].omega = x[j][i].omega - yo;
637:           }
638:           if (i == 0) {
639:             fomega = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx - bjiomega;
640:             ftemp  = x[j][i].temp - bjitemp;
641:             x[j][i].omega = x[j][i].omega - fomega;
642:             x[j][i].temp  = x[j][i].temp - ftemp;
643:           }
644:           if (i == info.mx - 1) {
645:             fomega = x[j][i].omega - (x[j][i].v - x[j][i-1].v)*dhx - bjiomega;
646:             ftemp  = x[j][i].temp - (PetscReal)(grashof>0) - bjitemp;
647:             x[j][i].omega = x[j][i].omega - fomega;
648:             x[j][i].temp  = x[j][i].temp - ftemp;
649:           }
650:           if (j == 0) {
651:             fomega = x[j][i].omega + (x[j+1][i].u - x[j][i].u)*dhy - bjiomega;
652:             ftemp  = x[j][i].temp-x[j+1][i].temp - bjitemp;
653:             x[j][i].omega = x[j][i].omega - fomega;
654:             x[j][i].temp  = x[j][i].temp - ftemp;
655:           }
656:           if (j == info.my - 1) {
657:             fomega = x[j][i].omega + (x[j][i].u - x[j-1][i].u)*dhy - bjiomega;
658:             ftemp  = x[j][i].temp-x[j-1][i].temp - bjitemp;
659:             x[j][i].omega = x[j][i].omega - fomega;
660:             x[j][i].temp  = x[j][i].temp - ftemp;
661:           }
662:           pfnorm = fu*fu + fv*fv + fomega*fomega + ftemp*ftemp;
663:           pfnorm = PetscSqrtScalar(pfnorm);
664:           if (l == 0) pfnorm0 = pfnorm;
665:           if (1e-15*PetscRealPart(pfnorm0) > PetscRealPart(pfnorm)) ptconverged = PETSC_TRUE;
666:         }
667:       }
668:     }
669:   }
670:   DMDAVecRestoreArray(da,localX,&x);
671:   if (B) {
672:     DMDAVecRestoreArray(da,localB,&b);
673:   }
674:   DMLocalToGlobalBegin(da,localX,INSERT_VALUES,X);
675:   DMLocalToGlobalEnd(da,localX,INSERT_VALUES,X);
676:   PetscLogFlops(n_sweeps*n_pointwise*(84.0 + 41)*info.ym*info.xm);
677:   if (B) {
678:     DMLocalToGlobalBegin(da,localB,INSERT_VALUES,B);
679:     DMLocalToGlobalEnd(da,localB,INSERT_VALUES,B);
680:   }
681:   DMRestoreLocalVector(da,&localX);
682:   if (B) {
683:     DMRestoreLocalVector(da,&localB);
684:   }
685:   return(0);
686: }