Actual source code: ex27.c

  2: static char help[] = "Nonlinear driven cavity with multigrid and pusedo timestepping 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";

 11: /*T
 12:    Concepts: SNES^solving a system of nonlinear equations (parallel multicomponent example);
 13:    Concepts: DA^using distributed arrays;
 14:    Concepts: multicomponent
 15:    Processors: n
 16: T*/

 18: /* ------------------------------------------------------------------------

 20:     We thank David E. Keyes for contributing the driven cavity discretization
 21:     within this example code.

 23:     This problem is modeled by the partial differential equation system
 24:   
 25:         dU/dt          - Lap(U)     - Grad_y(Omega) = 0
 26:         dV/dt          - Lap(V)     + Grad_x(Omega) = 0
 27:         d(omega)/dt    - Lap(Omega) + Div([U*Omega,V*Omega]) - GR*Grad_x(T) = 0
 28:         dT/dt          - Lap(T)     + PR*Div([U*T,V*T]) = 0

 30:     in the unit square, which is uniformly discretized in each of x and
 31:     y in this simple encoding.

 33:     No-slip, rigid-wall Dirichlet conditions are used for [U,V].
 34:     Dirichlet conditions are used for Omega, based on the definition of
 35:     vorticity: Omega = - Grad_y(U) + Grad_x(V), where along each
 36:     constant coordinate boundary, the tangential derivative is zero.
 37:     Dirichlet conditions are used for T on the left and right walls,
 38:     and insulation homogeneous Neumann conditions are used for T on
 39:     the top and bottom walls. 

 41:     A finite difference approximation with the usual 5-point stencil 
 42:     is used to discretize the boundary value problem to obtain a 
 43:     nonlinear system of equations.  Upwinding is used for the divergence
 44:     (convective) terms and central for the gradient (source) terms.
 45:     
 46:     The Jacobian can be either
 47:       * formed via finite differencing using coloring (the default), or
 48:       * applied matrix-free via the option -snes_mf 
 49:         (for larger grid problems this variant may not converge 
 50:         without a preconditioner due to ill-conditioning).

 52:   ------------------------------------------------------------------------- */

 54: /* 
 55:    Include "petscda.h" so that we can use distributed arrays (DAs).
 56:    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
 57:    file automatically includes:
 58:      petsc.h       - base PETSc routines   petscvec.h - vectors
 59:      petscsys.h    - system routines       petscmat.h - matrices
 60:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 61:      petscviewer.h - viewers               petscpc.h  - preconditioners
 62:      petscksp.h   - linear solvers 
 63: */
 64:  #include petscsnes.h
 65:  #include petscda.h

 67: /* 
 68:    User-defined routines and data structures
 69: */

 71: typedef struct {
 72:   PassiveScalar  fnorm_ini,dt_ini,cfl_ini;
 73:   PassiveScalar  fnorm,dt,cfl;
 74:   PassiveScalar  ptime;
 75:   PassiveScalar  cfl_max,max_time;
 76:   PassiveScalar  fnorm_ratio;
 77:   PetscInt       ires,iramp,itstep;
 78:   PetscInt       max_steps,print_freq;
 79:   PetscInt       LocalTimeStepping;
 80:   PetscTruth     use_parabolic;
 81: } TstepCtx;

 83: /*
 84:     The next two structures are essentially the same. The difference is that
 85:   the first does not contain derivative information (as used by ADIC) while the
 86:   second does. The first is used only to contain the solution at the previous time-step
 87:   which is a constant in the computation of the current time step and hence passive 
 88:   (meaning not active in terms of derivatives).
 89: */
 90: typedef struct {
 91:   PassiveScalar u,v,omega,temp;
 92: } PassiveField;

 94: typedef struct {
 95:   PetscScalar u,v,omega,temp;
 96: } Field;


 99: typedef struct {
100:   PetscInt     mglevels;
101:   PetscInt     cycles;                       /* numbers of time steps for integration */
102:   PassiveReal  lidvelocity,prandtl,grashof;  /* physical parameters */
103:   PetscTruth   draw_contours;                /* indicates drawing contours of solution */
104:   PetscTruth   PreLoading;
105: } Parameter;

107: typedef struct {
108:   Vec          Xold,func;
109:   TstepCtx     *tsCtx;
110:   Parameter    *param;
111: } AppCtx;



