Actual source code: ex19.c

petsc-master 2016-12-02
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).

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

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

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

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

 91: extern PetscErrorCode FormInitialGuess(AppCtx*,DM,Vec);
 92: 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: /* ------------------------------------------------------------------- */


189: /*
190:    FormInitialGuess - Forms initial approximation.

192:    Input Parameters:
193:    user - user-defined application context
194:    X - vector

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

207:   grashof = user->grashof;

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

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

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

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

241:   /*
242:      Restore vector
243:   */
244:   DMDAVecRestoreArray(da,X,&x);
245:   return(0);
246: }

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

260:   grashof = user->grashof;
261:   prandtl = user->prandtl;
262:   lid     = user->lidvelocity;

264:   /*
265:      Define mesh intervals ratios for uniform grid.

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


271:   */
272:   dhx   = (PetscReal)(info->mx-1);  dhy = (PetscReal)(info->my-1);
273:   hx    = 1.0/dhx;                   hy = 1.0/dhy;
274:   hxdhy = hx*dhy;                 hydhx = hy*dhx;

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

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

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

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

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

330:   /* Compute over the interior points */
331:   for (j=yints; j<yinte; j++) {
332:     for (i=xints; i<xinte; i++) {

334:       /*
335:        convective coefficients for upwinding
336:       */
337:       vx  = x[j][i].u; avx = PetscAbsScalar(vx);
338:       vxp = .5*(vx+avx); vxm = .5*(vx-avx);
339:       vy  = x[j][i].v; avy = PetscAbsScalar(vy);
340:       vyp = .5*(vy+avy); vym = .5*(vy-avy);

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

348:       /* V velocity */
349:       u         = x[j][i].v;
350:       uxx       = (2.0*u - x[j][i-1].v - x[j][i+1].v)*hydhx;
351:       uyy       = (2.0*u - x[j-1][i].v - x[j+1][i].v)*hxdhy;
352:       f[j][i].v = uxx + uyy + .5*(x[j][i+1].omega-x[j][i-1].omega)*hy;

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

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

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


381: PetscErrorCode NonlinearGS(SNES snes, Vec X, Vec B, void *ctx)
382: {
383:   DMDALocalInfo  info;
384:   Field          **x,**b;
386:   Vec            localX, localB;
387:   DM             da;
388:   PetscInt       xints,xinte,yints,yinte,i,j,k,l;
389:   PetscInt       max_its,tot_its;
390:   PetscInt       sweeps;
391:   PetscReal      rtol,atol,stol;
392:   PetscReal      hx,hy,dhx,dhy,hxdhy,hydhx;
393:   PetscReal      grashof,prandtl,lid;
394:   PetscScalar    u,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;
395:   PetscScalar    fu, fv, fomega, ftemp;
396:   PetscScalar    dfudu;
397:   PetscScalar    dfvdv;
398:   PetscScalar    dfodu, dfodv, dfodo;
399:   PetscScalar    dftdu, dftdv, dftdt;
400:   PetscScalar    yu=0, yv=0, yo=0, yt=0;
401:   PetscScalar    bjiu, bjiv, bjiomega, bjitemp;
402:   PetscBool      ptconverged;
403:   PetscReal      pfnorm,pfnorm0,pynorm,pxnorm;
404:   AppCtx         *user = (AppCtx*)ctx;

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

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

440:   /* Set the boundary conditions on the momentum equations */
441:   /* Test whether we are on the bottom edge of the global array */
442:   if (yints == 0) {
443:     j     = 0;
444:     /* bottom edge */
445:     for (i=info.xs; i<info.xs+info.xm; i++) {

447:       if (B) {
448:         bjiu = b[j][i].u;
449:         bjiv = b[j][i].v;
450:       } else {
451:         bjiu = 0.0;
452:         bjiv = 0.0;
453:       }
454:       x[j][i].u = 0.0 + bjiu;
455:       x[j][i].v = 0.0 + bjiv;
456:     }
457:   }

459:   /* Test whether we are on the top edge of the global array */
460:   if (yinte == info.my) {
461:     j     = info.my - 1;
462:     /* top edge */
463:     for (i=info.xs; i<info.xs+info.xm; i++) {
464:       if (B) {
465:         bjiu = b[j][i].u;
466:         bjiv = b[j][i].v;
467:       } else {
468:         bjiu = 0.0;
469:         bjiv = 0.0;
470:       }
471:       x[j][i].u = lid + bjiu;
472:       x[j][i].v = bjiv;
473:     }
474:   }

476:   /* Test whether we are on the left edge of the global array */
477:   if (xints == 0) {
478:     i     = 0;
479:     /* left edge */
480:     for (j=info.ys; j<info.ys+info.ym; j++) {
481:       if (B) {
482:         bjiu = b[j][i].u;
483:         bjiv = b[j][i].v;
484:       } else {
485:         bjiu = 0.0;
486:         bjiv = 0.0;
487:       }
488:       x[j][i].u = 0.0 + bjiu;
489:       x[j][i].v = 0.0 + bjiv;
490:     }
491:   }

493:   /* Test whether we are on the right edge of the global array */
494:   if (xinte == info.mx) {
495:     i     = info.mx - 1;
496:     /* right edge */
497:     for (j=info.ys; j<info.ys+info.ym; j++) {
498:       if (B) {
499:         bjiu = b[j][i].u;
500:         bjiv = b[j][i].v;
501:       } else {
502:         bjiu = 0.0;
503:         bjiv = 0.0;
504:       }
505:       x[j][i].u = 0.0 + bjiu;
506:       x[j][i].v = 0.0 + bjiv;
507:     }
508:   }

510:   for (k=0; k < sweeps; k++) {
511:     for (j=info.ys; j<info.ys + info.ym; j++) {
512:       for (i=info.xs; i<info.xs + info.xm; i++) {
513:         ptconverged = PETSC_FALSE;
514:         pfnorm0     = 0.0;
515:         fu          = 0.0;
516:         fv          = 0.0;
517:         fomega      = 0.0;
518:         ftemp       = 0.0;
519:         for (l = 0; l < max_its && !ptconverged; l++) {
520:           if (B) {
521:             bjiu     = b[j][i].u;
522:             bjiv     = b[j][i].v;
523:             bjiomega = b[j][i].omega;
524:             bjitemp  = b[j][i].temp;
525:           } else {
526:             bjiu     = 0.0;
527:             bjiv     = 0.0;
528:             bjiomega = 0.0;
529:             bjitemp  = 0.0;
530:           }

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

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

567:             /* Temperature */
568:             u     = x[j][i].temp;
569:             uxx   = (2.0*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx;
570:             uyy   = (2.0*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy;
571:             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;
572:             dftdt = 2.0*(hydhx + hxdhy) + prandtl*((vxp - vxm)*hy + (vyp - vym)*hx);
573:             if (PetscRealPart(vx) > 0.0) dftdu = prandtl*(u - x[j][i-1].temp)*hy;
574:             else dftdu = prandtl*(x[j][i+1].temp - u)*hy;

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

579:             /* invert the system:
580:              [ dfu / du     0        0        0    ][yu] = [fu]
581:              [     0    dfv / dv     0        0    ][yv]   [fv]
582:              [ dfo / du dfo / dv dfo / do     0    ][yo]   [fo]
583:              [ dft / du dft / dv     0    dft / dt ][yt]   [ft]
584:              by simple back-substitution
585:            */
586:             yu = fu / dfudu;
587:             yv = fv / dfvdv;
588:             yo = (fomega - (dfodu*yu + dfodv*yv)) / dfodo;
589:             yt = (ftemp - (dftdu*yu + dftdv*yv)) / dftdt;

591:             x[j][i].u     = x[j][i].u - yu;
592:             x[j][i].v     = x[j][i].v - yv;
593:             x[j][i].temp  = x[j][i].temp - yt;
594:             x[j][i].omega = x[j][i].omega - yo;
595:           }
596:           if (i == 0) {
597:             fomega        = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx - bjiomega;
598:             ftemp         = x[j][i].temp - 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 (i == info.mx - 1) {
605:             fomega        = x[j][i].omega - (x[j][i].v - x[j][i-1].v)*dhx - bjiomega;
606:             ftemp         = x[j][i].temp - (PetscReal)(grashof>0) - 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 == 0) {
613:             fomega        = x[j][i].omega + (x[j+1][i].u - x[j][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:           if (j == info.my - 1) {
621:             fomega        = x[j][i].omega + (x[j][i].u - x[j-1][i].u)*dhy - bjiomega;
622:             ftemp         = x[j][i].temp-x[j-1][i].temp - bjitemp;
623:             yo            = fomega;
624:             yt            = ftemp;
625:             x[j][i].omega = x[j][i].omega - fomega;
626:             x[j][i].temp  = x[j][i].temp - ftemp;
627:           }
628:           tot_its++;
629:           pfnorm = PetscRealPart(fu*fu + fv*fv + fomega*fomega + ftemp*ftemp);
630:           pfnorm = PetscSqrtReal(pfnorm);
631:           pynorm = PetscRealPart(yu*yu + yv*yv + yo*yo + yt*yt);
632:           pynorm = PetscSqrtReal(pynorm);
633:           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);
634:           pxnorm = PetscSqrtReal(pxnorm);
635:           if (l == 0) pfnorm0 = pfnorm;
636:           if (rtol*pfnorm0 >pfnorm ||
637:               atol > pfnorm ||
638:               pxnorm*stol > pynorm) ptconverged = PETSC_TRUE;
639:         }
640:       }
641:     }
642:   }
643:   DMDAVecRestoreArray(da,localX,&x);
644:   if (B) {
645:     DMDAVecRestoreArrayRead(da,localB,&b);
646:   }
647:   DMLocalToGlobalBegin(da,localX,INSERT_VALUES,X);
648:   DMLocalToGlobalEnd(da,localX,INSERT_VALUES,X);
649:   PetscLogFlops(tot_its*(84.0 + 41.0 + 26.0));
650:   DMRestoreLocalVector(da,&localX);
651:   if (B) {
652:     DMRestoreLocalVector(da,&localB);
653:   }
654:   return(0);
655: }