Actual source code: ex19.c

  1: /*$Id: ex19.c,v 1.30 2001/08/07 21:31:17 bsmith Exp $*/

  3: static char help[] = "Nonlinear driven cavity with multigrid in 2d.n
  4:   n
  5: The 2D driven cavity problem is solved in a velocity-vorticity formulation.n
  6: The flow can be driven with the lid or with bouyancy or both:n
  7:   -lidvelocity <lid>, where <lid> = dimensionless velocity of lidn
  8:   -grashof <gr>, where <gr> = dimensionless temperature gradentn
  9:   -prandtl <pr>, where <pr> = dimensionless thermal/momentum diffusity ration
 10:   -contours : draw contour plots of solutionnn";

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

 19: /* ------------------------------------------------------------------------

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

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

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

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

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

 53:   ------------------------------------------------------------------------- */

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

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

 75: extern int FormInitialGuess(SNES,Vec,void*);
 76: extern int FormFunction(SNES,Vec,Vec,void*);
 77: extern int FormFunctionLocal(DALocalInfo*,Field**,Field**,void*);
 78: extern int FormFunctionLocali(DALocalInfo*,MatStencil*,Field**,PetscScalar*,void*);

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

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

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


100:   PreLoadBegin(PETSC_TRUE,"SetUp");
101:     DMMGCreate(comm,2,&user,&dmmg);


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

112:     DAGetInfo(DMMGGetDA(dmmg),0,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
113:                      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:     DASetFieldName(DMMGGetDA(dmmg),0,"x-velocity");
126:     DASetFieldName(DMMGGetDA(dmmg),1,"y-velocity");
127:     DASetFieldName(DMMGGetDA(dmmg),2,"Omega");
128:     DASetFieldName(DMMGGetDA(dmmg),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:        Process adiC: FormFunctionLocal FormFunctionLocali
139:        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140:     PetscOptionsGetLogical(PETSC_NULL,"-localfunction",&localfunction,PETSC_IGNORE);
141:     if (localfunction) {
142:       DMMGSetSNESLocal(dmmg,FormFunctionLocal,0,ad_FormFunctionLocal,admf_FormFunctionLocal);
143:       DMMGSetSNESLocali(dmmg,FormFunctionLocali,ad_FormFunctionLocali,admf_FormFunctionLocali);
144:     } else {
145:       DMMGSetSNES(dmmg,FormFunction,0);
146:     }

148:     PetscPrintf(comm,"lid velocity = %g, prandtl # = %g, grashof # = %gn",
149:                        user.lidvelocity,user.prandtl,user.grashof);


152:     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153:        Solve the nonlinear system
154:        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
155:     DMMGSetInitialGuess(dmmg,FormInitialGuess);

157:   PreLoadStage("Solve");
158:     DMMGSolve(dmmg);

160:     snes = DMMGGetSNES(dmmg);
161:     SNESGetIterationNumber(snes,&its);
162:     PetscPrintf(comm,"Number of Newton iterations = %dn", its);

164:     /*
165:       Visualize solution
166:     */

168:     if (user.draw_contours) {
169:       VecView(DMMGGetx(dmmg),PETSC_VIEWER_DRAW_WORLD);
170:     }

172:     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
173:        Free work space.  All PETSc objects should be destroyed when they
174:        are no longer needed.
175:        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

177:     DMMGDestroy(dmmg);
178:   PreLoadEnd();

180:   PetscFinalize();
181:   return 0;
182: }

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


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

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

194:    Output Parameter:
195:    X - vector
196:  */
197: int FormInitialGuess(SNES snes,Vec X,void *ptr)
198: {
199:   DMMG      dmmg = (DMMG)ptr;
200:   AppCtx    *user = (AppCtx*)dmmg->user;
201:   DA        da = (DA)dmmg->dm;
202:   int       i,j,mx,ierr,xs,ys,xm,ym;
203:   PetscReal grashof,dx;
204:   Field     **x;

206:   grashof = user->grashof;

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

211:   /*
212:      Get local grid boundaries (for 2-dimensional DA):
213:        xs, ys   - starting grid indices (no ghost points)
214:        xm, ym   - widths of local grid (no ghost points)
215:   */
216:   DAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);

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

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

240:   /*
241:      Restore vector
242:   */
243:   DAVecRestoreArray(da,X,(void**)&x);
244:   return 0;
245: }
246: /* ------------------------------------------------------------------- */
247: /* 
248:    FormFunction - Evaluates the nonlinear function, F(x).

250:    Input Parameters:
251: .  snes - the SNES context
252: .  X - input vector
253: .  ptr - optional user-defined context, as set by SNESSetFunction()

255:    Output Parameter:
256: .  F - function vector

258:    Notes:
259:    We process the boundary nodes before handling the interior
260:    nodes, so that no conditional statements are needed within the
261:    PetscReal loop over the local grid indices. 
262:  */
263: int FormFunction(SNES snes,Vec X,Vec F,void *ptr)
264: {
265:   DMMG         dmmg = (DMMG)ptr;
266:   AppCtx       *user = (AppCtx*)dmmg->user;
267:   int          ierr,i,j,mx,my,xs,ys,xm,ym;
268:   int          xints,xinte,yints,yinte;
269:   PetscReal    two = 2.0,one = 1.0,p5 = 0.5,hx,hy,dhx,dhy,hxdhy,hydhx;
270:   PetscReal    grashof,prandtl,lid;
271:   PetscScalar  u,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;
272:   Field        **x,**f;
273:   Vec          localX;
274:   DA           da = (DA)dmmg->dm;

276:   DAGetLocalVector((DA)dmmg->dm,&localX);
277:   DAGetInfo(da,0,&mx,&my,0,0,0,0,0,0,0,0);

279:   grashof = user->grashof;
280:   prandtl = user->prandtl;
281:   lid     = user->lidvelocity;

283:   /* 
284:      Define mesh intervals ratios for uniform grid.
285:      [Note: FD formulae below are normalized by multiplying through by
286:      local volume element to obtain coefficients O(1) in two dimensions.]
287:   */
288:   dhx = (PetscReal)(mx-1);     dhy = (PetscReal)(my-1);
289:   hx = one/dhx;             hy = one/dhy;
290:   hxdhy = hx*dhy;           hydhx = hy*dhx;

292:   /*
293:      Scatter ghost points to local vector, using the 2-step process
294:         DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
295:      By placing code between these two statements, computations can be
296:      done while messages are in transition.
297:   */
298:   DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
299:   DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);