124: int main(int argc,char **argv)
125: {
126:   DMMG           *dmmg;               /* multilevel grid structure */
127:   AppCtx         *user;                /* user-defined work context */
128:   TstepCtx       tsCtx;
129:   Parameter      param;
130:   PetscInt       mx,my,i;
132:   MPI_Comm       comm;
133:   DA             da;

135:   PetscInitialize(&argc,&argv,(char *)0,help);
136:   comm = PETSC_COMM_WORLD;


139:   PreLoadBegin(PETSC_TRUE,"SetUp");

141:     param.PreLoading = PreLoading;
142:     DMMGCreate(comm,2,&user,&dmmg);
143:     param.mglevels = DMMGGetLevels(dmmg);


146:     /*
147:       Create distributed array multigrid object (DMMG) to manage parallel grid and vectors
148:       for principal unknowns (x) and governing residuals (f)
149:     */
150:     DACreate2d(comm,DA_NONPERIODIC,DA_STENCIL_STAR,-4,-4,PETSC_DECIDE,PETSC_DECIDE,4,1,0,0,&da);
151:     DMMGSetDM(dmmg,(DM)da);
152:     DADestroy(da);

154:     DAGetInfo(DMMGGetDA(dmmg),0,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
155:                      PETSC_IGNORE,PETSC_IGNORE);
156:     /* 
157:      Problem parameters (velocity of lid, prandtl, and grashof numbers)
158:     */
159:     param.lidvelocity = 1.0/(mx*my);
160:     param.prandtl     = 1.0;
161:     param.grashof     = 1.0;
162:     PetscOptionsGetReal(PETSC_NULL,"-lidvelocity",&param.lidvelocity,PETSC_NULL);
163:     PetscOptionsGetReal(PETSC_NULL,"-prandtl",&param.prandtl,PETSC_NULL);
164:     PetscOptionsGetReal(PETSC_NULL,"-grashof",&param.grashof,PETSC_NULL);
165:     PetscOptionsHasName(PETSC_NULL,"-contours",&param.draw_contours);

167:     DASetFieldName(DMMGGetDA(dmmg),0,"x-velocity");
168:     DASetFieldName(DMMGGetDA(dmmg),1,"y-velocity");
169:     DASetFieldName(DMMGGetDA(dmmg),2,"Omega");
170:     DASetFieldName(DMMGGetDA(dmmg),3,"temperature");

172:     /*======================================================================*/
173:     /* Initilize stuff related to time stepping */
174:     /*======================================================================*/
175:     tsCtx.fnorm_ini = 0.0;  tsCtx.cfl_ini     = 50.0;    tsCtx.cfl_max = 1.0e+06;
176:     tsCtx.max_steps = 50;   tsCtx.max_time    = 1.0e+12; tsCtx.iramp   = -50;
177:     tsCtx.dt        = 0.01; tsCtx.fnorm_ratio = 1.0e+10;
178:     tsCtx.LocalTimeStepping = 0;
179:     tsCtx.use_parabolic     = PETSC_FALSE;
180:     PetscOptionsGetLogical(PETSC_NULL,"-use_parabolic",&tsCtx.use_parabolic,PETSC_IGNORE);
181:     PetscOptionsGetInt(PETSC_NULL,"-max_st",&tsCtx.max_steps,PETSC_NULL);
182:     PetscOptionsGetReal(PETSC_NULL,"-ts_rtol",&tsCtx.fnorm_ratio,PETSC_NULL);
183:     PetscOptionsGetReal(PETSC_NULL,"-cfl_ini",&tsCtx.cfl_ini,PETSC_NULL);
184:     PetscOptionsGetReal(PETSC_NULL,"-cfl_max",&tsCtx.cfl_max,PETSC_NULL);
185:     tsCtx.print_freq = tsCtx.max_steps;
186:     PetscOptionsGetInt(PETSC_NULL,"-print_freq",&tsCtx.print_freq,PETSC_NULL);
187: 
188:     PetscMalloc(param.mglevels*sizeof(AppCtx),&user);
189:     for (i=0; i<param.mglevels; i++) {
190:       VecDuplicate(dmmg[i]->x, &(user[i].Xold));
191:       VecDuplicate(dmmg[i]->x, &(user[i].func));
192:       user[i].tsCtx = &tsCtx;
193:       user[i].param = &param;
194:       dmmg[i]->user = &user[i];
195:     }
196:     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
197:        Create user context, set problem data, create vector data structures.
198:        Also, compute the initial guess.
199:        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
200: 
201:     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
202:        Create nonlinear solver context
203:        
204:        Process adiC(20):  AddTSTermLocal FormFunctionLocal FormFunctionLocali
205:        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
206:     DMMGSetSNESLocal(dmmg,FormFunctionLocal,0,ad_FormFunctionLocal,admf_FormFunctionLocal);
207:     DMMGSetSNESLocali(dmmg,FormFunctionLocali,ad_FormFunctionLocali,admf_FormFunctionLocali);
208: 
209:     PetscPrintf(comm,"lid velocity = %g, prandtl # = %g, grashof # = %g\n",
210:                        param.lidvelocity,param.prandtl,param.grashof);
211: 
212: 
213:     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
214:        Solve the nonlinear system
215:        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
216:     DMMGSetInitialGuess(dmmg,FormInitialGuess);
217: 
218:     PreLoadStage("Solve");
219:     Initialize(dmmg);
220:     Update(dmmg);
221:     /*
222:       Visualize solution
223:     */
224: 
225:     if (param.draw_contours) {
226:       VecView(DMMGGetx(dmmg),PETSC_VIEWER_DRAW_WORLD);
227:     }
228:     /*VecView(DMMGGetx(dmmg),PETSC_VIEWER_STDOUT_WORLD);*/
229:     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
230:        Free work space.  All PETSc objects should be destroyed when they
231:        are no longer needed.
232:        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
233: 
234:     for (i=0; i<param.mglevels; i++) {
235:       VecDestroy(user[i].Xold);
236:       VecDestroy(user[i].func);
237:     }
238:     PetscFree(user);
239:     DMMGDestroy(dmmg);
240:     PreLoadEnd();
241: 
242:   PetscFinalize();
243:   return 0;
244: }

