Actual source code: ex19.c

petsc-3.5.4 2015-05-23
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 <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 <petscdm.h>
 72: #include <petscdmda.h>
 73: #endif

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

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

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

 89: extern PetscErrorCode FormInitialGuess(AppCtx*,DM,Vec);
 90: 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:   SNESSetDM(snes,(DM)da);
116:   SNESSetNGS(snes, NonlinearGS, (void*)&user);

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

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

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

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

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

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


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

158:   SNESSolve(snes,NULL,x);

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

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

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

181: /* ------------------------------------------------------------------- */


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

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

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

203:   grashof = user->grashof;

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

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

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

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

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

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

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

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

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


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

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

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

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

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

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

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

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

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

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

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

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

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


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

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

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

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

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

456:   /* Test whether we are on the top edge of the global array */
457:   if (yinte == info.my) {
458:     j     = info.my - 1;
459:     yinte = yinte - 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:     xints = xints + 1;
478:     /* left edge */
479:     for (j=info.ys; j<info.ys+info.ym; j++) {
480:       if (B) {
481:         bjiu = b[j][i].u;
482:         bjiv = b[j][i].v;
483:       } else {
484:         bjiu = 0.0;
485:         bjiv = 0.0;
486:       }
487:       x[j][i].u = 0.0 + bjiu;
488:       x[j][i].v = 0.0 + bjiv;
489:     }
490:   }

492:   /* Test whether we are on the right edge of the global array */
493:   if (xinte == info.mx) {
494:     i     = info.mx - 1;
495:     xinte = xinte - 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:         pfnorm      = 0.0;
516:         fu          = 0.0;
517:         fv          = 0.0;
518:         fomega      = 0.0;
519:         ftemp       = 0.0;
520:         for (l = 0; l < max_its && !ptconverged; l++) {
521:           if (B) {
522:             bjiu     = b[j][i].u;
523:             bjiv     = b[j][i].v;
524:             bjiomega = b[j][i].omega;
525:             bjitemp  = b[j][i].temp;
526:           } else {
527:             bjiu     = 0.0;
528:             bjiv     = 0.0;
529:             bjiomega = 0.0;
530:             bjitemp  = 0.0;
531:           }

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

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

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

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

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

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