Actual source code: ex19.c

  1: /*$Id: ex19.c,v 1.30 2001/08/07 21:31:17 bsmith Exp $*/

  3: static char help[] = "Nonlinear driven cavity with multigrid in 2d.n
  4:   n
  5: The 2D driven cavity problem is solved in a velocity-vorticity formulation.n
  6: The flow can be driven with the lid or with bouyancy or both:n
  7:   -lidvelocity <lid>, where <lid> = dimensionless velocity of lidn
  8:   -grashof <gr>, where <gr> = dimensionless temperature gradentn
  9:   -prandtl <pr>, where <pr> = dimensionless thermal/momentum diffusity ration
 10:   -contours : draw contour plots of solutionnn";

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

 19: /* ------------------------------------------------------------------------

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

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

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

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

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

 53:   ------------------------------------------------------------------------- */

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

 68: /* 
 69:    User-defined routines and data structures
 70: */
 71: typedef struct {
 72:   PetscScalar u,v,omega,temp;
 73: } Field;

 75: extern int FormInitialGuess(SNES,Vec,void*);
 76: extern int FormFunctionLocal(DALocalInfo*,Field**,Field**,void*);
 77: extern int FormFunctionLocali(DALocalInfo*,MatStencil*,Field**,PetscScalar*,void*);

 79: typedef struct {
 80:    PassiveReal  lidvelocity,prandtl,grashof;  /* physical parameters */
 81:    PetscTruth     draw_contours;                /* flag - 1 indicates drawing contours */
 82: } AppCtx;

 84: int main(int argc,char **argv)
 85: {
 86:   DMMG       *dmmg;               /* multilevel grid structure */
 87:   AppCtx     user;                /* user-defined work context */
 88:   int        mx,my,its;
 89:   int        ierr;
 90:   MPI_Comm   comm;
 91:   SNES       snes;
 92:   DA         da;

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


 98:   PreLoadBegin(PETSC_TRUE,"SetUp");
 99:     DMMGCreate(comm,2,&user,&dmmg);


102:     /*
103:       Create distributed array multigrid object (DMMG) to manage parallel grid and vectors
104:       for principal unknowns (x) and governing residuals (f)
105:     */
106:     DACreate2d(comm,DA_NONPERIODIC,DA_STENCIL_STAR,-4,-4,PETSC_DECIDE,PETSC_DECIDE,4,1,0,0,&da);
107:     DMMGSetDM(dmmg,(DM)da);
108:     DADestroy(da);

110:     DAGetInfo(DMMGGetDA(dmmg),0,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
111:                      PETSC_IGNORE,PETSC_IGNORE);
112:     /* 
113:      Problem parameters (velocity of lid, prandtl, and grashof numbers)
114:     */
115:     user.lidvelocity = 1.0/(mx*my);
116:     user.prandtl     = 1.0;
117:     user.grashof     = 1.0;
118:     PetscOptionsGetReal(PETSC_NULL,"-lidvelocity",&user.lidvelocity,PETSC_NULL);
119:     PetscOptionsGetReal(PETSC_NULL,"-prandtl",&user.prandtl,PETSC_NULL);
120:     PetscOptionsGetReal(PETSC_NULL,"-grashof",&user.grashof,PETSC_NULL);
121:     PetscOptionsHasName(PETSC_NULL,"-contours",&user.draw_contours);

123:     DASetFieldName(DMMGGetDA(dmmg),0,"x-velocity");
124:     DASetFieldName(DMMGGetDA(dmmg),1,"y-velocity");
125:     DASetFieldName(DMMGGetDA(dmmg),2,"Omega");
126:     DASetFieldName(DMMGGetDA(dmmg),3,"temperature");

128:     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129:        Create user context, set problem data, create vector data structures.
130:        Also, compute the initial guess.
131:        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

133:     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134:        Create nonlinear solver context

136:        Process adiC: FormFunctionLocal FormFunctionLocali
137:        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138:     DMMGSetSNESLocal(dmmg,FormFunctionLocal,0,ad_FormFunctionLocal,admf_FormFunctionLocal);
139:     DMMGSetSNESLocali(dmmg,FormFunctionLocali,ad_FormFunctionLocali,admf_FormFunctionLocali);

141:     PetscPrintf(comm,"lid velocity = %g, prandtl # = %g, grashof # = %gn",
142:                        user.lidvelocity,user.prandtl,user.grashof);


145:     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146:        Solve the nonlinear system
147:        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
148:     DMMGSetInitialGuess(dmmg,FormInitialGuess);

150:   PreLoadStage("Solve");
151:     DMMGSolve(dmmg);

153:     snes = DMMGGetSNES(dmmg);
154:     SNESGetIterationNumber(snes,&its);
155:     PetscPrintf(comm,"Number of Newton iterations = %dn", its);

157:     /*
158:       Visualize solution
159:     */

161:     if (user.draw_contours) {
162:       VecView(DMMGGetx(dmmg),PETSC_VIEWER_DRAW_WORLD);
163:     }

165:     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
166:        Free work space.  All PETSc objects should be destroyed when they
167:        are no longer needed.
168:        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

170:     DMMGDestroy(dmmg);
171:   PreLoadEnd();

173:   PetscFinalize();
174:   return 0;
175: }