246: /* ------------------------------------------------------------------- */
249: /* ------------------------------------------------------------------- */
250: PetscErrorCode Initialize(DMMG *dmmg)
251: {
252:   AppCtx         *user = (AppCtx*)dmmg[0]->user;
253:   DA             da;
254:   Parameter      *param = user->param;
255:   PetscInt       i,j,mx,xs,ys,xm,ym,mglevel;
257:   PetscReal      grashof,dx;
258:   Field          **x;

260:   da = (DA)(dmmg[param->mglevels-1]->dm);
261:   grashof = user->param->grashof;

263:   DAGetInfo(da,0,&mx,0,0,0,0,0,0,0,0,0);
264:   dx  = 1.0/(mx-1);

266:   /*
267:      Get local grid boundaries (for 2-dimensional DA):
268:        xs, ys   - starting grid indices (no ghost points)
269:        xm, ym   - widths of local grid (no ghost points)
270:   */
271:   DAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);

273:   /*
274:      Get a pointer to vector data.
275:        - For default PETSc vectors, VecGetArray() returns a pointer to
276:          the data array.  Otherwise, the routine is implementation dependent.
277:        - You MUST call VecRestoreArray() when you no longer need access to
278:          the array.
279:   */
280:   DAVecGetArray(da,((AppCtx*)dmmg[param->mglevels-1]->user)->Xold,&x);

282:   /*
283:      Compute initial guess over the locally owned part of the grid
284:      Initial condition is motionless fluid and equilibrium temperature
285:   */
286:   for (j=ys; j<ys+ym; j++) {
287:     for (i=xs; i<xs+xm; i++) {
288:       x[j][i].u     = 0.0;
289:       x[j][i].v     = 0.0;
290:       x[j][i].omega = 0.0;
291:       x[j][i].temp  = (grashof>0)*i*dx;
292:     }
293:   }

