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

 11: /*T
 12:    Concepts: SNES^solving a system of nonlinear equations (parallel multicomponent example);
 13:    Concepts: DA^using distributed arrays;
 14:    Concepts: multicomponent
 15:    Processors: n
 16: T*/
 17: /* ------------------------------------------------------------------------

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

 22:     This problem is modeled by the partial differential equation system
 23:   
 24:         - Lap(U) - Grad_y(Omega) = 0
 25:         - Lap(V) + Grad_x(Omega) = 0
 26:         - Lap(Omega) + Div([U*Omega,V*Omega]) - GR*Grad_x(T) = 0
 27:         - Lap(T) + PR*Div([U*T,V*T]) = 0

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

 32:     No-slip, rigid-wall Dirichlet conditions are used for [U,V].
 33:     Dirichlet conditions are used for Omega, based on the definition of
 34:     vorticity: Omega = - Grad_y(U) + Grad_x(V), where along each
 35:     constant coordinate boundary, the tangential derivative is zero.
 36:     Dirichlet conditions are used for T on the left and right walls,
 37:     and insulation homogeneous Neumann conditions are used for T on
 38:     the top and bottom walls. 

 40:     A finite difference approximation with the usual 5-point stencil 
 41:     is used to discretize the boundary value problem to obtain a 
 42:     nonlinear system of equations.  Upwinding is used for the divergence
 43:     (convective) terms and central for the gradient (source) terms.
 44:     
 45:     The Jacobian can be either
 46:       * formed via finite differencing using coloring (the default), or
 47:       * applied matrix-free via the option -snes_mf 
 48:         (for larger grid problems this variant may not converge 
 49:         without a preconditioner due to ill-conditioning).

 51:   ------------------------------------------------------------------------- */

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

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


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

 86: int main(int argc,char **argv)
 87: {
 88:   DMMG           *dmmg;               /* multilevel grid structure */
 89:   AppCtx         user;                /* user-defined work context */
 90:   PetscInt       mx,my,its,nlevels=2;
 92:   MPI_Comm       comm;
 93:   SNES           snes;
 94:   DA             da;

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

 99:   PetscOptionsGetInt(PETSC_NULL,"-nlevels",&nlevels,PETSC_NULL);
100:   PreLoadBegin(PETSC_TRUE,"SetUp");
101:     DMMGCreate(comm,nlevels,&user,&dmmg);

103:     /*
104:       Create distributed array multigrid object (DMMG) to manage parallel grid and vectors
105:       for principal unknowns (x) and governing residuals (f)
106:     */
107:     DACreate2d(comm,DA_NONPERIODIC,DA_STENCIL_STAR,-4,-4,PETSC_DECIDE,PETSC_DECIDE,4,1,0,0,&da);
108:     DMMGSetDM(dmmg,(DM)da);
109:     DADestroy(da);

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

124:     DASetFieldName(DMMGGetDA(dmmg),0,"x-velocity");
125:     DASetFieldName(DMMGGetDA(dmmg),1,"y-velocity");
126:     DASetFieldName(DMMGGetDA(dmmg),2,"Omega");
127:     DASetFieldName(DMMGGetDA(dmmg),3,"temperature");

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

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

137:        Process adiC(36): FormFunctionLocal FormFunctionLocali
138:        Process blockadiC(4): FormFunctionLocali4
139:        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140:     DMMGSetSNESLocal(dmmg,FormFunctionLocal,0,ad_FormFunctionLocal,admf_FormFunctionLocal);
141:     DMMGSetSNESLocali(dmmg,FormFunctionLocali,0,admf_FormFunctionLocali);
142:     DMMGSetSNESLocalib(dmmg,FormFunctionLocali4,0,admfb_FormFunctionLocali4);

144:     PetscPrintf(comm,"lid velocity = %G, prandtl # = %G, grashof # = %G\n",
145:                        user.lidvelocity,user.prandtl,user.grashof);


148:     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149:        Solve the nonlinear system
150:        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
151:     DMMGSetInitialGuess(dmmg,FormInitialGuess);

153:   PreLoadStage("Solve");
154:     DMMGSolve(dmmg);

156:     snes = DMMGGetSNES(dmmg);
157:     SNESGetIterationNumber(snes,&its);
158:     PetscPrintf(comm,"Number of Newton iterations = %D\n", its);

160:     /*
161:       Visualize solution
162:     */

164:     if (user.draw_contours) {
165:       VecView(DMMGGetx(dmmg),PETSC_VIEWER_DRAW_WORLD);
166:     }

168:     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
169:        Free work space.  All PETSc objects should be destroyed when they
170:        are no longer needed.
171:        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

173:     DMMGDestroy(dmmg);
174:   PreLoadEnd();

176:   PetscFinalize();
177:   return 0;
178: }