177: /* ------------------------------------------------------------------- */


180: /* 
181:    FormInitialGuess - Forms initial approximation.

183:    Input Parameters:
184:    user - user-defined application context
185:    X - vector

187:    Output Parameter:
188:    X - vector
189:  */
190: int FormInitialGuess(SNES snes,Vec X,void *ptr)
191: {
192:   DMMG      dmmg = (DMMG)ptr;
193:   AppCtx    *user = (AppCtx*)dmmg->user;
194:   DA        da = (DA)dmmg->dm;
195:   int       i,j,mx,ierr,xs,ys,xm,ym;
196:   PetscReal grashof,dx;
197:   Field     **x;

199:   grashof = user->grashof;

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

204:   /*
205:      Get local grid boundaries (for 2-dimensional DA):
206:        xs, ys   - starting grid indices (no ghost points)
207:        xm, ym   - widths of local grid (no ghost points)
208:   */
209:   DAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);

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

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

233:   /*
234:      Restore vector
235:   */
236:   DAVecRestoreArray(da,X,(void**)&x);
237:   return 0;
238: }
239: int FormFunctionLocal(DALocalInfo *info,Field **x,Field **f,void *ptr)
240:  {
241:   AppCtx       *user = (AppCtx*)ptr;
242:   int          ierr,i,j;
243:   int          xints,xinte,yints,yinte;
244:   PetscReal    hx,hy,dhx,dhy,hxdhy,hydhx;
245:   PetscReal    grashof,prandtl,lid;
246:   PetscScalar  u,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;

249:   grashof = user->grashof;
250:   prandtl = user->prandtl;
251:   lid     = user->lidvelocity;

253:   /* 
254:      Define mesh intervals ratios for uniform grid.
255:      [Note: FD formulae below are normalized by multiplying through by
256:      local volume element to obtain coefficients O(1) in two dimensions.]
257:   */
258:   dhx = (PetscReal)(info->mx-1);  dhy = (PetscReal)(info->my-1);
259:   hx = 1.0/dhx;                   hy = 1.0/dhy;
260:   hxdhy = hx*dhy;                 hydhx = hy*dhx;

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

264:   /* Test whether we are on the bottom edge of the global array */
265:   if (yints == 0) {
266:     j = 0;
267:     yints = yints + 1;
268:     /* bottom edge */
269:     for (i=info->xs; i<info->xs+info->xm; i++) {
270:         f[j][i].u     = x[j][i].u;
271:         f[j][i].v     = x[j][i].v;
272:         f[j][i].omega = x[j][i].omega + (x[j+1][i].u - x[j][i].u)*dhy;
273:         f[j][i].temp  = x[j][i].temp-x[j+1][i].temp;
274:     }
275:   }

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

290:   /* Test whether we are on the left edge of the global array */
291:   if (xints == 0) {
292:     i = 0;
293:     xints = xints + 1;
294:     /* left edge */
295:     for (j=info->ys; j<info->ys+info->ym; j++) {
296:       f[j][i].u     = x[j][i].u;
297:       f[j][i].v     = x[j][i].v;
298:       f[j][i].omega = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx;
299:       f[j][i].temp  = x[j][i].temp;
300:     }
301:   }

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

316:   /* Compute over the interior points */
317:   for (j=yints; j<yinte; j++) {
318:     for (i=xints; i<xinte; i++) {

320:         /*
321:           convective coefficients for upwinding
322:         */
323:         vx = x[j][i].u; avx = PetscAbsScalar(vx);
324:         vxp = .5*(vx+avx); vxm = .5*(vx-avx);
325:         vy = x[j][i].v; avy = PetscAbsScalar(vy);
326:         vyp = .5*(vy+avy); vym = .5*(vy-avy);

328:         /* U velocity */
329:         u          = x[j][i].u;
330:         uxx        = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
331:         uyy        = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
332:         f[j][i].u  = uxx + uyy - .5*(x[j+1][i].omega-x[j-1][i].omega)*hx;

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

340:         /* Omega */
341:         u          = x[j][i].omega;
342:         uxx        = (2.0*u - x[j][i-1].omega - x[j][i+1].omega)*hydhx;
343:         uyy        = (2.0*u - x[j-1][i].omega - x[j+1][i].omega)*hxdhy;
344:         f[j][i].omega = uxx + uyy +
345:                         (vxp*(u - x[j][i-1].omega) +
346:                           vxm*(x[j][i+1].omega - u)) * hy +
347:                         (vyp*(u - x[j-1][i].omega) +
348:                           vym*(x[j+1][i].omega - u)) * hx -
349:                         .5 * grashof * (x[j][i+1].temp - x[j][i-1].temp) * hy;

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

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

370: /*
371:     This is an experimental function and can be safely ignored.
372: */
373: int FormFunctionLocali(DALocalInfo *info,MatStencil *st,Field **x,PetscScalar *f,void *ptr)
374:  {
375:   AppCtx      *user = (AppCtx*)ptr;
376:   int         i,j,c;
377:   PassiveReal hx,hy,dhx,dhy,hxdhy,hydhx;
378:   PassiveReal grashof,prandtl,lid;
379:   PetscScalar u,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;

382:   grashof = user->grashof;
383:   prandtl = user->prandtl;
384:   lid     = user->lidvelocity;

386:   /* 
387:      Define mesh intervals ratios for uniform grid.
388:      [Note: FD formulae below are normalized by multiplying through by
389:      local volume element to obtain coefficients O(1) in two dimensions.]
390:   */
391:   dhx = (PetscReal)(info->mx-1);     dhy = (PetscReal)(info->my-1);
392:   hx = 1.0/dhx;                   hy = 1.0/dhy;
393:   hxdhy = hx*dhy;                 hydhx = hy*dhx;

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

397:   /* Test whether we are on the right edge of the global array */
398:   if (i == info->mx-1) {
399:     if (c == 0) *f     = x[j][i].u;
400:     else if (c == 1) *f     = x[j][i].v;
401:     else if (c == 2) *f = x[j][i].omega - (x[j][i].v - x[j][i-1].v)*dhx;
402:     else *f  = x[j][i].temp - (PetscReal)(grashof>0);

404:   /* Test whether we are on the left edge of the global array */
405:   } else if (i == 0) {
406:     if (c == 0) *f     = x[j][i].u;
407:     else if (c == 1) *f     = x[j][i].v;
408:     else if (c == 2) *f = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx;
409:     else *f  = x[j][i].temp;

411:   /* Test whether we are on the top edge of the global array */
412:   } else if (j == info->my-1) {
413:     if (c == 0) *f     = x[j][i].u - lid;
414:     else if (c == 1) *f     = x[j][i].v;
415:     else if (c == 2) *f = x[j][i].omega + (x[j][i].u - x[j-1][i].u)*dhy;
416:     else *f  = x[j][i].temp-x[j-1][i].temp;

418:   /* Test whether we are on the bottom edge of the global array */
419:   } else if (j == 0) {
420:     if (c == 0) *f     = x[j][i].u;
421:     else if (c == 1) *f     = x[j][i].v;
422:     else if (c == 2) *f = x[j][i].omega + (x[j+1][i].u - x[j][i].u)*dhy;
423:     else *f  = x[j][i].temp-x[j+1][i].temp;

425:   /* Compute over the interior points */
426:   } else {
427:     /*
428:       convective coefficients for upwinding
429:     */
430:     vx = x[j][i].u; avx = PetscAbsScalar(vx);
431:     vxp = .5*(vx+avx); vxm = .5*(vx-avx);
432:     vy = x[j][i].v; avy = PetscAbsScalar(vy);
433:     vyp = .5*(vy+avy); vym = .5*(vy-avy);

435:     /* U velocity */
436:     if (c == 0) {
437:       u          = x[j][i].u;
438:       uxx        = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
439:       uyy        = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
440:       *f         = uxx + uyy - .5*(x[j+1][i].omega-x[j-1][i].omega)*hx;

442:     /* V velocity */
443:     } else if (c == 1) {
444:       u          = x[j][i].v;
445:       uxx        = (2.0*u - x[j][i-1].v - x[j][i+1].v)*hydhx;
446:       uyy        = (2.0*u - x[j-1][i].v - x[j+1][i].v)*hxdhy;
447:       *f         = uxx + uyy + .5*(x[j][i+1].omega-x[j][i-1].omega)*hy;
448: 
449:     /* Omega */
450:     } else if (c == 2) {
451:       u          = x[j][i].omega;
452:       uxx        = (2.0*u - x[j][i-1].omega - x[j][i+1].omega)*hydhx;
453:       uyy        = (2.0*u - x[j-1][i].omega - x[j+1][i].omega)*hxdhy;
454:       *f         = uxx + uyy +
455:         (vxp*(u - x[j][i-1].omega) +
456:          vxm*(x[j][i+1].omega - u)) * hy +
457:         (vyp*(u - x[j-1][i].omega) +
458:          vym*(x[j+1][i].omega - u)) * hx -
459:         .5 * grashof * (x[j][i+1].temp - x[j][i-1].temp) * hy;
460: 
461:     /* Temperature */
462:     } else {
463:       u           = x[j][i].temp;
464:       uxx         = (2.0*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx;
465:       uyy         = (2.0*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy;
466:       *f          =  uxx + uyy  + prandtl * (
467:                                              (vxp*(u - x[j][i-1].temp) +
468:                                               vxm*(x[j][i+1].temp - u)) * hy +
469:                                              (vyp*(u - x[j-1][i].temp) +
470:                                               vym*(x[j+1][i].temp - u)) * hx);
471:     }
472:   }

474:   return(0);
475: }