Actual source code: ex19.c

petsc-master 2016-06-29
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:   SNESSetDM(snes,(DM)da);
118:   SNESSetNGS(snes, NonlinearGS, (void*)&user);

120:   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);
121:   /*
122:      Problem parameters (velocity of lid, prandtl, and grashof numbers)
123:   */
124:   user.lidvelocity = 1.0/(mx*my);
125:   user.prandtl     = 1.0;
126:   user.grashof     = 1.0;

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

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

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

143:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144:      Create nonlinear solver context

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


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

159:   SNESSolve(snes,NULL,x);

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

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

171:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
172:      Free work space.  All PETSc objects should be destroyed when they
173:      are no longer needed.
174:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
175:   VecDestroy(&x);
176:   DMDestroy(&da);
177:   SNESDestroy(&snes);
178:   PetscFinalize();
179:   return ierr;
180: }

182: /* ------------------------------------------------------------------- */


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: PetscErrorCode FormInitialGuess(AppCtx *user,DM da,Vec X)
198: {
199:   PetscInt       i,j,mx,xs,ys,xm,ym;
201:   PetscReal      grashof,dx;
202:   Field          **x;

205:   grashof = user->grashof;

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

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

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

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

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

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

258:   grashof = user->grashof;
259:   prandtl = user->prandtl;
260:   lid     = user->lidvelocity;

262:   /*
263:      Define mesh intervals ratios for uniform grid.

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


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

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

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

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

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

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

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

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

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

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

352:       /* Omega */
353:       u             = x[j][i].omega;
354:       uxx           = (2.0*u - x[j][i-1].omega - x[j][i+1].omega)*hydhx;
355:       uyy           = (2.0*u - x[j-1][i].omega - x[j+1][i].omega)*hxdhy;
356:       f[j][i].omega = uxx + uyy + (vxp*(u - x[j][i-1].omega) + vxm*(x[j][i+1].omega - u))*hy +
357:                       (vyp*(u - x[j-1][i].omega) + 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*((vxp*(u - x[j][i-1].temp) + vxm*(x[j][i+1].temp - u))*hy +
365:                                             (vyp*(u - x[j-1][i].temp) + vym*(x[j+1][i].temp - u))*hx);
366:     }
367:   }

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


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

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

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

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

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

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

474:   /* Test whether we are on the left edge of the global array */
475:   if (xints == 0) {
476:     i     = 0;
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:     /* right edge */
495:     for (j=info.ys; j<info.ys+info.ym; j++) {
496:       if (B) {
497:         bjiu = b[j][i].u;
498:         bjiv = b[j][i].v;
499:       } else {
500:         bjiu = 0.0;
501:         bjiv = 0.0;
502:       }
503:       x[j][i].u = 0.0 + bjiu;
504:       x[j][i].v = 0.0 + bjiv;
505:     }
506:   }

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

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

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

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

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

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

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