Actual source code: ex26.c

petsc-3.3-p7 2013-05-11
  2: static char help[] = "Transient nonlinear driven cavity 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/snes/examples/tutorials/ex50.c for the steady-state version.

 13:       See src/snes/examples/tutorials/ex27.c for a version written specifically
 14:       for pseudo-transient continuation, without using TS.
 15: */

 17: /*T
 18:    Concepts: TS^solving a system of nonlinear equations (parallel multicomponent example);
 19:    Concepts: DMDA^using distributed arrays;
 20:    Concepts: multicomponent;differential-algebraic equation
 21:    Processors: n
 22: T*/
 23: /* ------------------------------------------------------------------------

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

 28:     This problem is modeled by the partial differential equation system

 30:         - Lap(U) - Grad_y(Omega) = 0
 31:         - Lap(V) + Grad_x(Omega) = 0
 32:         Omega_t - Lap(Omega) + Div([U*Omega,V*Omega]) - GR*Grad_x(T) = 0
 33:         T_t - Lap(T) + PR*Div([U*T,V*T]) = 0

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

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

 46:     A finite difference approximation with the usual 5-point stencil
 47:     is used to discretize the boundary value problem to obtain a
 48:     nonlinear system of equations.  Upwinding is used for the divergence
 49:     (convective) terms and central for the gradient (source) terms.

 51:     The Jacobian can be either
 52:       * formed via finite differencing using coloring (the default), or
 53:       * applied matrix-free via the option -snes_mf
 54:         (for larger grid problems this variant may not converge
 55:         without a preconditioner due to ill-conditioning).

 57:   ------------------------------------------------------------------------- */

 59: /*
 60:    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
 61:    Include "petscts.h" so that we can use TS solvers.  Note that this
 62:    file automatically includes:
 63:      petscsys.h       - base PETSc routines   petscvec.h - vectors
 64:      petscmat.h - matrices
 65:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 66:      petscviewer.h - viewers               petscpc.h  - preconditioners
 67:      petscksp.h   - linear solvers         petscsnes.h - nonlinear solvers
 68: */
 69: #include <petscts.h>
 70: #include <petscdmda.h>

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

 79: PetscErrorCode FormIFunctionLocal(DMDALocalInfo*,PetscReal,Field**,Field**,Field**,void*);
 80: PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);

 82: typedef struct {
 83:   PassiveReal  lidvelocity,prandtl,grashof;  /* physical parameters */
 84:   PetscBool    parabolic;                    /* allow a transient term corresponding roughly to artificial compressibility */
 85:   PetscBool    draw_contours;                /* flag - 1 indicates drawing contours */
 86:   PetscReal    cfl_initial;                  /* CFL for first time step */
 87: } AppCtx;

 89: PetscErrorCode FormInitialSolution(TS,Vec,AppCtx*);

 93: int main(int argc,char **argv)
 94: {
 95:   AppCtx         user;                /* user-defined work context */
 96:   PetscInt       mx,my,steps;
 98:   TS             ts;
 99:   DM             da;
100:   Vec            X;
101:   PetscReal      ftime;
102:   PetscBool      fdcoloring;
103:   TSConvergedReason reason;
104:   MatFDColoring matfdcoloring = PETSC_NULL;

106:   PetscInitialize(&argc,&argv,(char *)0,help);if (ierr) return(1);

108:   TSCreate(PETSC_COMM_WORLD,&ts);

110:   DMDACreate2d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_STENCIL_STAR,-4,-4,PETSC_DECIDE,PETSC_DECIDE,4,1,0,0,&da);
111:   TSSetDM(ts,(DM)da);

