1 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/tutorials/ex19.c for the steady-state version. 12 13 There used to be a SNES example (src/snes/tutorials/ex27.c) that 14 implemented this algorithm without using TS and was used for the numerical 15 results in the paper 16 17 Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient 18 Continuation and Differential-Algebraic Equations, 2003. 19 20 That example was removed because it used obsolete interfaces, but the 21 algorithms from the paper can be reproduced using this example. 22 23 Note: The paper describes the algorithm as being linearly implicit but the 24 numerical results were created using nonlinearly implicit Euler. The 25 algorithm as described (linearly implicit) is more efficient and is the 26 default when using TSPSEUDO. If you want to reproduce the numerical 27 results from the paper, you'll have to change the SNES to converge the 28 nonlinear solve (e.g., -snes_type newtonls). The DAE versus ODE variants 29 are controlled using the -parabolic option. 30 31 Comment preserved from snes/tutorials/ex27.c, since removed: 32 33 [H]owever Figure 3.1 was generated with a slightly different algorithm 34 (see targets runex27 and runex27_p) than described in the paper. In 35 particular, the described algorithm is linearly implicit, advancing to 36 the next step after one Newton step, so that the steady state residual 37 is always used, but the figure was generated by converging each step to 38 a relative tolerance of 1.e-3. On the example problem, setting 39 -snes_type ksponly has only minor impact on number of steps, but 40 significantly reduces the required number of linear solves. 41 42 See also https://lists.mcs.anl.gov/pipermail/petsc-dev/2010-March/002362.html 43 */ 44 45 /* ------------------------------------------------------------------------ 46 47 We thank David E. Keyes for contributing the driven cavity discretization 48 within this example code. 49 50 This problem is modeled by the partial differential equation system 51 52 - Lap(U) - Grad_y(Omega) = 0 53 - Lap(V) + Grad_x(Omega) = 0 54 Omega_t - Lap(Omega) + Div([U*Omega,V*Omega]) - GR*Grad_x(T) = 0 55 T_t - Lap(T) + PR*Div([U*T,V*T]) = 0 56 57 in the unit square, which is uniformly discretized in each of x and 58 y in this simple encoding. 59 60 No-slip, rigid-wall Dirichlet conditions are used for [U,V]. 61 Dirichlet conditions are used for Omega, based on the definition of 62 vorticity: Omega = - Grad_y(U) + Grad_x(V), where along each 63 constant coordinate boundary, the tangential derivative is zero. 64 Dirichlet conditions are used for T on the left and right walls, 65 and insulation homogeneous Neumann conditions are used for T on 66 the top and bottom walls. 67 68 A finite difference approximation with the usual 5-point stencil 69 is used to discretize the boundary value problem to obtain a 70 nonlinear system of equations. Upwinding is used for the divergence 71 (convective) terms and central for the gradient (source) terms. 72 73 The Jacobian can be either 74 * formed via finite differencing using coloring (the default), or 75 * applied matrix-free via the option -snes_mf 76 (for larger grid problems this variant may not converge 77 without a preconditioner due to ill-conditioning). 78 79 ------------------------------------------------------------------------- */ 80 81 /* 82 Include "petscdmda.h" so that we can use distributed arrays (DMDAs). 83 Include "petscts.h" so that we can use TS solvers. Note that this 84 file automatically includes: 85 petscsys.h - base PETSc routines petscvec.h - vectors 86 petscmat.h - matrices 87 petscis.h - index sets petscksp.h - Krylov subspace methods 88 petscviewer.h - viewers petscpc.h - preconditioners 89 petscksp.h - linear solvers petscsnes.h - nonlinear solvers 90 */ 91 #include <petscts.h> 92 #include <petscdm.h> 93 #include <petscdmda.h> 94 95 /* 96 User-defined routines and data structures 97 */ 98 typedef struct { 99 PetscScalar u,v,omega,temp; 100 } Field; 101 102 PetscErrorCode FormIFunctionLocal(DMDALocalInfo*,PetscReal,Field**,Field**,Field**,void*); 103 104 typedef struct { 105 PetscReal lidvelocity,prandtl,grashof; /* physical parameters */ 106 PetscBool parabolic; /* allow a transient term corresponding roughly to artificial compressibility */ 107 PetscReal cfl_initial; /* CFL for first time step */ 108 } AppCtx; 109 110 PetscErrorCode FormInitialSolution(TS,Vec,AppCtx*); 111 112 int main(int argc,char **argv) 113 { 114 AppCtx user; /* user-defined work context */ 115 PetscInt mx,my,steps; 116 TS ts; 117 DM da; 118 Vec X; 119 PetscReal ftime; 120 TSConvergedReason reason; 121 122 PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 123 PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 124 PetscCall(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,4,4,PETSC_DECIDE,PETSC_DECIDE,4,1,0,0,&da)); 125 PetscCall(DMSetFromOptions(da)); 126 PetscCall(DMSetUp(da)); 127 PetscCall(TSSetDM(ts,(DM)da)); 128 129 PetscCall(DMDAGetInfo(da,0,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE)); 130 /* 131 Problem parameters (velocity of lid, prandtl, and grashof numbers) 132 */ 133 user.lidvelocity = 1.0/(mx*my); 134 user.prandtl = 1.0; 135 user.grashof = 1.0; 136 user.parabolic = PETSC_FALSE; 137 user.cfl_initial = 50.; 138 139 PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Driven cavity/natural convection options",""); 140 PetscCall(PetscOptionsReal("-lidvelocity","Lid velocity, related to Reynolds number","",user.lidvelocity,&user.lidvelocity,NULL)); 141 PetscCall(PetscOptionsReal("-prandtl","Ratio of viscous to thermal diffusivity","",user.prandtl,&user.prandtl,NULL)); 142 PetscCall(PetscOptionsReal("-grashof","Ratio of bouyant to viscous forces","",user.grashof,&user.grashof,NULL)); 143 PetscCall(PetscOptionsBool("-parabolic","Relax incompressibility to make the system parabolic instead of differential-algebraic","",user.parabolic,&user.parabolic,NULL)); 144 PetscCall(PetscOptionsReal("-cfl_initial","Advective CFL for the first time step","",user.cfl_initial,&user.cfl_initial,NULL)); 145 PetscOptionsEnd(); 146 147 PetscCall(DMDASetFieldName(da,0,"x-velocity")); 148 PetscCall(DMDASetFieldName(da,1,"y-velocity")); 149 PetscCall(DMDASetFieldName(da,2,"Omega")); 150 PetscCall(DMDASetFieldName(da,3,"temperature")); 151 152 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 153 Create user context, set problem data, create vector data structures. 154 Also, compute the initial guess. 155 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 156 157 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 158 Create time integration context 159 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 160 PetscCall(DMSetApplicationContext(da,&user)); 161 PetscCall(DMDATSSetIFunctionLocal(da,INSERT_VALUES,(DMDATSIFunctionLocal)FormIFunctionLocal,&user)); 162 PetscCall(TSSetMaxSteps(ts,10000)); 163 PetscCall(TSSetMaxTime(ts,1e12)); 164 PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 165 PetscCall(TSSetTimeStep(ts,user.cfl_initial/(user.lidvelocity*mx))); 166 PetscCall(TSSetFromOptions(ts)); 167 168 PetscCall(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)); 169 170 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 171 Solve the nonlinear system 172 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 173 174 PetscCall(DMCreateGlobalVector(da,&X)); 175 PetscCall(FormInitialSolution(ts,X,&user)); 176 177 PetscCall(TSSolve(ts,X)); 178 PetscCall(TSGetSolveTime(ts,&ftime)); 179 PetscCall(TSGetStepNumber(ts,&steps)); 180 PetscCall(TSGetConvergedReason(ts,&reason)); 181 182 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps)); 183 184 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 185 Free work space. All PETSc objects should be destroyed when they 186 are no longer needed. 187 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 188 PetscCall(VecDestroy(&X)); 189 PetscCall(DMDestroy(&da)); 190 PetscCall(TSDestroy(&ts)); 191 192 PetscCall(PetscFinalize()); 193 return 0; 194 } 195 196 /* ------------------------------------------------------------------- */ 197 198 /* 199 FormInitialSolution - Forms initial approximation. 200 201 Input Parameters: 202 user - user-defined application context 203 X - vector 204 205 Output Parameter: 206 X - vector 207 */ 208 PetscErrorCode FormInitialSolution(TS ts,Vec X,AppCtx *user) 209 { 210 DM da; 211 PetscInt i,j,mx,xs,ys,xm,ym; 212 PetscReal grashof,dx; 213 Field **x; 214 215 grashof = user->grashof; 216 PetscCall(TSGetDM(ts,&da)); 217 PetscCall(DMDAGetInfo(da,0,&mx,0,0,0,0,0,0,0,0,0,0,0)); 218 dx = 1.0/(mx-1); 219 220 /* 221 Get local grid boundaries (for 2-dimensional DMDA): 222 xs, ys - starting grid indices (no ghost points) 223 xm, ym - widths of local grid (no ghost points) 224 */ 225 PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); 226 227 /* 228 Get a pointer to vector data. 229 - For default PETSc vectors, VecGetArray() returns a pointer to 230 the data array. Otherwise, the routine is implementation dependent. 231 - You MUST call VecRestoreArray() when you no longer need access to 232 the array. 233 */ 234 PetscCall(DMDAVecGetArray(da,X,&x)); 235 236 /* 237 Compute initial guess over the locally owned part of the grid 238 Initial condition is motionless fluid and equilibrium temperature 239 */ 240 for (j=ys; j<ys+ym; j++) { 241 for (i=xs; i<xs+xm; i++) { 242 x[j][i].u = 0.0; 243 x[j][i].v = 0.0; 244 x[j][i].omega = 0.0; 245 x[j][i].temp = (grashof>0)*i*dx; 246 } 247 } 248 249 /* 250 Restore vector 251 */ 252 PetscCall(DMDAVecRestoreArray(da,X,&x)); 253 return 0; 254 } 255 256 PetscErrorCode FormIFunctionLocal(DMDALocalInfo *info,PetscReal ptime,Field **x,Field **xdot,Field **f,void *ptr) 257 { 258 AppCtx *user = (AppCtx*)ptr; 259 PetscInt xints,xinte,yints,yinte,i,j; 260 PetscReal hx,hy,dhx,dhy,hxdhy,hydhx; 261 PetscReal grashof,prandtl,lid; 262 PetscScalar u,udot,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym; 263 264 PetscFunctionBeginUser; 265 grashof = user->grashof; 266 prandtl = user->prandtl; 267 lid = user->lidvelocity; 268 269 /* 270 Define mesh intervals ratios for uniform grid. 271 272 Note: FD formulae below are normalized by multiplying through by 273 local volume element (i.e. hx*hy) to obtain coefficients O(1) in two dimensions. 274 275 */ 276 dhx = (PetscReal)(info->mx-1); dhy = (PetscReal)(info->my-1); 277 hx = 1.0/dhx; hy = 1.0/dhy; 278 hxdhy = hx*dhy; hydhx = hy*dhx; 279 280 xints = info->xs; xinte = info->xs+info->xm; yints = info->ys; yinte = info->ys+info->ym; 281 282 /* Test whether we are on the bottom edge of the global array */ 283 if (yints == 0) { 284 j = 0; 285 yints = yints + 1; 286 /* bottom edge */ 287 for (i=info->xs; i<info->xs+info->xm; i++) { 288 f[j][i].u = x[j][i].u; 289 f[j][i].v = x[j][i].v; 290 f[j][i].omega = x[j][i].omega + (x[j+1][i].u - x[j][i].u)*dhy; 291 f[j][i].temp = x[j][i].temp-x[j+1][i].temp; 292 } 293 } 294 295 /* Test whether we are on the top edge of the global array */ 296 if (yinte == info->my) { 297 j = info->my - 1; 298 yinte = yinte - 1; 299 /* top edge */ 300 for (i=info->xs; i<info->xs+info->xm; i++) { 301 f[j][i].u = x[j][i].u - lid; 302 f[j][i].v = x[j][i].v; 303 f[j][i].omega = x[j][i].omega + (x[j][i].u - x[j-1][i].u)*dhy; 304 f[j][i].temp = x[j][i].temp-x[j-1][i].temp; 305 } 306 } 307 308 /* Test whether we are on the left edge of the global array */ 309 if (xints == 0) { 310 i = 0; 311 xints = xints + 1; 312 /* left edge */ 313 for (j=info->ys; j<info->ys+info->ym; j++) { 314 f[j][i].u = x[j][i].u; 315 f[j][i].v = x[j][i].v; 316 f[j][i].omega = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx; 317 f[j][i].temp = x[j][i].temp; 318 } 319 } 320 321 /* Test whether we are on the right edge of the global array */ 322 if (xinte == info->mx) { 323 i = info->mx - 1; 324 xinte = xinte - 1; 325 /* right edge */ 326 for (j=info->ys; j<info->ys+info->ym; j++) { 327 f[j][i].u = x[j][i].u; 328 f[j][i].v = x[j][i].v; 329 f[j][i].omega = x[j][i].omega - (x[j][i].v - x[j][i-1].v)*dhx; 330 f[j][i].temp = x[j][i].temp - (PetscReal)(grashof>0); 331 } 332 } 333 334 /* Compute over the interior points */ 335 for (j=yints; j<yinte; j++) { 336 for (i=xints; i<xinte; i++) { 337 338 /* 339 convective coefficients for upwinding 340 */ 341 vx = x[j][i].u; avx = PetscAbsScalar(vx); 342 vxp = .5*(vx+avx); vxm = .5*(vx-avx); 343 vy = x[j][i].v; avy = PetscAbsScalar(vy); 344 vyp = .5*(vy+avy); vym = .5*(vy-avy); 345 346 /* U velocity */ 347 u = x[j][i].u; 348 udot = user->parabolic ? xdot[j][i].u : 0.; 349 uxx = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx; 350 uyy = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy; 351 f[j][i].u = udot + uxx + uyy - .5*(x[j+1][i].omega-x[j-1][i].omega)*hx; 352 353 /* V velocity */ 354 u = x[j][i].v; 355 udot = user->parabolic ? xdot[j][i].v : 0.; 356 uxx = (2.0*u - x[j][i-1].v - x[j][i+1].v)*hydhx; 357 uyy = (2.0*u - x[j-1][i].v - x[j+1][i].v)*hxdhy; 358 f[j][i].v = udot + uxx + uyy + .5*(x[j][i+1].omega-x[j][i-1].omega)*hy; 359 360 /* Omega */ 361 u = x[j][i].omega; 362 uxx = (2.0*u - x[j][i-1].omega - x[j][i+1].omega)*hydhx; 363 uyy = (2.0*u - x[j-1][i].omega - x[j+1][i].omega)*hxdhy; 364 f[j][i].omega = (xdot[j][i].omega + uxx + uyy 365 + (vxp*(u - x[j][i-1].omega) 366 + vxm*(x[j][i+1].omega - u)) * hy 367 + (vyp*(u - x[j-1][i].omega) 368 + vym*(x[j+1][i].omega - u)) * hx 369 - .5 * grashof * (x[j][i+1].temp - x[j][i-1].temp) * hy); 370 371 /* Temperature */ 372 u = x[j][i].temp; 373 uxx = (2.0*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx; 374 uyy = (2.0*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy; 375 f[j][i].temp = (xdot[j][i].temp + uxx + uyy 376 + prandtl * ((vxp*(u - x[j][i-1].temp) 377 + vxm*(x[j][i+1].temp - u)) * hy 378 + (vyp*(u - x[j-1][i].temp) 379 + vym*(x[j+1][i].temp - u)) * hx)); 380 } 381 } 382 383 /* 384 Flop count (multiply-adds are counted as 2 operations) 385 */ 386 PetscCall(PetscLogFlops(84.0*info->ym*info->xm)); 387 PetscFunctionReturn(0); 388 } 389 390 /*TEST 391 392 test: 393 args: -da_grid_x 20 -da_grid_y 20 -lidvelocity 100 -grashof 1e3 -ts_max_steps 100 -ts_rtol 1e-3 -ts_atol 1e-3 -ts_type rosw -ts_rosw_type ra3pw -ts_monitor -ts_monitor_solution_vtk 'foo-%03D.vts' 394 requires: !complex !single 395 396 test: 397 suffix: 2 398 nsize: 4 399 args: -da_grid_x 20 -da_grid_y 20 -lidvelocity 100 -grashof 1e3 -ts_max_steps 100 -ts_rtol 1e-3 -ts_atol 1e-3 -ts_type rosw -ts_rosw_type ra3pw -ts_monitor -ts_monitor_solution_vtk 'foo-%03D.vts' 400 requires: !complex !single 401 402 test: 403 suffix: 3 404 nsize: 4 405 args: -da_refine 2 -lidvelocity 100 -grashof 1e3 -ts_max_steps 10 -ts_rtol 1e-3 -ts_atol 1e-3 -pc_type none -ts_type beuler -ts_monitor -snes_monitor_short -snes_type aspin -da_overlap 4 406 requires: !complex !single 407 408 test: 409 suffix: 4 410 nsize: 2 411 args: -da_refine 1 -lidvelocity 100 -grashof 1e3 -ts_max_steps 10 -ts_rtol 1e-3 -ts_atol 1e-3 412 requires: !complex !single 413 414 test: 415 suffix: asm 416 nsize: 4 417 args: -da_refine 1 -lidvelocity 100 -grashof 1e3 -ts_max_steps 10 -ts_rtol 1e-3 -ts_atol 1e-3 418 requires: !complex !single 419 420 TEST*/ 421