295:   /*
296:      Restore vector
297:   */
298:   DAVecRestoreArray(da,((AppCtx*)dmmg[param->mglevels-1]->user)->Xold,&x);
299:   /* Restrict Xold to coarser levels */
300:   for (mglevel=param->mglevels-1; mglevel>0; mglevel--) {
301:     MatRestrict(dmmg[mglevel]->R, ((AppCtx*)dmmg[mglevel]->user)->Xold, ((AppCtx*)dmmg[mglevel-1]->user)->Xold);
302:     VecPointwiseMult(dmmg[mglevel]->Rscale,((AppCtx*)dmmg[mglevel-1]->user)->Xold,((AppCtx*)dmmg[mglevel-1]->user)->Xold);
303:   }
304: 
305:   /* Store X in the qold for time stepping */
306:   /*VecDuplicate(X,&tsCtx->qold);
307:   VecDuplicate(X,&tsCtx->func);
308:   VecCopy(X,tsCtx->Xold);
309:   tsCtx->ires = 0;
310:   SNESComputeFunction(snes,tsCtx->Xold,tsCtx->func);
311:   VecNorm(tsCtx->func,NORM_2,&tsCtx->fnorm_ini);
312:   tsCtx->ires = 1;
313:   PetscPrintf(PETSC_COMM_WORLD,"Initial function norm is %g\n",tsCtx->fnorm_ini);*/
314:   return 0;
315: }
316: /* ------------------------------------------------------------------- */
319: /* 
320:    FormInitialGuess - Forms initial approximation.

322:    Input Parameters:
323:    user - user-defined application context
324:    X - vector

326:    Output Parameter:
327:    X - vector
328:  */
329: PetscErrorCode FormInitialGuess(DMMG dmmg,Vec X)
330: {
331:   AppCtx         *user = (AppCtx*)dmmg->user;
332:   TstepCtx       *tsCtx = user->tsCtx;

335:   VecCopy(user->Xold, X);
336:   /* calculate the residual on fine mesh */
337:   if (user->tsCtx->fnorm_ini == 0.0) {
338:     tsCtx->ires = 0;
339:     SNESComputeFunction(dmmg->snes,user->Xold,user->func);
340:     VecNorm(user->func,NORM_2,&tsCtx->fnorm_ini);
341:     tsCtx->ires = 1;
342:     tsCtx->cfl = tsCtx->cfl_ini;
343:   }
344:   return 0;
345: }

347: /*---------------------------------------------------------------------*/
350: PetscErrorCode AddTSTermLocal(DALocalInfo* info,Field **x,Field **f,void *ptr)
351: /*---------------------------------------------------------------------*/
352: {
353:   AppCtx         *user = (AppCtx*)ptr;
354:   TstepCtx       *tsCtx = user->tsCtx;
355:   DA             da = info->da;
357:   PetscInt       i,j, xints,xinte,yints,yinte;
358:   PetscReal      hx,hy,dhx,dhy,hxhy;
359:   PassiveScalar  dtinv;
360:   PassiveField   **xold;
361:   PetscTruth     use_parab = tsCtx->use_parabolic;

364:   xints = info->xs; xinte = info->xs+info->xm; yints = info->ys; yinte = info->ys+info->ym;
365:   dhx = (PetscReal)(info->mx-1);  dhy = (PetscReal)(info->my-1);
366:   hx = 1.0/dhx;                   hy = 1.0/dhy;
367:   hxhy = hx*hy;
368:   DAVecGetArray(da,user->Xold,&xold);
369:   dtinv = hxhy/(tsCtx->cfl*tsCtx->dt);
370:   /* 
371:      use_parab = PETSC_TRUE for parabolic equations; all the four equations have temporal term.
372:                = PETSC_FALSE for differential algebraic equtions (DAE); 
373:                  velocity equations do not have temporal term.
374:   */
375:   for (j=yints; j<yinte; j++) {
376:     for (i=xints; i<xinte; i++) {
377:       if (use_parab) {
378:         f[j][i].u     += dtinv*(x[j][i].u-xold[j][i].u);
379:         f[j][i].v     += dtinv*(x[j][i].v-xold[j][i].v);
380:       }
381:       f[j][i].omega += dtinv*(x[j][i].omega-xold[j][i].omega);
382:       f[j][i].temp  += dtinv*(x[j][i].temp-xold[j][i].temp);
383:     }
384:   }
385:   DAVecRestoreArray(da,user->Xold,&xold);
386:   return(0);
387: }

