Actual source code: ex19.c

petsc-3.4.5 2014-06-29
  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*/

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

 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: #if defined(PETSC_APPLE_FRAMEWORK)
 67: #import <PETSc/petscsnes.h>
 68: #import <PETSc/petscdmda.h>
 69: #else
 70: #include <petscsnes.h>
 71: #include <petscdmda.h>
 72: #endif

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

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

 83: typedef struct {
 84:   PassiveReal lidvelocity,prandtl,grashof;  /* physical parameters */
 85:   PetscBool   draw_contours;                /* flag - 1 indicates drawing contours */
 86: } AppCtx;

 88: extern PetscErrorCode FormInitialGuess(AppCtx*,DM,Vec);
 89: extern PetscErrorCode NonlinearGS(SNES,Vec,Vec,void*);

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

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

106:   comm = PETSC_COMM_WORLD;
107:   SNESCreate(comm,&snes);

109:   /*
110:       Create distributed array object to manage parallel grid and vectors
111:       for principal unknowns (x) and governing residuals (f)
112:   */
113:   DMDACreate2d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_STENCIL_STAR,-4,-4,PETSC_DECIDE,PETSC_DECIDE,4,1,0,0,&da);
114:   SNESSetDM(snes,(DM)da);
115:   SNESSetGS(snes, NonlinearGS, (void*)&user);

117:   DMDAGetInfo(da,0,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
118:                      PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
119:   /*
120:      Problem parameters (velocity of lid, prandtl, and grashof numbers)
121:   */
122:   user.lidvelocity = 1.0/(mx*my);
123:   user.prandtl     = 1.0;
124:   user.grashof     = 1.0;

126:   PetscOptionsGetReal(NULL,"-lidvelocity",&user.lidvelocity,NULL);
127:   PetscOptionsGetReal(NULL,"-prandtl",&user.prandtl,NULL);
128:   PetscOptionsGetReal(NULL,"-grashof",&user.grashof,NULL);
129:   PetscOptionsHasName(NULL,"-contours",&user.draw_contours);

131:   DMDASetFieldName(da,0,"x_velocity");
132:   DMDASetFieldName(da,1,"y_velocity");
133:   DMDASetFieldName(da,2,"Omega");
134:   DMDASetFieldName(da,3,"temperature");

136:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137:      Create user context, set problem data, create vector data structures.
138:      Also, compute the initial guess.
139:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

141:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142:      Create nonlinear solver context

144:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145:   DMSetApplicationContext(da,&user);
146:   DMDASNESSetFunctionLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))FormFunctionLocal,&user);
147:   SNESSetFromOptions(snes);
148:   PetscPrintf(comm,"lid velocity = %G, prandtl # = %G, grashof # = %G\n",user.lidvelocity,user.prandtl,user.grashof);


151:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152:      Solve the nonlinear system
153:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
154:   DMCreateGlobalVector(da,&x);
155:   FormInitialGuess(&user,da,x);

157:   SNESSolve(snes,NULL,x);

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

162:   /*
163:      Visualize solution
164:   */
165:   if (user.draw_contours) {
166:     VecView(x,PETSC_VIEWER_DRAW_WORLD);
167:   }

169:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
170:      Free work space.  All PETSc objects should be destroyed when they
171:      are no longer needed.
172:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
173:   VecDestroy(&x);
174:   DMDestroy(&da);
175:   SNESDestroy(&snes);
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(AppCtx *user,DM da,Vec X)
196: {
197:   PetscInt       i,j,mx,xs,ys,xm,ym;
199:   PetscReal      grashof,dx;
200:   Field          **x;

202:   grashof = user->grashof;

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

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

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

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

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

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

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

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

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


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

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

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

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

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

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

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

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

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

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

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

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

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


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

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

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

435:   /* Set the boundary conditions on the momentum equations */
436:   /* Test whether we are on the bottom edge of the global array */
437:   if (yints == 0) {
438:     j     = 0;
439:     yints = yints + 1;
440:     /* bottom edge */
441:     for (i=info.xs; i<info.xs+info.xm; i++) {

443:       if (B) {
444:         bjiu = b[j][i].u;
445:         bjiv = b[j][i].v;
446:       } else {
447:         bjiu = 0.0;
448:         bjiv = 0.0;
449:       }
450:       x[j][i].u = 0.0 + bjiu;
451:       x[j][i].v = 0.0 + bjiv;
452:     }
453:   }

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

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

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

509:   for (k=0; k < sweeps; k++) {
510:     for (j=info.ys; j<info.ys + info.ym; j++) {
511:       for (i=info.xs; i<info.xs + info.xm; i++) {
512:         ptconverged = PETSC_FALSE;
513:         pfnorm0     = 0.0;
514:         pfnorm      = 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 +
572:                                            (vyp*(u - x[j-1][i].temp) + vym*(x[j+1][i].temp - u))*hx) - bjitemp;
573:             dftdt = 2.0*(hydhx + hxdhy) + prandtl*((vxp - vxm)*hy + (vyp - vym)*hx);
574:             if (PetscRealPart(vx) > 0.0) dftdu = prandtl*(u - x[j][i-1].temp)*hy;
575:             else dftdu = prandtl*(x[j][i+1].temp - u)*hy;

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

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

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