180: /* ------------------------------------------------------------------- */


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

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

192:    Output Parameter:
193:    X - vector
194:  */
195: PetscErrorCode FormInitialGuess(DMMG dmmg,Vec X)
196: {
197:   AppCtx         *user = (AppCtx*)dmmg->user;
198:   DA             da = (DA)dmmg->dm;
199:   PetscInt       i,j,mx,xs,ys,xm,ym;
201:   PetscReal      grashof,dx;
202:   Field          **x;

204:   grashof = user->grashof;

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

209:   /*
210:      Get local grid boundaries (for 2-dimensional DA):
211:        xs, ys   - starting grid indices (no ghost points)
212:        xm, ym   - widths of local grid (no ghost points)
213:   */
214:   DAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_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:   DAVecGetArray(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:   DAVecRestoreArray(da,X,&x);
242:   return 0;
243: }
244: 
245: PetscErrorCode FormFunctionLocal(DALocalInfo *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.

265:      
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 +
354:                         (vxp*(u - x[j][i-1].omega) +
355:                           vxm*(x[j][i+1].omega - u)) * hy +
356:                         (vyp*(u - x[j-1][i].omega) +
357:                           vym*(x[j+1][i].omega - u)) * hx -
358:                         .5 * grashof * (x[j][i+1].temp - x[j][i-1].temp) * hy;

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

372:   /*
373:      Flop count (multiply-adds are counted as 2 operations)
374:   */
375:   PetscLogFlops(84*info->ym*info->xm);
376:   return(0);
377: }

379: /*
380:     This function that evaluates the function for a single 
381:     degree of freedom. It is used by the -dmmg_fas solver
382: */
383: PetscErrorCode FormFunctionLocali(DALocalInfo *info,MatStencil *st,Field **x,PetscScalar *f,void *ptr)
384: {
385:   AppCtx      *user = (AppCtx*)ptr;
386:   PetscInt    i,j,c;
387:   PassiveReal hx,hy,dhx,dhy,hxdhy,hydhx;
388:   PassiveReal grashof,prandtl,lid;
389:   PetscScalar u,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;

392:   grashof = user->grashof;
393:   prandtl = user->prandtl;
394:   lid     = user->lidvelocity;

396:   /* 
397:      Define mesh intervals ratios for uniform grid.
398:      [Note: FD formulae below are normalized by multiplying through by
399:      local volume element to obtain coefficients O(1) in two dimensions.]
400:   */
401:   dhx = (PetscReal)(info->mx-1);     dhy = (PetscReal)(info->my-1);
402:   hx = 1.0/dhx;                   hy = 1.0/dhy;
403:   hxdhy = hx*dhy;                 hydhx = hy*dhx;

405:   i = st->i; j = st->j; c = st->c;

407:   /* Test whether we are on the right edge of the global array */
408:   if (i == info->mx-1) {
409:     if (c == 0)      *f = x[j][i].u;
410:     else if (c == 1) *f = x[j][i].v;
411:     else if (c == 2) *f = x[j][i].omega - (x[j][i].v - x[j][i-1].v)*dhx;
412:     else             *f = x[j][i].temp - (PetscReal)(grashof>0);

414:   /* Test whether we are on the left edge of the global array */
415:   } else if (i == 0) {
416:     if (c == 0)      *f = x[j][i].u;
417:     else if (c == 1) *f = x[j][i].v;
418:     else if (c == 2) *f = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx;
419:     else             *f = x[j][i].temp;

421:   /* Test whether we are on the top edge of the global array */
422:   } else if (j == info->my-1) {
423:     if (c == 0)      *f = x[j][i].u - lid;
424:     else if (c == 1) *f = x[j][i].v;
425:     else if (c == 2) *f = x[j][i].omega + (x[j][i].u - x[j-1][i].u)*dhy;
426:     else             *f = x[j][i].temp-x[j-1][i].temp;

428:   /* Test whether we are on the bottom edge of the global array */
429:   } else if (j == 0) {
430:     if (c == 0)      *f = x[j][i].u;
431:     else if (c == 1) *f = x[j][i].v;
432:     else if (c == 2) *f = x[j][i].omega + (x[j+1][i].u - x[j][i].u)*dhy;
433:     else             *f = x[j][i].temp - x[j+1][i].temp;

435:   /* Compute over the interior points */
436:   } else {
437:     /*
438:       convective coefficients for upwinding
439:     */
440:     vx = x[j][i].u; avx = PetscAbsScalar(vx);
441:     vxp = .5*(vx+avx); vxm = .5*(vx-avx);
442:     vy = x[j][i].v; avy = PetscAbsScalar(vy);
443:     vyp = .5*(vy+avy); vym = .5*(vy-avy);

445:     /* U velocity */
446:     if (c == 0) {
447:       u          = x[j][i].u;
448:       uxx        = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
449:       uyy        = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
450:       *f         = uxx + uyy - .5*(x[j+1][i].omega-x[j-1][i].omega)*hx;

452:     /* V velocity */
453:     } else if (c == 1) {
454:       u          = x[j][i].v;
455:       uxx        = (2.0*u - x[j][i-1].v - x[j][i+1].v)*hydhx;
456:       uyy        = (2.0*u - x[j-1][i].v - x[j+1][i].v)*hxdhy;
457:       *f         = uxx + uyy + .5*(x[j][i+1].omega-x[j][i-1].omega)*hy;
458: 
459:     /* Omega */
460:     } else if (c == 2) {
461:       u          = x[j][i].omega;
462:       uxx        = (2.0*u - x[j][i-1].omega - x[j][i+1].omega)*hydhx;
463:       uyy        = (2.0*u - x[j-1][i].omega - x[j+1][i].omega)*hxdhy;
464:       *f         = uxx + uyy +
465:         (vxp*(u - x[j][i-1].omega) +
466:          vxm*(x[j][i+1].omega - u)) * hy +
467:         (vyp*(u - x[j-1][i].omega) +
468:          vym*(x[j+1][i].omega - u)) * hx -
469:          .5 * grashof * (x[j][i+1].temp - x[j][i-1].temp) * hy;
470: 
471:     /* Temperature */
472:     } else {
473:       u           = x[j][i].temp;
474:       uxx         = (2.0*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx;
475:       uyy         = (2.0*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy;
476:       *f          =  uxx + uyy  + prandtl * (
477:                                              (vxp*(u - x[j][i-1].temp) +
478:                                               vxm*(x[j][i+1].temp - u)) * hy +
479:                                              (vyp*(u - x[j-1][i].temp) +
480:                                              vym*(x[j+1][i].temp - u)) * hx);
481:     }
482:   }

484:   return(0);
485: }

487: /*
488:     This function that evaluates the function for a single 
489:     grid point. It is used by the -dmmg_fas -dmmg_fas_block solver
490: */
491: PetscErrorCode FormFunctionLocali4(DALocalInfo *info,MatStencil *st,Field **x,PetscScalar *ff,void *ptr)
492: {
493:   Field       *f = (Field*)ff;
494:   AppCtx      *user = (AppCtx*)ptr;
495:   PetscInt    i,j;
496:   PassiveReal hx,hy,dhx,dhy,hxdhy,hydhx;
497:   PassiveReal grashof,prandtl,lid;
498:   PetscScalar u,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;

501:   grashof = user->grashof;
502:   prandtl = user->prandtl;
503:   lid     = user->lidvelocity;

505:   /* 
506:      Define mesh intervals ratios for uniform grid.
507:      [Note: FD formulae below are normalized by multiplying through by
508:      local volume element to obtain coefficients O(1) in two dimensions.]
509:   */
510:   dhx = (PetscReal)(info->mx-1);     dhy = (PetscReal)(info->my-1);
511:   hx = 1.0/dhx;                   hy = 1.0/dhy;
512:   hxdhy = hx*dhy;                 hydhx = hy*dhx;

514:   i = st->i; j = st->j;

516:   /* Test whether we are on the right edge of the global array */
517:   if (i == info->mx-1) {
518:     f->u = x[j][i].u;
519:     f->v = x[j][i].v;
520:     f->omega = x[j][i].omega - (x[j][i].v - x[j][i-1].v)*dhx;
521:     f->temp = x[j][i].temp - (PetscReal)(grashof>0);

523:     /* Test whether we are on the left edge of the global array */
524:   } else if (i == 0) {
525:     f->u = x[j][i].u;
526:     f->v = x[j][i].v;
527:     f->omega = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx;
528:     f->temp = x[j][i].temp;
529: 
530:     /* Test whether we are on the top edge of the global array */
531:   } else if (j == info->my-1) {
532:     f->u = x[j][i].u - lid;
533:     f->v = x[j][i].v;
534:     f->omega = x[j][i].omega + (x[j][i].u - x[j-1][i].u)*dhy;
535:     f->temp = x[j][i].temp-x[j-1][i].temp;

537:     /* Test whether we are on the bottom edge of the global array */
538:   } else if (j == 0) {
539:     f->u = x[j][i].u;
540:     f->v = x[j][i].v;
541:     f->omega = x[j][i].omega + (x[j+1][i].u - x[j][i].u)*dhy;
542:     f->temp = x[j][i].temp - x[j+1][i].temp;

544:     /* Compute over the interior points */
545:   } else {
546:     /*
547:       convective coefficients for upwinding
548:     */
549:     vx = x[j][i].u; avx = PetscAbsScalar(vx);
550:     vxp = .5*(vx+avx); vxm = .5*(vx-avx);
551:     vy = x[j][i].v; avy = PetscAbsScalar(vy);
552:     vyp = .5*(vy+avy); vym = .5*(vy-avy);

554:     /* U velocity */
555:     u          = x[j][i].u;
556:     uxx        = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
557:     uyy        = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
558:     f->u        = uxx + uyy - .5*(x[j+1][i].omega-x[j-1][i].omega)*hx;

560:     /* V velocity */
561:     u          = x[j][i].v;
562:     uxx        = (2.0*u - x[j][i-1].v - x[j][i+1].v)*hydhx;
563:     uyy        = (2.0*u - x[j-1][i].v - x[j+1][i].v)*hxdhy;
564:     f->v        = uxx + uyy + .5*(x[j][i+1].omega-x[j][i-1].omega)*hy;
565: 
566:     /* Omega */
567:     u          = x[j][i].omega;
568:     uxx        = (2.0*u - x[j][i-1].omega - x[j][i+1].omega)*hydhx;
569:     uyy        = (2.0*u - x[j-1][i].omega - x[j+1][i].omega)*hxdhy;
570:     f->omega    = uxx + uyy +
571:       (vxp*(u - x[j][i-1].omega) +
572:        vxm*(x[j][i+1].omega - u)) * hy +
573:       (vyp*(u - x[j-1][i].omega) +
574:        vym*(x[j+1][i].omega - u)) * hx -
575:        .5 * grashof * (x[j][i+1].temp - x[j][i-1].temp) * hy;
576: 
577:     /* Temperature */
578: 
579:     u           = x[j][i].temp;
580:     uxx         = (2.0*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx;
581:     uyy         = (2.0*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy;
582:     f->temp     =  uxx + uyy  + prandtl * (
583:                                            (vxp*(u - x[j][i-1].temp) +
584:                                             vxm*(x[j][i+1].temp - u)) * hy +
585:                                            (vyp*(u - x[j-1][i].temp) +
586:                                             vym*(x[j+1][i].temp - u)) * hx);
587:   }
588:   return(0);
589: }