Actual source code: ex26.c

petsc-master 2017-04-29
Report Typos and Errors

  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/ex19.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: TS^multicomponent
 21:    Concepts: TS^differential-algebraic equation
 22:    Processors: n
 23: T*/
 24: /* ------------------------------------------------------------------------

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

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

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

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

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

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

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

 58:   ------------------------------------------------------------------------- */

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

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

 81: PetscErrorCode FormIFunctionLocal(DMDALocalInfo*,PetscReal,Field**,Field**,Field**,void*);

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

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

 91: int main(int argc,char **argv)
 92: {
 93:   AppCtx            user;             /* user-defined work context */
 94:   PetscInt          mx,my,steps;
 95:   PetscErrorCode    ierr;
 96:   TS                ts;
 97:   DM                da;
 98:   Vec               X;
 99:   PetscReal         ftime;
100:   TSConvergedReason reason;

102:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
103:   TSCreate(PETSC_COMM_WORLD,&ts);
104:   DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,4,4,PETSC_DECIDE,PETSC_DECIDE,4,1,0,0,&da);
105:   DMSetFromOptions(da);
106:   DMSetUp(da);
107:   TSSetDM(ts,(DM)da);

109:   DMDAGetInfo(da,0,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
110:                      PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
111:   /*
112:      Problem parameters (velocity of lid, prandtl, and grashof numbers)
113:   */
114:   user.lidvelocity = 1.0/(mx*my);
115:   user.prandtl     = 1.0;
116:   user.grashof     = 1.0;
117:   user.parabolic   = PETSC_FALSE;
118:   user.cfl_initial = 50.;

120:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Driven cavity/natural convection options","");
121:   PetscOptionsReal("-lidvelocity","Lid velocity, related to Reynolds number","",user.lidvelocity,&user.lidvelocity,NULL);
122:   PetscOptionsReal("-prandtl","Ratio of viscous to thermal diffusivity","",user.prandtl,&user.prandtl,NULL);
123:   PetscOptionsReal("-grashof","Ratio of bouyant to viscous forces","",user.grashof,&user.grashof,NULL);
124:   PetscOptionsBool("-parabolic","Relax incompressibility to make the system parabolic instead of differential-algebraic","",user.parabolic,&user.parabolic,NULL);
125:   PetscOptionsReal("-cfl_initial","Advective CFL for the first time step","",user.cfl_initial,&user.cfl_initial,NULL);
126:   PetscOptionsEnd();

128:   DMDASetFieldName(da,0,"x-velocity");
129:   DMDASetFieldName(da,1,"y-velocity");
130:   DMDASetFieldName(da,2,"Omega");
131:   DMDASetFieldName(da,3,"temperature");

133:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134:      Create user context, set problem data, create vector data structures.
135:      Also, compute the initial guess.
136:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

138:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
139:      Create time integration context
140:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
141:   DMSetApplicationContext(da,&user);
142:   DMDATSSetIFunctionLocal(da,INSERT_VALUES,(DMDATSIFunctionLocal)FormIFunctionLocal,&user);
143:   TSSetDuration(ts,10000,1e12);
144:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
145:   TSSetInitialTimeStep(ts,0.0,user.cfl_initial/(user.lidvelocity*mx));
146:   TSSetFromOptions(ts);

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


151:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152:      Solve the nonlinear system
153:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

155:   DMCreateGlobalVector(da,&X);
156:   FormInitialSolution(ts,X,&user);

158:   TSSolve(ts,X);
159:   TSGetSolveTime(ts,&ftime);
160:   TSGetTimeStepNumber(ts,&steps);
161:   TSGetConvergedReason(ts,&reason);

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

165:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
166:      Free work space.  All PETSc objects should be destroyed when they
167:      are no longer needed.
168:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
169:   VecDestroy(&X);
170:   DMDestroy(&da);
171:   TSDestroy(&ts);

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

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


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

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

187:    Output Parameter:
188:    X - vector
189:  */
190: PetscErrorCode FormInitialSolution(TS ts,Vec X,AppCtx *user)
191: {
192:   DM             da;
193:   PetscInt       i,j,mx,xs,ys,xm,ym;
195:   PetscReal      grashof,dx;
196:   Field          **x;

198:   grashof = user->grashof;
199:   TSGetDM(ts,&da);
200:   DMDAGetInfo(da,0,&mx,0,0,0,0,0,0,0,0,0,0,0);
201:   dx      = 1.0/(mx-1);

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

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

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

232:   /*
233:      Restore vector
234:   */
235:   DMDAVecRestoreArray(da,X,&x);
236:   return 0;
237: }

239: PetscErrorCode FormIFunctionLocal(DMDALocalInfo *info,PetscReal ptime,Field **x,Field **xdot,Field **f,void *ptr)
240: {
241:   AppCtx         *user = (AppCtx*)ptr;
243:   PetscInt       xints,xinte,yints,yinte,i,j;
244:   PetscReal      hx,hy,dhx,dhy,hxdhy,hydhx;
245:   PetscReal      grashof,prandtl,lid;
246:   PetscScalar    u,udot,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.

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


260:   */
261:   dhx   = (PetscReal)(info->mx-1);  dhy = (PetscReal)(info->my-1);
262:   hx    = 1.0/dhx;                   hy = 1.0/dhy;
263:   hxdhy = hx*dhy;                 hydhx = hy*dhx;

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

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

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

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

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

319:   /* Compute over the interior points */
320:   for (j=yints; j<yinte; j++) {
321:     for (i=xints; i<xinte; i++) {

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

331:       /* U velocity */
332:       u         = x[j][i].u;
333:       udot      = user->parabolic ? xdot[j][i].u : 0.;
334:       uxx       = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
335:       uyy       = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
336:       f[j][i].u = udot + uxx + uyy - .5*(x[j+1][i].omega-x[j-1][i].omega)*hx;

338:       /* V velocity */
339:       u         = x[j][i].v;
340:       udot      = user->parabolic ? xdot[j][i].v : 0.;
341:       uxx       = (2.0*u - x[j][i-1].v - x[j][i+1].v)*hydhx;
342:       uyy       = (2.0*u - x[j-1][i].v - x[j+1][i].v)*hxdhy;
343:       f[j][i].v = udot + uxx + uyy + .5*(x[j][i+1].omega-x[j][i-1].omega)*hy;

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

356:       /* Temperature */
357:       u            = x[j][i].temp;
358:       uxx          = (2.0*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx;
359:       uyy          = (2.0*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy;
360:       f[j][i].temp =  (xdot[j][i].temp + uxx + uyy
361:                        + prandtl * ((vxp*(u - x[j][i-1].temp)
362:                                      + vxm*(x[j][i+1].temp - u)) * hy
363:                                     + (vyp*(u - x[j-1][i].temp)
364:                                        + 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: }