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