113:   DMDAGetInfo(da,0,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
114:                    PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
115:   /*
116:      Problem parameters (velocity of lid, prandtl, and grashof numbers)
117:   */
118:   user.lidvelocity   = 1.0/(mx*my);
119:   user.prandtl       = 1.0;
120:   user.grashof       = 1.0;
121:   user.draw_contours = PETSC_FALSE;
122:   user.parabolic     = PETSC_FALSE;
123:   user.cfl_initial   = 50.;
124:   fdcoloring         = PETSC_TRUE;
125:   PetscOptionsBegin(PETSC_COMM_WORLD,PETSC_NULL,"Driven cavity/natural convection options","");
126:   PetscOptionsReal("-lidvelocity","Lid velocity, related to Reynolds number","",user.lidvelocity,&user.lidvelocity,PETSC_NULL);
127:   PetscOptionsReal("-prandtl","Ratio of viscous to thermal diffusivity","",user.prandtl,&user.prandtl,PETSC_NULL);
128:   PetscOptionsReal("-grashof","Ratio of bouyant to viscous forces","",user.grashof,&user.grashof,PETSC_NULL);
129:   PetscOptionsBool("-draw_contours","Plot the solution with contours","",user.draw_contours,&user.draw_contours,PETSC_NULL);
130:   PetscOptionsBool("-parabolic","Relax incompressibility to make the system parabolic instead of differential-algebraic","",user.parabolic,&user.parabolic,PETSC_NULL);
131:   PetscOptionsReal("-cfl_initial","Advective CFL for the first time step","",user.cfl_initial,&user.cfl_initial,PETSC_NULL);
132:   PetscOptionsBool("-fdcoloring","Use finite difference coloring to compute Jacobian","",fdcoloring,&fdcoloring,PETSC_NULL);
133:   PetscOptionsEnd();

135:   DMDASetFieldName(da,0,"x-velocity");
136:   DMDASetFieldName(da,1,"y-velocity");
137:   DMDASetFieldName(da,2,"Omega");
138:   DMDASetFieldName(da,3,"temperature");

140:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141:      Create user context, set problem data, create vector data structures.
142:      Also, compute the initial guess.
143:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

145:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146:      Create time integration context
147:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
148:   DMSetApplicationContext(da,&user);
149:   TSSetIFunction(ts,PETSC_NULL,FormIFunction,&user);
150:   if (fdcoloring) {
151:     Mat B;
152:     SNES snes;
153:     ISColoring iscoloring;
154:     TSGetSNES(ts,&snes);
155:     DMCreateMatrix(da,MATAIJ,&B);
156:     DMCreateColoring(da,IS_COLORING_GLOBAL,MATAIJ,&iscoloring);
157:     MatFDColoringCreate(B,iscoloring,&matfdcoloring);
158:     ISColoringDestroy(&iscoloring);
159:     MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode(*)(void))SNESTSFormFunction,ts);
160:     MatFDColoringSetFromOptions(matfdcoloring);
161:     SNESSetJacobian(snes,B,B,SNESDefaultComputeJacobianColor,matfdcoloring);
162:     MatDestroy(&B);
163:   }
164:   TSSetDuration(ts,10000,1e12);
165:   TSSetInitialTimeStep(ts,0.0,user.cfl_initial/(user.lidvelocity*mx));
166:   TSSetFromOptions(ts);

168:   PetscPrintf(PETSC_COMM_WORLD,"%Dx%D grid, lid velocity = %G, prandtl # = %G, grashof # = %G\n",mx,my,user.lidvelocity,user.prandtl,user.grashof);


171:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
172:      Solve the nonlinear system
173:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

175:   DMCreateGlobalVector(da,&X);
176:   FormInitialSolution(ts,X,&user);

178:   TSSolve(ts,X,&ftime);
179:   TSGetTimeStepNumber(ts,&steps);
180:   TSGetConvergedReason(ts,&reason);

182:   PetscPrintf(PETSC_COMM_WORLD,"%s at time %G after %D steps\n",TSConvergedReasons[reason],ftime,steps);

184:   /*
185:      Visualize solution
186:   */
187:   if (user.draw_contours) {
188:     VecView(X,PETSC_VIEWER_DRAW_WORLD);
189:   }

191:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
192:      Free work space.  All PETSc objects should be destroyed when they
193:      are no longer needed.
194:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
195:   VecDestroy(&X);
196:   MatFDColoringDestroy(&matfdcoloring);
197:   DMDestroy(&da);
198:   TSDestroy(&ts);

200:   PetscFinalize();
201:   return 0;
202: }

204: /* ------------------------------------------------------------------- */


