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