301:   /*
302:      Get pointers to vector data
303:   */
304:   DAVecGetArray((DA)dmmg->dm,localX,(void**)&x);
305:   DAVecGetArray((DA)dmmg->dm,F,(void**)&f);

307:   /*
308:      Get local grid boundaries
309:   */
310:   DAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);

312:   /*
313:      Compute function over the locally owned part of the grid
314:      (physical corner points are set twice to avoid more conditionals).
315:   */
316:   xints = xs; xinte = xs+xm; yints = ys; yinte = ys+ym;

318:   /* Test whether we are on the bottom edge of the global array */
319:   if (yints == 0) {
320:     j = 0;
321:     yints = yints + 1;
322:     /* bottom edge */
323:     for (i=xs; i<xs+xm; i++) {
324:         f[j][i].u     = x[j][i].u;
325:         f[j][i].v     = x[j][i].v;
326:         f[j][i].omega = x[j][i].omega + (x[j+1][i].u - x[j][i].u)*dhy;
327:         f[j][i].temp  = x[j][i].temp-x[j+1][i].temp;
328:     }
329:   }

331:   /* Test whether we are on the top edge of the global array */
332:   if (yinte == my) {
333:     j = my - 1;
334:     yinte = yinte - 1;
335:     /* top edge */
336:     for (i=xs; i<xs+xm; i++) {
337:         f[j][i].u     = x[j][i].u - lid;
338:         f[j][i].v     = x[j][i].v;
339:         f[j][i].omega = x[j][i].omega + (x[j][i].u - x[j-1][i].u)*dhy;
340:         f[j][i].temp  = x[j][i].temp-x[j-1][i].temp;
341:     }
342:   }

344:   /* Test whether we are on the left edge of the global array */
345:   if (xints == 0) {
346:     i = 0;
347:     xints = xints + 1;
348:     /* left edge */
349:     for (j=ys; j<ys+ym; j++) {
350:       f[j][i].u     = x[j][i].u;
351:       f[j][i].v     = x[j][i].v;
352:       f[j][i].omega = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx;
353:       f[j][i].temp  = x[j][i].temp;
354:     }
355:   }

357:   /* Test whether we are on the right edge of the global array */
358:   if (xinte == mx) {
359:     i = mx - 1;
360:     xinte = xinte - 1;
361:     /* right edge */
362:     for (j=ys; j<ys+ym; j++) {
363:       f[j][i].u     = x[j][i].u;
364:       f[j][i].v     = x[j][i].v;
365:       f[j][i].omega = x[j][i].omega - (x[j][i].v - x[j][i-1].v)*dhx;
366:       f[j][i].temp  = x[j][i].temp - (PetscReal)(grashof>0);
367:     }
368:   }