209: /*
210:    FormInitialSolution - Forms initial approximation.

212:    Input Parameters:
213:    user - user-defined application context
214:    X - vector

216:    Output Parameter:
217:    X - vector
218:  */
219: PetscErrorCode FormInitialSolution(TS ts,Vec X,AppCtx *user)
220: {
221:   DM             da;
222:   PetscInt       i,j,mx,xs,ys,xm,ym;
224:   PetscReal      grashof,dx;
225:   Field          **x;

227:   grashof = user->grashof;
228:   TSGetDM(ts,&da);
229:   DMDAGetInfo(da,0,&mx,0,0,0,0,0,0,0,0,0,0,0);
230:   dx  = 1.0/(mx-1);

232:   /*
233:      Get local grid boundaries (for 2-dimensional DMDA):
234:        xs, ys   - starting grid indices (no ghost points)
235:        xm, ym   - widths of local grid (no ghost points)
236:   */
237:   DMDAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);

239:   /*
240:      Get a pointer to vector data.
241:        - For default PETSc vectors, VecGetArray() returns a pointer to
242:          the data array.  Otherwise, the routine is implementation dependent.
243:        - You MUST call VecRestoreArray() when you no longer need access to
244:          the array.
245:   */
246:   DMDAVecGetArray(da,X,&x);

248:   /*
249:      Compute initial guess over the locally owned part of the grid
250:      Initial condition is motionless fluid and equilibrium temperature
251:   */
252:   for (j=ys; j<ys+ym; j++) {
253:     for (i=xs; i<xs+xm; i++) {
254:       x[j][i].u     = 0.0;
255:       x[j][i].v     = 0.0;
256:       x[j][i].omega = 0.0;
257:       x[j][i].temp  = (grashof>0)*i*dx;
258:     }
259:   }

261:   /*
262:      Restore vector
263:   */
264:   DMDAVecRestoreArray(da,X,&x);
265:   return 0;
266: }