391: PetscErrorCode FormFunctionLocal(DALocalInfo *info,Field **x,Field **f,void *ptr)
392:  {
393:   AppCtx         *user = (AppCtx*)ptr;
394:   TstepCtx       *tsCtx = user->tsCtx;
396:   PetscInt       xints,xinte,yints,yinte,i,j;
397:   PetscReal      hx,hy,dhx,dhy,hxdhy,hydhx;
398:   PetscReal      grashof,prandtl,lid;
399:   PetscScalar    u,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;

402:   grashof = user->param->grashof;
403:   prandtl = user->param->prandtl;
404:   lid     = user->param->lidvelocity;

406:   /* 
407:      Define mesh intervals ratios for uniform grid.
408:      [Note: FD formulae below are normalized by multiplying through by
409:      local volume element to obtain coefficients O(1) in two dimensions.]
410:   */
411:   dhx = (PetscReal)(info->mx-1);  dhy = (PetscReal)(info->my-1);
412:   hx = 1.0/dhx;                   hy = 1.0/dhy;
413:   hxdhy = hx*dhy;                 hydhx = hy*dhx;

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

417:   /* Test whether we are on the bottom edge of the global array */
418:   if (yints == 0) {
419:     j = 0;
420:     yints = yints + 1;
421:     /* bottom edge */
422:     for (i=info->xs; i<info->xs+info->xm; i++) {
423:         f[j][i].u     = x[j][i].u;
424:         f[j][i].v     = x[j][i].v;
425:         f[j][i].omega = x[j][i].omega + (x[j+1][i].u - x[j][i].u)*dhy;
426:         f[j][i].temp  = x[j][i].temp-x[j+1][i].temp;
427:     }
428:   }

430:   /* Test whether we are on the top edge of the global array */
431:   if (yinte == info->my) {
432:     j = info->my - 1;
433:     yinte = yinte - 1;
434:     /* top edge */
435:     for (i=info->xs; i<info->xs+info->xm; i++) {
436:         f[j][i].u     = x[j][i].u - lid;
437:         f[j][i].v     = x[j][i].v;
438:         f[j][i].omega = x[j][i].omega + (x[j][i].u - x[j-1][i].u)*dhy;
439:         f[j][i].temp  = x[j][i].temp-x[j-1][i].temp;
440:     }
441:   }

443:   /* Test whether we are on the left edge of the global array */
444:   if (xints == 0) {
445:     i = 0;
446:     xints = xints + 1;
447:     /* left edge */
448:     for (j=info->ys; j<info->ys+info->ym; j++) {
449:       f[j][i].u     = x[j][i].u;
450:       f[j][i].v     = x[j][i].v;
451:       f[j][i].omega = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx;
452:       f[j][i].temp  = x[j][i].temp;
453:     }
454:   }

456:   /* Test whether we are on the right edge of the global array */
457:   if (xinte == info->mx) {
458:     i = info->mx - 1;
459:     xinte = xinte - 1;
460:     /* right edge */
461:     for (j=info->ys; j<info->ys+info->ym; j++) {
462:       f[j][i].u     = x[j][i].u;
463:       f[j][i].v     = x[j][i].v;
464:       f[j][i].omega = x[j][i].omega - (x[j][i].v - x[j][i-1].v)*dhx;
465:       f[j][i].temp  = x[j][i].temp - (PetscReal)(grashof>0);
466:     }
467:   }

469:   /* Compute over the interior points */
470:   for (j=yints; j<yinte; j++) {
471:     for (i=xints; i<xinte; i++) {

473:         /*
474:           convective coefficients for upwinding
475:         */
476:         vx = x[j][i].u; avx = PetscAbsScalar(vx);
477:         vxp = .5*(vx+avx); vxm = .5*(vx-avx);
478:         vy = x[j][i].v; avy = PetscAbsScalar(vy);
479:         vyp = .5*(vy+avy); vym = .5*(vy-avy);

481:         /* U velocity */
482:         u          = x[j][i].u;
483:         uxx        = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
484:         uyy        = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
485:         f[j][i].u  = uxx + uyy - .5*(x[j+1][i].omega-x[j-1][i].omega)*hx;

487:         /* V velocity */
488:         u          = x[j][i].v;
489:         uxx        = (2.0*u - x[j][i-1].v - x[j][i+1].v)*hydhx;
490:         uyy        = (2.0*u - x[j-1][i].v - x[j+1][i].v)*hxdhy;
491:         f[j][i].v  = uxx + uyy + .5*(x[j][i+1].omega-x[j][i-1].omega)*hy;

493:         /* Omega */
494:         u          = x[j][i].omega;
495:         uxx        = (2.0*u - x[j][i-1].omega - x[j][i+1].omega)*hydhx;
496:         uyy        = (2.0*u - x[j-1][i].omega - x[j+1][i].omega)*hxdhy;
497:         f[j][i].omega = uxx + uyy +
498:                         (vxp*(u - x[j][i-1].omega) +
499:                           vxm*(x[j][i+1].omega - u)) * hy +
500:                         (vyp*(u - x[j-1][i].omega) +
501:                           vym*(x[j+1][i].omega - u)) * hx -
502:                         .5 * grashof * (x[j][i+1].temp - x[j][i-1].temp) * hy;

504:         /* Temperature */
505:         u             = x[j][i].temp;
506:         uxx           = (2.0*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx;
507:         uyy           = (2.0*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy;
508:         f[j][i].temp =  uxx + uyy  + prandtl * (
509:                         (vxp*(u - x[j][i-1].temp) +
510:                           vxm*(x[j][i+1].temp - u)) * hy +
511:                         (vyp*(u - x[j-1][i].temp) +
512:                                  vym*(x[j+1][i].temp - u)) * hx);
513:     }
514:   }

516:   /* Add time step contribution */
517:   if (tsCtx->ires) {
518:     AddTSTermLocal(info,x,f,ptr);
519:   }
520:   /*
521:      Flop count (multiply-adds are counted as 2 operations)
522:   */
523:   PetscLogFlops(84*info->ym*info->xm);
524:   return(0);
525: }

529: PetscErrorCode FormFunctionLocali(DALocalInfo *info,MatStencil *st,Field **x,PetscScalar *f,void *ptr)
530:  {
531:   AppCtx      *user = (AppCtx*)ptr;
532:   PetscInt    i,j,c;
533:   PassiveReal hx,hy,dhx,dhy,hxdhy,hydhx;
534:   PassiveReal grashof,prandtl,lid;
535:   PetscScalar u,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;

538:   grashof = user->param->grashof;
539:   prandtl = user->param->prandtl;
540:   lid     = user->param->lidvelocity;

542:   /* 
543:      Define mesh intervals ratios for uniform grid.
544:      [Note: FD formulae below are normalized by multiplying through by
545:      local volume element to obtain coefficients O(1) in two dimensions.]
546:   */
547:   dhx = (PetscReal)(info->mx-1);     dhy = (PetscReal)(info->my-1);
548:   hx = 1.0/dhx;                   hy = 1.0/dhy;
549:   hxdhy = hx*dhy;                 hydhx = hy*dhx;

551:   i = st->i; j = st->j; c = st->c;

553:   /* Test whether we are on the right edge of the global array */
554:   if (i == info->mx-1) {
555:     if (c == 0) *f     = x[j][i].u;
556:     else if (c == 1) *f     = x[j][i].v;
557:     else if (c == 2) *f = x[j][i].omega - (x[j][i].v - x[j][i-1].v)*dhx;
558:     else *f  = x[j][i].temp - (PetscReal)(grashof>0);

560:   /* Test whether we are on the left edge of the global array */
561:   } else if (i == 0) {
562:     if (c == 0) *f     = x[j][i].u;
563:     else if (c == 1) *f     = x[j][i].v;
564:     else if (c == 2) *f = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx;
565:     else *f  = x[j][i].temp;

567:   /* Test whether we are on the top edge of the global array */
568:   } else if (j == info->my-1) {
569:     if (c == 0) *f     = x[j][i].u - lid;
570:     else if (c == 1) *f     = x[j][i].v;
571:     else if (c == 2) *f = x[j][i].omega + (x[j][i].u - x[j-1][i].u)*dhy;
572:     else *f  = x[j][i].temp-x[j-1][i].temp;

574:   /* Test whether we are on the bottom edge of the global array */
575:   } else if (j == 0) {
576:     if (c == 0) *f     = x[j][i].u;
577:     else if (c == 1) *f     = x[j][i].v;
578:     else if (c == 2) *f = x[j][i].omega + (x[j+1][i].u - x[j][i].u)*dhy;
579:     else *f  = x[j][i].temp-x[j+1][i].temp;

581:   /* Compute over the interior points */
582:   } else {
583:     /*
584:       convective coefficients for upwinding
585:     */
586:     vx = x[j][i].u; avx = PetscAbsScalar(vx);
587:     vxp = .5*(vx+avx); vxm = .5*(vx-avx);
588:     vy = x[j][i].v; avy = PetscAbsScalar(vy);
589:     vyp = .5*(vy+avy); vym = .5*(vy-avy);

591:     /* U velocity */
592:     if (c == 0) {
593:       u          = x[j][i].u;
594:       uxx        = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
595:       uyy        = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
596:       *f         = uxx + uyy - .5*(x[j+1][i].omega-x[j-1][i].omega)*hx;

598:     /* V velocity */
599:     } else if (c == 1) {
600:       u          = x[j][i].v;
601:       uxx        = (2.0*u - x[j][i-1].v - x[j][i+1].v)*hydhx;
602:       uyy        = (2.0*u - x[j-1][i].v - x[j+1][i].v)*hxdhy;
603:       *f         = uxx + uyy + .5*(x[j][i+1].omega-x[j][i-1].omega)*hy;
604: 
605:     /* Omega */
606:     } else if (c == 2) {
607:       u          = x[j][i].omega;
608:       uxx        = (2.0*u - x[j][i-1].omega - x[j][i+1].omega)*hydhx;
609:       uyy        = (2.0*u - x[j-1][i].omega - x[j+1][i].omega)*hxdhy;
610:       *f         = uxx + uyy +
611:         (vxp*(u - x[j][i-1].omega) +
612:          vxm*(x[j][i+1].omega - u)) * hy +
613:         (vyp*(u - x[j-1][i].omega) +
614:          vym*(x[j+1][i].omega - u)) * hx -
615:         .5 * grashof * (x[j][i+1].temp - x[j][i-1].temp) * hy;
616: 
617:     /* Temperature */
618:     } else {
619:       u           = x[j][i].temp;
620:       uxx         = (2.0*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx;
621:       uyy         = (2.0*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy;
622:       *f          =  uxx + uyy  + prandtl * (
623:                                              (vxp*(u - x[j][i-1].temp) +
624:                                               vxm*(x[j][i+1].temp - u)) * hy +
625:                                              (vyp*(u - x[j-1][i].temp) +
626:                                               vym*(x[j+1][i].temp - u)) * hx);
627:     }
628:   }

630:   return(0);
631: }


634: /*---------------------------------------------------------------------*/
637: PetscErrorCode Update(DMMG *dmmg)
638: /*---------------------------------------------------------------------*/
639: {
640:  AppCtx         *user = (AppCtx *) ((dmmg[0])->user);
641:  TstepCtx         *tsCtx = user->tsCtx;
642:  Parameter      *param = user->param;
643:  SNES           snes;
645:  PetscScalar         fratio;
646:  PetscScalar         time1,time2,cpuloc;
647:  PetscInt       max_steps,its;
648:  PetscTruth     print_flag = PETSC_FALSE;
649:  PetscInt       nfailsCum = 0,nfails = 0;


653:   PetscOptionsHasName(PETSC_NULL,"-print",&print_flag);
654:   if (user->param->PreLoading)
655:    max_steps = 1;
656:   else
657:    max_steps = tsCtx->max_steps;
658:   fratio = 1.0;
659: 
660:   PetscGetTime(&time1);
661:   for (tsCtx->itstep = 0; (tsCtx->itstep < max_steps) &&
662:          (fratio <= tsCtx->fnorm_ratio); tsCtx->itstep++) {
663:     DMMGSolve(dmmg);
664:     snes = DMMGGetSNES(dmmg);
665:     SNESGetIterationNumber(snes,&its);
666:     PetscPrintf(PETSC_COMM_WORLD,"Number of Newton iterations = %D\n", its);
667:     SNESGetNumberUnsuccessfulSteps(snes,&nfails);
668:     nfailsCum += nfails; nfails = 0;
669:     if (nfailsCum >= 2) SETERRQ(1,"Unable to find a Newton Step");
670:     /*tsCtx->qcur = DMMGGetx(dmmg);
671:       VecCopy(tsCtx->qcur,tsCtx->qold);*/

673:     VecCopy(dmmg[param->mglevels-1]->x, ((AppCtx*)dmmg[param->mglevels-1]->user)->Xold);
674:     for (its=param->mglevels-1; its>0 ;its--) {
675:       MatRestrict(dmmg[its]->R, ((AppCtx*)dmmg[its]->user)->Xold, ((AppCtx*)dmmg[its-1]->user)->Xold);
676:       VecPointwiseMult(dmmg[its]->Rscale,((AppCtx*)dmmg[its-1]->user)->Xold,((AppCtx*)dmmg[its-1]->user)->Xold);
677:     }

679: 
680:     ComputeTimeStep(snes,((AppCtx*)dmmg[param->mglevels-1]->user));
681:     if (print_flag) {
682:       PetscPrintf(PETSC_COMM_WORLD,"At Time Step %D cfl = %g and fnorm = %g\n",
683:                          tsCtx->itstep,tsCtx->cfl,tsCtx->fnorm);
684:       PetscPrintf(PETSC_COMM_WORLD,"Wall clock time needed %g seconds for %D time steps\n",
685:                          cpuloc,tsCtx->itstep);
686:     }
687:     fratio = tsCtx->fnorm_ini/tsCtx->fnorm;
688:     PetscGetTime(&time2);
689:     cpuloc = time2-time1;
690:     MPI_Barrier(PETSC_COMM_WORLD);
691:   } /* End of time step loop */
692: 
693:   PetscPrintf(PETSC_COMM_WORLD,"Total wall clock time needed %g seconds for %D time steps\n",
694:                      cpuloc,tsCtx->itstep);
695:   PetscPrintf(PETSC_COMM_WORLD,"cfl = %g fnorm = %g\n",tsCtx->cfl,tsCtx->fnorm);
696:   if (user->param->PreLoading) {
697:     tsCtx->fnorm_ini = 0.0;
698:     PetscPrintf(PETSC_COMM_WORLD,"Preloading done ...\n");
699:   }
700:   /*
701:   {
702:     Vec xx,yy;
703:     PetscScalar fnorm,fnorm1;
704:     SNESGetFunctionNorm(snes,&fnorm);
705:     xx = DMMGGetx(dmmg);
706:     VecDuplicate(xx,&yy);
707:     SNESComputeFunction(snes,xx,yy);
708:     VecNorm(yy,NORM_2,&fnorm1);
709:     PetscPrintf(PETSC_COMM_WORLD,"fnorm = %g, fnorm1 = %g\n",fnorm,fnorm1);
710:     
711:   }
712:   */

714:   return(0);
715: }

717: /*---------------------------------------------------------------------*/
720: PetscErrorCode ComputeTimeStep(SNES snes,void *ptr)
721: /*---------------------------------------------------------------------*/
722: {
723:   AppCtx         *user = (AppCtx*)ptr;
724:   TstepCtx       *tsCtx = user->tsCtx;
725:   Vec                 func = user->func;
726:   PetscScalar    inc = 1.1,  newcfl;
728:   /*int               iramp = tsCtx->iramp;*/
729: 
731:   tsCtx->dt = 0.01;
732:   PetscOptionsGetScalar(PETSC_NULL,"-deltat",&tsCtx->dt,PETSC_NULL);
733:   tsCtx->ires = 0;
734:   SNESComputeFunction(snes,user->Xold,user->func);
735:   tsCtx->ires = 1;
736:   VecNorm(func,NORM_2,&tsCtx->fnorm);
737:   newcfl     = inc*tsCtx->cfl_ini*tsCtx->fnorm_ini/tsCtx->fnorm;
738:   tsCtx->cfl = PetscMin(newcfl,tsCtx->cfl_max);
739:   /* first time through so compute initial function norm */
740:   /*if (tsCtx->fnorm_ini == 0.0) {
741:     tsCtx->fnorm_ini = tsCtx->fnorm;
742:     tsCtx->cfl       = tsCtx->cfl_ini;
743:   } else {
744:     newcfl     = inc*tsCtx->cfl_ini*tsCtx->fnorm_ini/tsCtx->fnorm;
745:     tsCtx->cfl = PetscMin(newcfl,tsCtx->cfl_max);
746:     }*/
747:   return(0);
748: }