370:   /* Compute over the interior points */
371:   for (j=yints; j<yinte; j++) {
372:     for (i=xints; i<xinte; i++) {

374:         /*
375:           convective coefficients for upwinding
376:         */
377:         vx = x[j][i].u; avx = PetscAbsScalar(vx);
378:         vxp = p5*(vx+avx); vxm = p5*(vx-avx);
379:         vy = x[j][i].v; avy = PetscAbsScalar(vy);
380:         vyp = p5*(vy+avy); vym = p5*(vy-avy);

382:         /* U velocity */
383:         u          = x[j][i].u;
384:         uxx        = (two*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
385:         uyy        = (two*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
386:         f[j][i].u  = uxx + uyy - p5*(x[j+1][i].omega-x[j-1][i].omega)*hx;

388:         /* V velocity */
389:         u          = x[j][i].v;
390:         uxx        = (two*u - x[j][i-1].v - x[j][i+1].v)*hydhx;
391:         uyy        = (two*u - x[j-1][i].v - x[j+1][i].v)*hxdhy;
392:         f[j][i].v  = uxx + uyy + p5*(x[j][i+1].omega-x[j][i-1].omega)*hy;

394:         /* Omega */
395:         u          = x[j][i].omega;
396:         uxx        = (two*u - x[j][i-1].omega - x[j][i+1].omega)*hydhx;
397:         uyy        = (two*u - x[j-1][i].omega - x[j+1][i].omega)*hxdhy;
398:         f[j][i].omega = uxx + uyy +
399:                         (vxp*(u - x[j][i-1].omega) +
400:                           vxm*(x[j][i+1].omega - u)) * hy +
401:                         (vyp*(u - x[j-1][i].omega) +
402:                           vym*(x[j+1][i].omega - u)) * hx -
403:                         p5 * grashof * (x[j][i+1].temp - x[j][i-1].temp) * hy;

405:         /* Temperature */
406:         u             = x[j][i].temp;
407:         uxx           = (two*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx;
408:         uyy           = (two*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy;
409:         f[j][i].temp =  uxx + uyy  + prandtl * (
410:                         (vxp*(u - x[j][i-1].temp) +
411:                           vxm*(x[j][i+1].temp - u)) * hy +
412:                         (vyp*(u - x[j-1][i].temp) +
413:                                  vym*(x[j+1][i].temp - u)) * hx);
414:     }
415:   }

417:   /*
418:      Restore vectors
419:   */
420:   DAVecRestoreArray((DA)dmmg->dm,localX,(void**)&x);
421:   DAVecRestoreArray((DA)dmmg->dm,F,(void**)&f);

423:   DARestoreLocalVector((DA)dmmg->dm,&localX);

425:   /*
426:      Flop count (multiply-adds are counted as 2 operations)
427:   */
428:   PetscLogFlops(84*ym*xm);

430:   return 0;
431: }

433: int FormFunctionLocal(DALocalInfo *info,Field **x,Field **f,void *ptr)
434:  {
435:   AppCtx       *user = (AppCtx*)ptr;
436:   int          ierr,i,j;
437:   int          xints,xinte,yints,yinte;
438:   PetscReal    hx,hy,dhx,dhy,hxdhy,hydhx;
439:   PetscReal    grashof,prandtl,lid;
440:   PetscScalar  u,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;

443:   grashof = user->grashof;
444:   prandtl = user->prandtl;
445:   lid     = user->lidvelocity;

447:   /* 
448:      Define mesh intervals ratios for uniform grid.
449:      [Note: FD formulae below are normalized by multiplying through by
450:      local volume element to obtain coefficients O(1) in two dimensions.]
451:   */
452:   dhx = (PetscReal)(info->mx-1);     dhy = (PetscReal)(info->my-1);
453:   hx = 1.0/dhx;                   hy = 1.0/dhy;
454:   hxdhy = hx*dhy;                 hydhx = hy*dhx;

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

458:   /* Test whether we are on the bottom edge of the global array */
459:   if (yints == 0) {
460:     j = 0;
461:     yints = yints + 1;
462:     /* bottom edge */
463:     for (i=info->xs; i<info->xs+info->xm; i++) {
464:         f[j][i].u     = x[j][i].u;
465:         f[j][i].v     = x[j][i].v;
466:         f[j][i].omega = x[j][i].omega + (x[j+1][i].u - x[j][i].u)*dhy;
467:         f[j][i].temp  = x[j][i].temp-x[j+1][i].temp;
468:     }
469:   }

471:   /* Test whether we are on the top edge of the global array */
472:   if (yinte == info->my) {
473:     j = info->my - 1;
474:     yinte = yinte - 1;
475:     /* top edge */
476:     for (i=info->xs; i<info->xs+info->xm; i++) {
477:         f[j][i].u     = x[j][i].u - lid;
478:         f[j][i].v     = x[j][i].v;
479:         f[j][i].omega = x[j][i].omega + (x[j][i].u - x[j-1][i].u)*dhy;
480:         f[j][i].temp  = x[j][i].temp-x[j-1][i].temp;
481:     }
482:   }

484:   /* Test whether we are on the left edge of the global array */
485:   if (xints == 0) {
486:     i = 0;
487:     xints = xints + 1;
488:     /* left edge */
489:     for (j=info->ys; j<info->ys+info->ym; j++) {
490:       f[j][i].u     = x[j][i].u;
491:       f[j][i].v     = x[j][i].v;
492:       f[j][i].omega = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx;
493:       f[j][i].temp  = x[j][i].temp;
494:     }
495:   }

497:   /* Test whether we are on the right edge of the global array */
498:   if (xinte == info->mx) {
499:     i = info->mx - 1;
500:     xinte = xinte - 1;
501:     /* right edge */
502:     for (j=info->ys; j<info->ys+info->ym; j++) {
503:       f[j][i].u     = x[j][i].u;
504:       f[j][i].v     = x[j][i].v;
505:       f[j][i].omega = x[j][i].omega - (x[j][i].v - x[j][i-1].v)*dhx;
506:       f[j][i].temp  = x[j][i].temp - (PetscReal)(grashof>0);
507:     }
508:   }

510:   /* Compute over the interior points */
511:   for (j=yints; j<yinte; j++) {
512:     for (i=xints; i<xinte; i++) {

514:         /*
515:           convective coefficients for upwinding
516:         */
517:         vx = x[j][i].u; avx = PetscAbsScalar(vx);
518:         vxp = .5*(vx+avx); vxm = .5*(vx-avx);
519:         vy = x[j][i].v; avy = PetscAbsScalar(vy);
520:         vyp = .5*(vy+avy); vym = .5*(vy-avy);

522:         /* U velocity */
523:         u          = x[j][i].u;
524:         uxx        = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
525:         uyy        = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
526:         f[j][i].u  = uxx + uyy - .5*(x[j+1][i].omega-x[j-1][i].omega)*hx;

528:         /* V velocity */
529:         u          = x[j][i].v;
530:         uxx        = (2.0*u - x[j][i-1].v - x[j][i+1].v)*hydhx;
531:         uyy        = (2.0*u - x[j-1][i].v - x[j+1][i].v)*hxdhy;
532:         f[j][i].v  = uxx + uyy + .5*(x[j][i+1].omega-x[j][i-1].omega)*hy;

534:         /* Omega */
535:         u          = x[j][i].omega;
536:         uxx        = (2.0*u - x[j][i-1].omega - x[j][i+1].omega)*hydhx;
537:         uyy        = (2.0*u - x[j-1][i].omega - x[j+1][i].omega)*hxdhy;
538:         f[j][i].omega = uxx + uyy +
539:                         (vxp*(u - x[j][i-1].omega) +
540:                           vxm*(x[j][i+1].omega - u)) * hy +
541:                         (vyp*(u - x[j-1][i].omega) +
542:                           vym*(x[j+1][i].omega - u)) * hx -
543:                         .5 * grashof * (x[j][i+1].temp - x[j][i-1].temp) * hy;

545:         /* Temperature */
546:         u             = x[j][i].temp;
547:         uxx           = (2.0*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx;
548:         uyy           = (2.0*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy;
549:         f[j][i].temp =  uxx + uyy  + prandtl * (
550:                         (vxp*(u - x[j][i-1].temp) +
551:                           vxm*(x[j][i+1].temp - u)) * hy +
552:                         (vyp*(u - x[j-1][i].temp) +
553:                                  vym*(x[j+1][i].temp - u)) * hx);
554:     }
555:   }

557:   /*
558:      Flop count (multiply-adds are counted as 2 operations)
559:   */
560:   PetscLogFlops(84*info->ym*info->xm);
561:   return(0);
562: }

564: int FormFunctionLocali(DALocalInfo *info,MatStencil *st,Field **x,PetscScalar *f,void *ptr)
565:  {
566:   AppCtx      *user = (AppCtx*)ptr;
567:   int         i,j,c;
568:   PassiveReal hx,hy,dhx,dhy,hxdhy,hydhx;
569:   PassiveReal grashof,prandtl,lid;
570:   PetscScalar u,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;

573:   grashof = user->grashof;
574:   prandtl = user->prandtl;
575:   lid     = user->lidvelocity;

577:   /* 
578:      Define mesh intervals ratios for uniform grid.
579:      [Note: FD formulae below are normalized by multiplying through by
580:      local volume element to obtain coefficients O(1) in two dimensions.]
581:   */
582:   dhx = (PetscReal)(info->mx-1);     dhy = (PetscReal)(info->my-1);
583:   hx = 1.0/dhx;                   hy = 1.0/dhy;
584:   hxdhy = hx*dhy;                 hydhx = hy*dhx;

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

588:   /* Test whether we are on the right edge of the global array */
589:   if (i == info->mx-1) {
590:     if (c == 0) *f     = x[j][i].u;
591:     else if (c == 1) *f     = x[j][i].v;
592:     else if (c == 2) *f = x[j][i].omega - (x[j][i].v - x[j][i-1].v)*dhx;
593:     else *f  = x[j][i].temp - (PetscReal)(grashof>0);

595:   /* Test whether we are on the left edge of the global array */
596:   } else if (i == 0) {
597:     if (c == 0) *f     = x[j][i].u;
598:     else if (c == 1) *f     = x[j][i].v;
599:     else if (c == 2) *f = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx;
600:     else *f  = x[j][i].temp;

602:   /* Test whether we are on the top edge of the global array */
603:   } else if (j == info->my-1) {
604:     if (c == 0) *f     = x[j][i].u - lid;
605:     else if (c == 1) *f     = x[j][i].v;
606:     else if (c == 2) *f = x[j][i].omega + (x[j][i].u - x[j-1][i].u)*dhy;
607:     else *f  = x[j][i].temp-x[j-1][i].temp;

609:   /* Test whether we are on the bottom edge of the global array */
610:   } else if (j == 0) {
611:     if (c == 0) *f     = x[j][i].u;
612:     else if (c == 1) *f     = x[j][i].v;
613:     else if (c == 2) *f = x[j][i].omega + (x[j+1][i].u - x[j][i].u)*dhy;
614:     else *f  = x[j][i].temp-x[j+1][i].temp;

616:   /* Compute over the interior points */
617:   } else {
618:     /*
619:       convective coefficients for upwinding
620:     */
621:     vx = x[j][i].u; avx = PetscAbsScalar(vx);
622:     vxp = .5*(vx+avx); vxm = .5*(vx-avx);
623:     vy = x[j][i].v; avy = PetscAbsScalar(vy);
624:     vyp = .5*(vy+avy); vym = .5*(vy-avy);

626:     /* U velocity */
627:     if (c == 0) {
628:       u          = x[j][i].u;
629:       uxx        = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
630:       uyy        = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
631:       *f         = uxx + uyy - .5*(x[j+1][i].omega-x[j-1][i].omega)*hx;

633:     /* V velocity */
634:     } else if (c == 1) {
635:       u          = x[j][i].v;
636:       uxx        = (2.0*u - x[j][i-1].v - x[j][i+1].v)*hydhx;
637:       uyy        = (2.0*u - x[j-1][i].v - x[j+1][i].v)*hxdhy;
638:       *f         = uxx + uyy + .5*(x[j][i+1].omega-x[j][i-1].omega)*hy;
639: 
640:     /* Omega */
641:     } else if (c == 2) {
642:       u          = x[j][i].omega;
643:       uxx        = (2.0*u - x[j][i-1].omega - x[j][i+1].omega)*hydhx;
644:       uyy        = (2.0*u - x[j-1][i].omega - x[j+1][i].omega)*hxdhy;
645:       *f         = uxx + uyy +
646:         (vxp*(u - x[j][i-1].omega) +
647:          vxm*(x[j][i+1].omega - u)) * hy +
648:         (vyp*(u - x[j-1][i].omega) +
649:          vym*(x[j+1][i].omega - u)) * hx -
650:         .5 * grashof * (x[j][i+1].temp - x[j][i-1].temp) * hy;
651: 
652:     /* Temperature */
653:     } else {
654:       u           = x[j][i].temp;
655:       uxx         = (2.0*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx;
656:       uyy         = (2.0*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy;
657:       *f          =  uxx + uyy  + prandtl * (
658:                                              (vxp*(u - x[j][i-1].temp) +
659:                                               vxm*(x[j][i+1].temp - u)) * hy +
660:                                              (vyp*(u - x[j-1][i].temp) +
661:                                               vym*(x[j+1][i].temp - u)) * hx);
662:     }
663:   }

665:   return(0);
666: }