xref: /petsc/src/ts/tutorials/advection-diffusion-reaction/ex3.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1 
2 static char help[] = "Model Equations for Advection-Diffusion\n";
3 
4 /*
5     Page 9, Section 1.2 Model Equations for Advection-Diffusion
6 
7           u_t = a u_x + d u_xx
8 
9    The initial conditions used here different then in the book.
10 
11 */
12 
13 /*
14      Helpful runtime linear solver options:
15            -pc_type mg -da_refine 2 -snes_monitor -ksp_monitor -ts_view   (geometric multigrid with three levels)
16 
17 */
18 
19 /*
20    Include "petscts.h" so that we can use TS solvers.  Note that this file
21    automatically includes:
22      petscsys.h       - base PETSc routines   petscvec.h  - vectors
23      petscmat.h  - matrices
24      petscis.h     - index sets            petscksp.h  - Krylov subspace methods
25      petscviewer.h - viewers               petscpc.h   - preconditioners
26      petscksp.h   - linear solvers        petscsnes.h - nonlinear solvers
27 */
28 
29 #include <petscts.h>
30 #include <petscdm.h>
31 #include <petscdmda.h>
32 
33 /*
34    User-defined application context - contains data needed by the
35    application-provided call-back routines.
36 */
37 typedef struct {
38   PetscScalar a, d; /* advection and diffusion strength */
39   PetscBool   upwind;
40 } AppCtx;
41 
42 /*
43    User-defined routines
44 */
45 extern PetscErrorCode InitialConditions(TS, Vec, AppCtx *);
46 extern PetscErrorCode RHSMatrixHeat(TS, PetscReal, Vec, Mat, Mat, void *);
47 extern PetscErrorCode Solution(TS, PetscReal, Vec, AppCtx *);
48 
49 int main(int argc, char **argv)
50 {
51   AppCtx    appctx; /* user-defined application context */
52   TS        ts;     /* timestepping context */
53   Vec       U;      /* approximate solution vector */
54   PetscReal dt;
55   DM        da;
56   PetscInt  M;
57 
58   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
59      Initialize program and set problem parameters
60      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61 
62   PetscFunctionBeginUser;
63   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
64   appctx.a = 1.0;
65   appctx.d = 0.0;
66   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-a", &appctx.a, NULL));
67   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-d", &appctx.d, NULL));
68   appctx.upwind = PETSC_TRUE;
69   PetscCall(PetscOptionsGetBool(NULL, NULL, "-upwind", &appctx.upwind, NULL));
70 
71   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 60, 1, 1, NULL, &da));
72   PetscCall(DMSetFromOptions(da));
73   PetscCall(DMSetUp(da));
74   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75      Create vector data structures
76      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
77 
78   /*
79      Create vector data structures for approximate and exact solutions
80   */
81   PetscCall(DMCreateGlobalVector(da, &U));
82 
83   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84      Create timestepping solver context
85      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86 
87   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
88   PetscCall(TSSetDM(ts, da));
89 
90   /*
91       For linear problems with a time-dependent f(U,t) in the equation
92      u_t = f(u,t), the user provides the discretized right-hand-side
93       as a time-dependent matrix.
94   */
95   PetscCall(TSSetRHSFunction(ts, NULL, TSComputeRHSFunctionLinear, &appctx));
96   PetscCall(TSSetRHSJacobian(ts, NULL, NULL, RHSMatrixHeat, &appctx));
97   PetscCall(TSSetSolutionFunction(ts, (PetscErrorCode(*)(TS, PetscReal, Vec, void *))Solution, &appctx));
98 
99   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100      Customize timestepping solver:
101        - Set timestepping duration info
102      Then set runtime options, which can override these defaults.
103      For example,
104           -ts_max_steps <maxsteps> -ts_max_time <maxtime>
105      to override the defaults set by TSSetMaxSteps()/TSSetMaxTime().
106      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
107 
108   PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
109   dt = .48 / (M * M);
110   PetscCall(TSSetTimeStep(ts, dt));
111   PetscCall(TSSetMaxSteps(ts, 1000));
112   PetscCall(TSSetMaxTime(ts, 100.0));
113   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
114   PetscCall(TSSetType(ts, TSARKIMEX));
115   PetscCall(TSSetFromOptions(ts));
116 
117   /*
118      Evaluate initial conditions
119   */
120   PetscCall(InitialConditions(ts, U, &appctx));
121 
122   /*
123      Run the timestepping solver
124   */
125   PetscCall(TSSolve(ts, U));
126 
127   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128      Free work space.  All PETSc objects should be destroyed when they
129      are no longer needed.
130      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131 
132   PetscCall(TSDestroy(&ts));
133   PetscCall(VecDestroy(&U));
134   PetscCall(DMDestroy(&da));
135 
136   /*
137      Always call PetscFinalize() before exiting a program.  This routine
138        - finalizes the PETSc libraries as well as MPI
139        - provides summary and diagnostic information if certain runtime
140          options are chosen (e.g., -log_view).
141   */
142   PetscCall(PetscFinalize());
143   return 0;
144 }
145 /* --------------------------------------------------------------------- */
146 /*
147    InitialConditions - Computes the solution at the initial time.
148 
149    Input Parameter:
150    u - uninitialized solution vector (global)
151    appctx - user-defined application context
152 
153    Output Parameter:
154    u - vector with solution at initial time (global)
155 */
156 PetscErrorCode InitialConditions(TS ts, Vec U, AppCtx *appctx)
157 {
158   PetscScalar *u, h;
159   PetscInt     i, mstart, mend, xm, M;
160   DM           da;
161 
162   PetscCall(TSGetDM(ts, &da));
163   PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &xm, 0, 0));
164   PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
165   h    = 1.0 / M;
166   mend = mstart + xm;
167   /*
168     Get a pointer to vector data.
169     - For default PETSc vectors, VecGetArray() returns a pointer to
170       the data array.  Otherwise, the routine is implementation dependent.
171     - You MUST call VecRestoreArray() when you no longer need access to
172       the array.
173     - Note that the Fortran interface to VecGetArray() differs from the
174       C version.  See the users manual for details.
175   */
176   PetscCall(DMDAVecGetArray(da, U, &u));
177 
178   /*
179      We initialize the solution array by simply writing the solution
180      directly into the array locations.  Alternatively, we could use
181      VecSetValues() or VecSetValuesLocal().
182   */
183   for (i = mstart; i < mend; i++) u[i] = PetscSinScalar(PETSC_PI * i * 6. * h) + 3. * PetscSinScalar(PETSC_PI * i * 2. * h);
184 
185   /*
186      Restore vector
187   */
188   PetscCall(DMDAVecRestoreArray(da, U, &u));
189   return 0;
190 }
191 /* --------------------------------------------------------------------- */
192 /*
193    Solution - Computes the exact solution at a given time.
194 
195    Input Parameters:
196    t - current time
197    solution - vector in which exact solution will be computed
198    appctx - user-defined application context
199 
200    Output Parameter:
201    solution - vector with the newly computed exact solution
202 */
203 PetscErrorCode Solution(TS ts, PetscReal t, Vec U, AppCtx *appctx)
204 {
205   PetscScalar *u, ex1, ex2, sc1, sc2, h;
206   PetscInt     i, mstart, mend, xm, M;
207   DM           da;
208 
209   PetscCall(TSGetDM(ts, &da));
210   PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &xm, 0, 0));
211   PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
212   h    = 1.0 / M;
213   mend = mstart + xm;
214   /*
215      Get a pointer to vector data.
216   */
217   PetscCall(DMDAVecGetArray(da, U, &u));
218 
219   /*
220      Simply write the solution directly into the array locations.
221      Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
222   */
223   ex1 = PetscExpScalar(-36. * PETSC_PI * PETSC_PI * appctx->d * t);
224   ex2 = PetscExpScalar(-4. * PETSC_PI * PETSC_PI * appctx->d * t);
225   sc1 = PETSC_PI * 6. * h;
226   sc2 = PETSC_PI * 2. * h;
227   for (i = mstart; i < mend; i++) u[i] = PetscSinScalar(sc1 * (PetscReal)i + appctx->a * PETSC_PI * 6. * t) * ex1 + 3. * PetscSinScalar(sc2 * (PetscReal)i + appctx->a * PETSC_PI * 2. * t) * ex2;
228 
229   /*
230      Restore vector
231   */
232   PetscCall(DMDAVecRestoreArray(da, U, &u));
233   return 0;
234 }
235 
236 /* --------------------------------------------------------------------- */
237 /*
238    RHSMatrixHeat - User-provided routine to compute the right-hand-side
239    matrix for the heat equation.
240 
241    Input Parameters:
242    ts - the TS context
243    t - current time
244    global_in - global input vector
245    dummy - optional user-defined context, as set by TSetRHSJacobian()
246 
247    Output Parameters:
248    AA - Jacobian matrix
249    BB - optionally different preconditioning matrix
250    str - flag indicating matrix structure
251 
252    Notes:
253    Recall that MatSetValues() uses 0-based row and column numbers
254    in Fortran as well as in C.
255 */
256 PetscErrorCode RHSMatrixHeat(TS ts, PetscReal t, Vec U, Mat AA, Mat BB, void *ctx)
257 {
258   Mat         A      = AA;            /* Jacobian matrix */
259   AppCtx     *appctx = (AppCtx *)ctx; /* user-defined application context */
260   PetscInt    mstart, mend;
261   PetscInt    i, idx[3], M, xm;
262   PetscScalar v[3], h;
263   DM          da;
264 
265   PetscCall(TSGetDM(ts, &da));
266   PetscCall(DMDAGetInfo(da, 0, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
267   PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &xm, 0, 0));
268   h    = 1.0 / M;
269   mend = mstart + xm;
270   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
271      Compute entries for the locally owned part of the matrix
272      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
273   /*
274      Set matrix rows corresponding to boundary data
275   */
276 
277   /* diffusion */
278   v[0] = appctx->d / (h * h);
279   v[1] = -2.0 * appctx->d / (h * h);
280   v[2] = appctx->d / (h * h);
281   if (!mstart) {
282     idx[0] = M - 1;
283     idx[1] = 0;
284     idx[2] = 1;
285     PetscCall(MatSetValues(A, 1, &mstart, 3, idx, v, INSERT_VALUES));
286     mstart++;
287   }
288 
289   if (mend == M) {
290     mend--;
291     idx[0] = M - 2;
292     idx[1] = M - 1;
293     idx[2] = 0;
294     PetscCall(MatSetValues(A, 1, &mend, 3, idx, v, INSERT_VALUES));
295   }
296 
297   /*
298      Set matrix rows corresponding to interior data.  We construct the
299      matrix one row at a time.
300   */
301   for (i = mstart; i < mend; i++) {
302     idx[0] = i - 1;
303     idx[1] = i;
304     idx[2] = i + 1;
305     PetscCall(MatSetValues(A, 1, &i, 3, idx, v, INSERT_VALUES));
306   }
307   PetscCall(MatAssemblyBegin(A, MAT_FLUSH_ASSEMBLY));
308   PetscCall(MatAssemblyEnd(A, MAT_FLUSH_ASSEMBLY));
309 
310   PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &xm, 0, 0));
311   mend = mstart + xm;
312   if (!appctx->upwind) {
313     /* advection -- centered differencing */
314     v[0] = -.5 * appctx->a / (h);
315     v[1] = .5 * appctx->a / (h);
316     if (!mstart) {
317       idx[0] = M - 1;
318       idx[1] = 1;
319       PetscCall(MatSetValues(A, 1, &mstart, 2, idx, v, ADD_VALUES));
320       mstart++;
321     }
322 
323     if (mend == M) {
324       mend--;
325       idx[0] = M - 2;
326       idx[1] = 0;
327       PetscCall(MatSetValues(A, 1, &mend, 2, idx, v, ADD_VALUES));
328     }
329 
330     for (i = mstart; i < mend; i++) {
331       idx[0] = i - 1;
332       idx[1] = i + 1;
333       PetscCall(MatSetValues(A, 1, &i, 2, idx, v, ADD_VALUES));
334     }
335   } else {
336     /* advection -- upwinding */
337     v[0] = -appctx->a / (h);
338     v[1] = appctx->a / (h);
339     if (!mstart) {
340       idx[0] = 0;
341       idx[1] = 1;
342       PetscCall(MatSetValues(A, 1, &mstart, 2, idx, v, ADD_VALUES));
343       mstart++;
344     }
345 
346     if (mend == M) {
347       mend--;
348       idx[0] = M - 1;
349       idx[1] = 0;
350       PetscCall(MatSetValues(A, 1, &mend, 2, idx, v, ADD_VALUES));
351     }
352 
353     for (i = mstart; i < mend; i++) {
354       idx[0] = i;
355       idx[1] = i + 1;
356       PetscCall(MatSetValues(A, 1, &i, 2, idx, v, ADD_VALUES));
357     }
358   }
359 
360   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
361      Complete the matrix assembly process and set some options
362      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
363   /*
364      Assemble matrix, using the 2-step process:
365        MatAssemblyBegin(), MatAssemblyEnd()
366      Computations can be done while messages are in transition
367      by placing code between these two statements.
368   */
369   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
370   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
371 
372   /*
373      Set and option to indicate that we will never add a new nonzero location
374      to the matrix. If we do, it will generate an error.
375   */
376   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
377   return 0;
378 }
379 
380 /*TEST
381 
382    test:
383       args: -pc_type mg -da_refine 2  -ts_view  -ts_monitor -ts_max_time .3 -mg_levels_ksp_max_it 3
384       requires: double
385       filter: grep -v "total number of"
386 
387    test:
388       suffix: 2
389       args:  -pc_type mg -da_refine 2  -ts_view  -ts_monitor_draw_solution -ts_monitor -ts_max_time .3 -mg_levels_ksp_max_it 3
390       requires: x
391       output_file: output/ex3_1.out
392       requires: double
393       filter: grep -v "total number of"
394 
395 TEST*/
396