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