Actual source code: ex49.c

  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:       The same as ex19.c except it uses DMMGSetSNES() instead of DMMGSetSNESLocal() 
 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*/
 20: /* ------------------------------------------------------------------------

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

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

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

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

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

 54:   ------------------------------------------------------------------------- */

 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: #include <petscsnes.h>
 67: #include <petscdmda.h>
 68: #include <petscdmmg.h>

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

 77: PetscErrorCode FormInitialGuess(DMMG,Vec);
 78: PetscErrorCode FormFunctionLocal(DMDALocalInfo*,Field**,Field**,void*);
 79: PetscErrorCode FormFunction(SNES,Vec,Vec,void*);

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

 88: int main(int argc,char **argv)
 89: {
 90:   DMMG           *dmmg;               /* multilevel grid structure */
 91:   AppCtx         user;                /* user-defined work context */
 92:   PetscInt       mx,my,its,nlevels=2;
 94:   MPI_Comm       comm;
 95:   SNES           snes;
 96:   DM             da;

 98:   PetscInitialize(&argc,&argv,(char *)0,help);if (ierr) return(1);
 99:   comm = PETSC_COMM_WORLD;

101:   PetscOptionsGetInt(PETSC_NULL,"-nlevels",&nlevels,PETSC_NULL);
102:     DMMGCreate(comm,nlevels,&user,&dmmg);

104:     /*
105:       Create distributed array multigrid object (DMMG) to manage parallel grid and vectors
106:       for principal unknowns (x) and governing residuals (f)
107:     */
108:     DMDACreate2d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_STENCIL_STAR,-4,-4,PETSC_DECIDE,PETSC_DECIDE,4,1,0,0,&da);
109:     DMMGSetDM(dmmg,(DM)da);
110:     DMDestroy(&da);

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

125:     DMDASetFieldName(DMMGGetDM(dmmg),0,"x-velocity");
126:     DMDASetFieldName(DMMGGetDM(dmmg),1,"y-velocity");
127:     DMDASetFieldName(DMMGGetDM(dmmg),2,"Omega");
128:     DMDASetFieldName(DMMGGetDM(dmmg),3,"temperature");

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

135:     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136:        Create nonlinear solver context

138:        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139:     DMMGSetSNES(dmmg,FormFunction,0);
140:     DMMGSetFromOptions(dmmg);

142:     PetscPrintf(comm,"lid velocity = %G, prandtl # = %G, grashof # = %G\n",
143:                        user.lidvelocity,user.prandtl,user.grashof);


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

151:     DMMGSolve(dmmg);

153:     snes = DMMGGetSNES(dmmg);
154:     SNESGetIterationNumber(snes,&its);
155:     PetscPrintf(comm,"Number of Newton iterations = %D\n", 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);

172: /********  PetscDraw draw;
173:   PetscViewerDrawGetDraw(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),0,&draw);
174:   PetscDrawSetPause(draw,-1);
175:   PetscDrawPause(draw); */
176:   PetscFinalize();
177:   return 0;
178: }

180: /* ------------------------------------------------------------------- */


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

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

192:    Output Parameter:
193:    X - vector
194:  */
195: PetscErrorCode FormInitialGuess(DMMG dmmg,Vec X)
196: {
197:   AppCtx         *user = (AppCtx*)dmmg->user;
198:   DM             da = dmmg->dm;
199:   PetscInt       i,j,mx,xs,ys,xm,ym;
201:   PetscReal      grashof,dx;
202:   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,PETSC_NULL,&xm,&ym,PETSC_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: }
244: 
245: PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,Field **x,Field **f,void *ptr)
246:  {
247:   AppCtx         *user = (AppCtx*)ptr;
249:   PetscInt       xints,xinte,yints,yinte,i,j;
250:   PetscReal      hx,hy,dhx,dhy,hxdhy,hydhx;
251:   PetscReal      grashof,prandtl,lid;
252:   PetscScalar    u,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

372:   /*
373:      Flop count (multiply-adds are counted as 2 operations)
374:   */
375:   PetscLogFlops(84.0*info->ym*info->xm);
376:   return(0);
377: }

381: PetscErrorCode FormFunction(SNES snes, Vec X, Vec F, void *ctx)
382: {
383:   DMDALocalInfo    info;
384:   void          *u;
385:   void          *fu;
387:   DMMG           dmmg = (DMMG)ctx;
388:   Vec            localX;
389:   DM             da = dmmg->dm;

392:   DMGetLocalVector(da,&localX);
393:   /*
394:      Scatter ghost points to local vector, using the 2-step process
395:         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
396:   */
397:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
398:   DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
399:   DMDAGetLocalInfo(da,&info);
400:   DMDAVecGetArray(da,localX,&u);
401:   DMDAVecGetArray(da,F,&fu);
402:   FormFunctionLocal(&info,u,fu,dmmg->user);
403:   DMDAVecRestoreArray(da,localX,&u);
404:   DMDAVecRestoreArray(da,F,&fu);
405:   DMRestoreLocalVector(da,&localX);
406:   return(0);
407: }