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