Actual source code: ex19.c

petsc-master 2014-12-21
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;

204:   grashof = user->grashof;

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

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

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

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

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

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

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

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

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


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

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

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

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

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

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

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

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

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

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

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

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

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


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

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

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

437:   /* Set the boundary conditions on the momentum equations */
438:   /* Test whether we are on the bottom edge of the global array */
439:   if (yints == 0) {
440:     j     = 0;
441:     yints = yints + 1;
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:     yinte = yinte - 1;
461:     /* top edge */
462:     for (i=info.xs; i<info.xs+info.xm; i++) {
463:       if (B) {
464:         bjiu = b[j][i].u;
465:         bjiv = b[j][i].v;
466:       } else {
467:         bjiu = 0.0;
468:         bjiv = 0.0;
469:       }
470:       x[j][i].u = lid + bjiu;
471:       x[j][i].v = bjiv;
472:     }
473:   }

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

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

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

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

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

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

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

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

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