270: PetscErrorCode FormIFunctionLocal(DMDALocalInfo *info,PetscReal ptime,Field **x,Field **xdot,Field **f,void *ptr)
271:  {
272:   AppCtx         *user = (AppCtx*)ptr;
274:   PetscInt       xints,xinte,yints,yinte,i,j;
275:   PetscReal      hx,hy,dhx,dhy,hxdhy,hydhx;
276:   PetscReal      grashof,prandtl,lid;
277:   PetscScalar    u,udot,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;

280:   grashof = user->grashof;
281:   prandtl = user->prandtl;
282:   lid     = user->lidvelocity;

284:   /*
285:      Define mesh intervals ratios for uniform grid.

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


291:   */
292:   dhx = (PetscReal)(info->mx-1);  dhy = (PetscReal)(info->my-1);
293:   hx = 1.0/dhx;                   hy = 1.0/dhy;
294:   hxdhy = hx*dhy;                 hydhx = hy*dhx;

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

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

311:   /* Test whether we are on the top edge of the global array */
312:   if (yinte == info->my) {
313:     j = info->my - 1;
314:     yinte = yinte - 1;
315:     /* top edge */
316:     for (i=info->xs; i<info->xs+info->xm; i++) {
317:         f[j][i].u     = x[j][i].u - lid;
318:         f[j][i].v     = x[j][i].v;
319:         f[j][i].omega = x[j][i].omega + (x[j][i].u - x[j-1][i].u)*dhy;
320:         f[j][i].temp  = x[j][i].temp-x[j-1][i].temp;
321:     }
322:   }

324:   /* Test whether we are on the left edge of the global array */
325:   if (xints == 0) {
326:     i = 0;
327:     xints = xints + 1;
328:     /* left edge */
329:     for (j=info->ys; j<info->ys+info->ym; j++) {
330:       f[j][i].u     = x[j][i].u;
331:       f[j][i].v     = x[j][i].v;
332:       f[j][i].omega = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx;
333:       f[j][i].temp  = x[j][i].temp;
334:     }
335:   }

337:   /* Test whether we are on the right edge of the global array */
338:   if (xinte == info->mx) {
339:     i = info->mx - 1;
340:     xinte = xinte - 1;
341:     /* right edge */
342:     for (j=info->ys; j<info->ys+info->ym; j++) {
343:       f[j][i].u     = x[j][i].u;
344:       f[j][i].v     = x[j][i].v;
345:       f[j][i].omega = x[j][i].omega - (x[j][i].v - x[j][i-1].v)*dhx;
346:       f[j][i].temp  = x[j][i].temp - (PetscReal)(grashof>0);
347:     }
348:   }

350:   /* Compute over the interior points */
351:   for (j=yints; j<yinte; j++) {
352:     for (i=xints; i<xinte; i++) {

354:         /*
355:           convective coefficients for upwinding
356:         */
357:         vx = x[j][i].u; avx = PetscAbsScalar(vx);
358:         vxp = .5*(vx+avx); vxm = .5*(vx-avx);
359:         vy = x[j][i].v; avy = PetscAbsScalar(vy);
360:         vyp = .5*(vy+avy); vym = .5*(vy-avy);

362:         /* U velocity */
363:         u          = x[j][i].u;
364:         udot       = user->parabolic ? xdot[j][i].u : 0.;
365:         uxx        = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
366:         uyy        = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
367:         f[j][i].u  = udot + uxx + uyy - .5*(x[j+1][i].omega-x[j-1][i].omega)*hx;

369:         /* V velocity */
370:         u          = x[j][i].v;
371:         udot       = user->parabolic ? xdot[j][i].v : 0.;
372:         uxx        = (2.0*u - x[j][i-1].v - x[j][i+1].v)*hydhx;
373:         uyy        = (2.0*u - x[j-1][i].v - x[j+1][i].v)*hxdhy;
374:         f[j][i].v  = udot + uxx + uyy + .5*(x[j][i+1].omega-x[j][i-1].omega)*hy;

376:         /* Omega */
377:         u          = x[j][i].omega;
378:         uxx        = (2.0*u - x[j][i-1].omega - x[j][i+1].omega)*hydhx;
379:         uyy        = (2.0*u - x[j-1][i].omega - x[j+1][i].omega)*hxdhy;
380:         f[j][i].omega = (xdot[j][i].omega + uxx + uyy
381:                          + (vxp*(u - x[j][i-1].omega)
382:                             + vxm*(x[j][i+1].omega - u)) * hy
383:                          + (vyp*(u - x[j-1][i].omega)
384:                             + vym*(x[j+1][i].omega - u)) * hx
385:                          - .5 * grashof * (x[j][i+1].temp - x[j][i-1].temp) * hy);

387:         /* Temperature */
388:         u             = x[j][i].temp;
389:         uxx           = (2.0*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx;
390:         uyy           = (2.0*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy;
391:         f[j][i].temp =  (xdot[j][i].temp + uxx + uyy
392:                          + prandtl * ((vxp*(u - x[j][i-1].temp)
393:                                        + vxm*(x[j][i+1].temp - u)) * hy
394:                                       + (vyp*(u - x[j-1][i].temp)
395:                                          + vym*(x[j+1][i].temp - u)) * hx));
396:     }
397:   }

399:   /*
400:      Flop count (multiply-adds are counted as 2 operations)
401:   */
402:   PetscLogFlops(84.0*info->ym*info->xm);
403:   return(0);
404: }

408: PetscErrorCode FormIFunction(TS ts,PetscReal ptime,Vec X,Vec Xdot,Vec F,void *user)
409: {
410:   DMDALocalInfo  info;
411:   Field          **u,**udot,**fu;
413:   Vec            localX;
414:   DM             da;

417:   TSGetDM(ts,(DM*)&da);
418:   DMGetLocalVector(da,&localX);
419:   /*
420:      Scatter ghost points to local vector, using the 2-step process
421:         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
422:   */
423:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
424:   DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
425:   DMDAGetLocalInfo(da,&info);
426:   DMDAVecGetArray(da,localX,&u);
427:   DMDAVecGetArray(da,Xdot,&udot);
428:   DMDAVecGetArray(da,F,&fu);
429:   FormIFunctionLocal(&info,ptime,u,udot,fu,user);
430:   DMDAVecRestoreArray(da,localX,&u);
431:   DMDAVecRestoreArray(da,Xdot,&udot);
432:   DMDAVecRestoreArray(da,F,&fu);
433:   DMRestoreLocalVector(da,&localX);
434:   return(0);
435: }