xref: /petsc/src/ts/tutorials/advection-diffusion-reaction/ex6.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Model Equations for Advection \n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*
5c4762a1bSJed Brown     Modified from ex3.c
6c4762a1bSJed Brown     Page 9, Section 1.2 Model Equations for Advection-Diffusion
7c4762a1bSJed Brown 
8c4762a1bSJed Brown           u_t + a u_x = 0, 0<= x <= 1.0
9c4762a1bSJed Brown 
10c4762a1bSJed Brown    The initial conditions used here different from the book.
11c4762a1bSJed Brown 
12c4762a1bSJed Brown    Example:
13c4762a1bSJed Brown      ./ex6 -ts_monitor -ts_view_solution -ts_max_steps 100 -ts_monitor_solution draw -draw_pause .1
14c4762a1bSJed Brown      ./ex6 -ts_monitor -ts_max_steps 100 -ts_monitor_lg_error -draw_pause .1
15c4762a1bSJed Brown */
16c4762a1bSJed Brown 
17c4762a1bSJed Brown #include <petscts.h>
18c4762a1bSJed Brown #include <petscdm.h>
19c4762a1bSJed Brown #include <petscdmda.h>
20c4762a1bSJed Brown 
21c4762a1bSJed Brown /*
22c4762a1bSJed Brown    User-defined application context - contains data needed by the
23c4762a1bSJed Brown    application-provided call-back routines.
24c4762a1bSJed Brown */
25c4762a1bSJed Brown typedef struct {
26c4762a1bSJed Brown   PetscReal a; /* advection strength */
27c4762a1bSJed Brown } AppCtx;
28c4762a1bSJed Brown 
29c4762a1bSJed Brown /* User-defined routines */
30c4762a1bSJed Brown extern PetscErrorCode InitialConditions(TS, Vec, AppCtx *);
31c4762a1bSJed Brown extern PetscErrorCode Solution(TS, PetscReal, Vec, AppCtx *);
32c4762a1bSJed Brown extern PetscErrorCode IFunction_LaxFriedrichs(TS, PetscReal, Vec, Vec, Vec, void *);
33c4762a1bSJed Brown extern PetscErrorCode IFunction_LaxWendroff(TS, PetscReal, Vec, Vec, Vec, void *);
34c4762a1bSJed Brown 
359371c9d4SSatish Balay int main(int argc, char **argv) {
36c4762a1bSJed Brown   AppCtx      appctx; /* user-defined application context */
37c4762a1bSJed Brown   TS          ts;     /* timestepping context */
38c4762a1bSJed Brown   Vec         U;      /* approximate solution vector */
39c4762a1bSJed Brown   PetscReal   dt;
40c4762a1bSJed Brown   DM          da;
41c4762a1bSJed Brown   PetscInt    M;
42c4762a1bSJed Brown   PetscMPIInt rank;
43c4762a1bSJed Brown   PetscBool   useLaxWendroff = PETSC_TRUE;
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   /* Initialize program and set problem parameters */
46327415f7SBarry Smith   PetscFunctionBeginUser;
479566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
489566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   appctx.a = -1.0;
519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-a", &appctx.a, NULL));
52c4762a1bSJed Brown 
539566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 60, 1, 1, NULL, &da));
549566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
559566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   /* Create vector data structures for approximate and exact solutions */
589566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &U));
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   /* Create timestepping solver context */
619566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
629566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, da));
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   /* Function evaluation */
659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-useLaxWendroff", &useLaxWendroff, NULL));
66c4762a1bSJed Brown   if (useLaxWendroff) {
67*48a46eb9SPierre Jolivet     if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "... Use Lax-Wendroff finite volume\n"));
689566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(ts, NULL, IFunction_LaxWendroff, &appctx));
69c4762a1bSJed Brown   } else {
70*48a46eb9SPierre Jolivet     if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "... Use Lax-LaxFriedrichs finite difference\n"));
719566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(ts, NULL, IFunction_LaxFriedrichs, &appctx));
72c4762a1bSJed Brown   }
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   /* Customize timestepping solver */
759566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
76c4762a1bSJed Brown   dt = 1.0 / (PetscAbsReal(appctx.a) * M);
779566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
789566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, 100));
799566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 100.0));
809566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
819566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSBEULER));
829566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   /* Evaluate initial conditions */
859566063dSJacob Faibussowitsch   PetscCall(InitialConditions(ts, U, &appctx));
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   /* For testing accuracy of TS with already known solution, e.g., '-ts_monitor_lg_error' */
889566063dSJacob Faibussowitsch   PetscCall(TSSetSolutionFunction(ts, (PetscErrorCode(*)(TS, PetscReal, Vec, void *))Solution, &appctx));
89c4762a1bSJed Brown 
90c4762a1bSJed Brown   /* Run the timestepping solver */
919566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, U));
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   /* Free work space */
949566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
969566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
97c4762a1bSJed Brown 
989566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
99b122ec5aSJacob Faibussowitsch   return 0;
100c4762a1bSJed Brown }
101c4762a1bSJed Brown /* --------------------------------------------------------------------- */
102c4762a1bSJed Brown /*
103c4762a1bSJed Brown    InitialConditions - Computes the solution at the initial time.
104c4762a1bSJed Brown 
105c4762a1bSJed Brown    Input Parameter:
106c4762a1bSJed Brown    u - uninitialized solution vector (global)
107c4762a1bSJed Brown    appctx - user-defined application context
108c4762a1bSJed Brown 
109c4762a1bSJed Brown    Output Parameter:
110c4762a1bSJed Brown    u - vector with solution at initial time (global)
111c4762a1bSJed Brown */
1129371c9d4SSatish Balay PetscErrorCode InitialConditions(TS ts, Vec U, AppCtx *appctx) {
113c4762a1bSJed Brown   PetscScalar *u;
114c4762a1bSJed Brown   PetscInt     i, mstart, mend, um, M;
115c4762a1bSJed Brown   DM           da;
116c4762a1bSJed Brown   PetscReal    h;
117c4762a1bSJed Brown 
1189566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
1199566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &um, 0, 0));
1209566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
121c4762a1bSJed Brown   h    = 1.0 / M;
122c4762a1bSJed Brown   mend = mstart + um;
123c4762a1bSJed Brown   /*
124c4762a1bSJed Brown     Get a pointer to vector data.
125c4762a1bSJed Brown     - For default PETSc vectors, VecGetArray() returns a pointer to
126c4762a1bSJed Brown       the data array.  Otherwise, the routine is implementation dependent.
127c4762a1bSJed Brown     - You MUST call VecRestoreArray() when you no longer need access to
128c4762a1bSJed Brown       the array.
129c4762a1bSJed Brown     - Note that the Fortran interface to VecGetArray() differs from the
130c4762a1bSJed Brown       C version.  See the users manual for details.
131c4762a1bSJed Brown   */
1329566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, U, &u));
133c4762a1bSJed Brown 
134c4762a1bSJed Brown   /*
135c4762a1bSJed Brown      We initialize the solution array by simply writing the solution
136c4762a1bSJed Brown      directly into the array locations.  Alternatively, we could use
137c4762a1bSJed Brown      VecSetValues() or VecSetValuesLocal().
138c4762a1bSJed Brown   */
139c4762a1bSJed Brown   for (i = mstart; i < mend; i++) u[i] = PetscSinReal(PETSC_PI * i * 6. * h) + 3. * PetscSinReal(PETSC_PI * i * 2. * h);
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   /* Restore vector */
1429566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, U, &u));
143c4762a1bSJed Brown   return 0;
144c4762a1bSJed Brown }
145c4762a1bSJed Brown /* --------------------------------------------------------------------- */
146c4762a1bSJed Brown /*
147c4762a1bSJed Brown    Solution - Computes the exact solution at a given time
148c4762a1bSJed Brown 
149c4762a1bSJed Brown    Input Parameters:
150c4762a1bSJed Brown    t - current time
151c4762a1bSJed Brown    solution - vector in which exact solution will be computed
152c4762a1bSJed Brown    appctx - user-defined application context
153c4762a1bSJed Brown 
154c4762a1bSJed Brown    Output Parameter:
155c4762a1bSJed Brown    solution - vector with the newly computed exact solution
156c4762a1bSJed Brown               u(x,t) = sin(6*PI*(x - a*t)) + 3 * sin(2*PI*(x - a*t))
157c4762a1bSJed Brown */
1589371c9d4SSatish Balay PetscErrorCode Solution(TS ts, PetscReal t, Vec U, AppCtx *appctx) {
159c4762a1bSJed Brown   PetscScalar *u;
160c4762a1bSJed Brown   PetscReal    a = appctx->a, h, PI6, PI2;
161c4762a1bSJed Brown   PetscInt     i, mstart, mend, um, M;
162c4762a1bSJed Brown   DM           da;
163c4762a1bSJed Brown 
1649566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
1659566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &um, 0, 0));
1669566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
167c4762a1bSJed Brown   h    = 1.0 / M;
168c4762a1bSJed Brown   mend = mstart + um;
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   /* Get a pointer to vector data. */
1719566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, U, &u));
172c4762a1bSJed Brown 
173c4762a1bSJed Brown   /* u[i] = sin(6*PI*(x[i] - a*t)) + 3 * sin(2*PI*(x[i] - a*t)) */
174c4762a1bSJed Brown   PI6 = PETSC_PI * 6.;
175c4762a1bSJed Brown   PI2 = PETSC_PI * 2.;
1769371c9d4SSatish Balay   for (i = mstart; i < mend; i++) { u[i] = PetscSinReal(PI6 * (i * h - a * t)) + 3. * PetscSinReal(PI2 * (i * h - a * t)); }
177c4762a1bSJed Brown 
178c4762a1bSJed Brown   /* Restore vector */
1799566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, U, &u));
180c4762a1bSJed Brown   return 0;
181c4762a1bSJed Brown }
182c4762a1bSJed Brown 
183c4762a1bSJed Brown /* --------------------------------------------------------------------- */
184c4762a1bSJed Brown /*
185c4762a1bSJed Brown  Use Lax-Friedrichs method to evaluate F(u,t) = du/dt + a *  du/dx
186c4762a1bSJed Brown 
187c4762a1bSJed Brown  See https://en.wikipedia.org/wiki/Lax%E2%80%93Friedrichs_method
188c4762a1bSJed Brown  */
1899371c9d4SSatish Balay PetscErrorCode IFunction_LaxFriedrichs(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx) {
190c4762a1bSJed Brown   AppCtx      *appctx = (AppCtx *)ctx;
191c4762a1bSJed Brown   PetscInt     mstart, mend, M, i, um;
192c4762a1bSJed Brown   DM           da;
193c4762a1bSJed Brown   Vec          Uold, localUold;
194c4762a1bSJed Brown   PetscScalar *uarray, *f, *uoldarray, h, uave, c;
195c4762a1bSJed Brown   PetscReal    dt;
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   PetscFunctionBegin;
1989566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
1999566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &Uold));
200c4762a1bSJed Brown 
2019566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
2029566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
2039566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &um, 0, 0));
204c4762a1bSJed Brown   h    = 1.0 / M;
205c4762a1bSJed Brown   mend = mstart + um;
206c4762a1bSJed Brown   /* printf(" mstart %d, um %d\n",mstart,um); */
207c4762a1bSJed Brown 
2089566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localUold));
2099566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, Uold, INSERT_VALUES, localUold));
2109566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, Uold, INSERT_VALUES, localUold));
211c4762a1bSJed Brown 
212c4762a1bSJed Brown   /* Get pointers to vector data */
2139566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, U, &uarray));
2149566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, localUold, &uoldarray));
2159566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, F, &f));
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   /* advection */
218c4762a1bSJed Brown   c = appctx->a * dt / h; /* Courant-Friedrichs-Lewy number (CFL number) */
219c4762a1bSJed Brown 
220c4762a1bSJed Brown   for (i = mstart; i < mend; i++) {
221c4762a1bSJed Brown     uave = 0.5 * (uoldarray[i - 1] + uoldarray[i + 1]);
222c4762a1bSJed Brown     f[i] = uarray[i] - uave + c * 0.5 * (uoldarray[i + 1] - uoldarray[i - 1]);
223c4762a1bSJed Brown   }
224c4762a1bSJed Brown 
225c4762a1bSJed Brown   /* Restore vectors */
2269566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, U, &uarray));
2279566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, localUold, &uoldarray));
2289566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, F, &f));
2299566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localUold));
230c4762a1bSJed Brown   PetscFunctionReturn(0);
231c4762a1bSJed Brown }
232c4762a1bSJed Brown 
233c4762a1bSJed Brown /*
234c4762a1bSJed Brown  Use Lax-Wendroff method to evaluate F(u,t) = du/dt + a *  du/dx
235c4762a1bSJed Brown */
2369371c9d4SSatish Balay PetscErrorCode IFunction_LaxWendroff(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx) {
237c4762a1bSJed Brown   AppCtx      *appctx = (AppCtx *)ctx;
238c4762a1bSJed Brown   PetscInt     mstart, mend, M, i, um;
239c4762a1bSJed Brown   DM           da;
240c4762a1bSJed Brown   Vec          Uold, localUold;
241c4762a1bSJed Brown   PetscScalar *uarray, *f, *uoldarray, h, RFlux, LFlux, lambda;
242c4762a1bSJed Brown   PetscReal    dt, a;
243c4762a1bSJed Brown 
244c4762a1bSJed Brown   PetscFunctionBegin;
2459566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
2469566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &Uold));
247c4762a1bSJed Brown 
2489566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
2499566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
2509566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &um, 0, 0));
251c4762a1bSJed Brown   h    = 1.0 / M;
252c4762a1bSJed Brown   mend = mstart + um;
253c4762a1bSJed Brown   /* printf(" mstart %d, um %d\n",mstart,um); */
254c4762a1bSJed Brown 
2559566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localUold));
2569566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, Uold, INSERT_VALUES, localUold));
2579566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, Uold, INSERT_VALUES, localUold));
258c4762a1bSJed Brown 
259c4762a1bSJed Brown   /* Get pointers to vector data */
2609566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, U, &uarray));
2619566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, localUold, &uoldarray));
2629566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, F, &f));
263c4762a1bSJed Brown 
264c4762a1bSJed Brown   /* advection -- finite volume (appctx->a < 0 -- can be relaxed?) */
265c4762a1bSJed Brown   lambda = dt / h;
266c4762a1bSJed Brown   a      = appctx->a;
267c4762a1bSJed Brown 
268c4762a1bSJed Brown   for (i = mstart; i < mend; i++) {
269c4762a1bSJed Brown     RFlux = 0.5 * a * (uoldarray[i + 1] + uoldarray[i]) - a * a * 0.5 * lambda * (uoldarray[i + 1] - uoldarray[i]);
270c4762a1bSJed Brown     LFlux = 0.5 * a * (uoldarray[i - 1] + uoldarray[i]) - a * a * 0.5 * lambda * (uoldarray[i] - uoldarray[i - 1]);
271c4762a1bSJed Brown     f[i]  = uarray[i] - uoldarray[i] + lambda * (RFlux - LFlux);
272c4762a1bSJed Brown   }
273c4762a1bSJed Brown 
274c4762a1bSJed Brown   /* Restore vectors */
2759566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, U, &uarray));
2769566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, localUold, &uoldarray));
2779566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, F, &f));
2789566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localUold));
279c4762a1bSJed Brown   PetscFunctionReturn(0);
280c4762a1bSJed Brown }
281c4762a1bSJed Brown 
282c4762a1bSJed Brown /*TEST
283c4762a1bSJed Brown 
284c4762a1bSJed Brown    test:
285c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_monitor
286c4762a1bSJed Brown 
287c4762a1bSJed Brown    test:
288c4762a1bSJed Brown       suffix: 2
289c4762a1bSJed Brown       nsize: 3
290c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_monitor
291c4762a1bSJed Brown       output_file: output/ex6_1.out
292c4762a1bSJed Brown 
293c4762a1bSJed Brown    test:
294c4762a1bSJed Brown       suffix: 3
295c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_monitor -useLaxWendroff false
296c4762a1bSJed Brown 
297c4762a1bSJed Brown    test:
298c4762a1bSJed Brown       suffix: 4
299c4762a1bSJed Brown       nsize: 3
300c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_monitor -useLaxWendroff false
301c4762a1bSJed Brown       output_file: output/ex6_3.out
302c4762a1bSJed Brown 
303c4762a1bSJed Brown TEST*/
304