Actual source code: ex19.c

petsc-3.8.4 2018-03-24
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*);

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

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

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

110:   /*
111:       Create distributed array object to manage parallel grid and vectors
112:       for principal unknowns (x) and governing residuals (f)
113:   */
114:   DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,4,4,PETSC_DECIDE,PETSC_DECIDE,4,1,0,0,&da);
115:   DMSetFromOptions(da);
116:   DMSetUp(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: /* ------------------------------------------------------------------- */

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

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

191:    Output Parameter:
192:    X - vector
193: */
194: PetscErrorCode FormInitialGuess(AppCtx *user,DM da,Vec X)
195: {
196:   PetscInt       i,j,mx,xs,ys,xm,ym;
198:   PetscReal      grashof,dx;
199:   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: }

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

253:   grashof = user->grashof;
254:   prandtl = user->prandtl;
255:   lid     = user->lidvelocity;

257:   /*
258:      Define mesh intervals ratios for uniform grid.

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


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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

430:   /* Set the boundary conditions on the momentum equations */
431:   /* Test whether we are on the bottom edge of the global array */
432:   if (yints == 0) {
433:     j     = 0;
434:     /* bottom edge */
435:     for (i=info.xs; i<info.xs+info.xm; i++) {

437:       if (B) {
438:         bjiu = b[j][i].u;
439:         bjiv = b[j][i].v;
440:       } else {
441:         bjiu = 0.0;
442:         bjiv = 0.0;
443:       }
444:       x[j][i].u = 0.0 + bjiu;
445:       x[j][i].v = 0.0 + bjiv;
446:     }
447:   }

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

466:   /* Test whether we are on the left edge of the global array */
467:   if (xints == 0) {
468:     i     = 0;
469:     /* left edge */
470:     for (j=info.ys; j<info.ys+info.ym; j++) {
471:       if (B) {
472:         bjiu = b[j][i].u;
473:         bjiv = b[j][i].v;
474:       } else {
475:         bjiu = 0.0;
476:         bjiv = 0.0;
477:       }
478:       x[j][i].u = 0.0 + bjiu;
479:       x[j][i].v = 0.0 + bjiv;
480:     }
481:   }

483:   /* Test whether we are on the right edge of the global array */
484:   if (xinte == info.mx) {
485:     i     = info.mx - 1;
486:     /* right edge */
487:     for (j=info.ys; j<info.ys+info.ym; j++) {
488:       if (B) {
489:         bjiu = b[j][i].u;
490:         bjiv = b[j][i].v;
491:       } else {
492:         bjiu = 0.0;
493:         bjiv = 0.0;
494:       }
495:       x[j][i].u = 0.0 + bjiu;
496:       x[j][i].v = 0.0 + bjiv;
497:     }
498:   }

500:   for (k=0; k < sweeps; k++) {
501:     for (j=info.ys; j<info.ys + info.ym; j++) {
502:       for (i=info.xs; i<info.xs + info.xm; i++) {
503:         ptconverged = PETSC_FALSE;
504:         pfnorm0     = 0.0;
505:         fu          = 0.0;
506:         fv          = 0.0;
507:         fomega      = 0.0;
508:         ftemp       = 0.0;
509:         for (l = 0; l < max_its && !ptconverged; l++) {
510:           if (B) {
511:             bjiu     = b[j][i].u;
512:             bjiv     = b[j][i].v;
513:             bjiomega = b[j][i].omega;
514:             bjitemp  = b[j][i].temp;
515:           } else {
516:             bjiu     = 0.0;
517:             bjiv     = 0.0;
518:             bjiomega = 0.0;
519:             bjitemp  = 0.0;
520:           }

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

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

557:             /* Temperature */
558:             u     = x[j][i].temp;
559:             uxx   = (2.0*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx;
560:             uyy   = (2.0*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy;
561:             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;
562:             dftdt = 2.0*(hydhx + hxdhy) + prandtl*((vxp - vxm)*hy + (vyp - vym)*hx);
563:             if (PetscRealPart(vx) > 0.0) dftdu = prandtl*(u - x[j][i-1].temp)*hy;
564:             else dftdu = prandtl*(x[j][i+1].temp - u)*hy;

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

569:             /* invert the system:
570:              [ dfu / du     0        0        0    ][yu] = [fu]
571:              [     0    dfv / dv     0        0    ][yv]   [fv]
572:              [ dfo / du dfo / dv dfo / do     0    ][yo]   [fo]
573:              [ dft / du dft / dv     0    dft / dt ][yt]   [ft]
574:              by simple back-substitution
575:            */
576:             yu = fu / dfudu;
577:             yv = fv / dfvdv;
578:             yo = (fomega - (dfodu*yu + dfodv*yv)) / dfodo;
579:             yt = (ftemp - (dftdu*yu + dftdv*yv)) / dftdt;

581:             x[j][i].u     = x[j][i].u - yu;
582:             x[j][i].v     = x[j][i].v - yv;
583:             x[j][i].temp  = x[j][i].temp - yt;
584:             x[j][i].omega = x[j][i].omega - yo;
585:           }
586:           if (i == 0) {
587:             fomega        = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx - bjiomega;
588:             ftemp         = x[j][i].temp - bjitemp;
589:             yo            = fomega;
590:             yt            = ftemp;
591:             x[j][i].omega = x[j][i].omega - fomega;
592:             x[j][i].temp  = x[j][i].temp - ftemp;
593:           }
594:           if (i == info.mx - 1) {
595:             fomega        = x[j][i].omega - (x[j][i].v - x[j][i-1].v)*dhx - bjiomega;
596:             ftemp         = x[j][i].temp - (PetscReal)(grashof>0) - 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 (j == 0) {
603:             fomega        = x[j][i].omega + (x[j+1][i].u - x[j][i].u)*dhy - bjiomega;
604:             ftemp         = x[j][i].temp-x[j+1][i].temp - 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 == info.my - 1) {
611:             fomega        = x[j][i].omega + (x[j][i].u - x[j-1][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:           tot_its++;
619:           pfnorm = PetscRealPart(fu*fu + fv*fv + fomega*fomega + ftemp*ftemp);
620:           pfnorm = PetscSqrtReal(pfnorm);
621:           pynorm = PetscRealPart(yu*yu + yv*yv + yo*yo + yt*yt);
622:           pynorm = PetscSqrtReal(pynorm);
623:           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);
624:           pxnorm = PetscSqrtReal(pxnorm);
625:           if (l == 0) pfnorm0 = pfnorm;
626:           if (rtol*pfnorm0 >pfnorm ||
627:               atol > pfnorm ||
628:               pxnorm*stol > pynorm) ptconverged = PETSC_TRUE;
629:         }
630:       }
631:     }
632:   }
633:   DMDAVecRestoreArray(da,localX,&x);
634:   if (B) {
635:     DMDAVecRestoreArrayRead(da,localB,&b);
636:   }
637:   DMLocalToGlobalBegin(da,localX,INSERT_VALUES,X);
638:   DMLocalToGlobalEnd(da,localX,INSERT_VALUES,X);
639:   PetscLogFlops(tot_its*(84.0 + 41.0 + 26.0));
640:   DMRestoreLocalVector(da,&localX);
641:   if (B) {
642:     DMRestoreLocalVector(da,&localB);
643:   }
644:   return(0);
645: }