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