xref: /petsc/src/ts/tutorials/advection-diffusion-reaction/ex6.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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 
35*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
36*d71ae5a4SJacob Faibussowitsch {
37c4762a1bSJed Brown   AppCtx      appctx; /* user-defined application context */
38c4762a1bSJed Brown   TS          ts;     /* timestepping context */
39c4762a1bSJed Brown   Vec         U;      /* approximate solution vector */
40c4762a1bSJed Brown   PetscReal   dt;
41c4762a1bSJed Brown   DM          da;
42c4762a1bSJed Brown   PetscInt    M;
43c4762a1bSJed Brown   PetscMPIInt rank;
44c4762a1bSJed Brown   PetscBool   useLaxWendroff = PETSC_TRUE;
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   /* Initialize program and set problem parameters */
47327415f7SBarry Smith   PetscFunctionBeginUser;
489566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
499566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   appctx.a = -1.0;
529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-a", &appctx.a, NULL));
53c4762a1bSJed Brown 
549566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 60, 1, 1, NULL, &da));
559566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
569566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   /* Create vector data structures for approximate and exact solutions */
599566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &U));
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /* Create timestepping solver context */
629566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
639566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, da));
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   /* Function evaluation */
669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-useLaxWendroff", &useLaxWendroff, NULL));
67c4762a1bSJed Brown   if (useLaxWendroff) {
6848a46eb9SPierre Jolivet     if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "... Use Lax-Wendroff finite volume\n"));
699566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(ts, NULL, IFunction_LaxWendroff, &appctx));
70c4762a1bSJed Brown   } else {
7148a46eb9SPierre Jolivet     if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "... Use Lax-LaxFriedrichs finite difference\n"));
729566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(ts, NULL, IFunction_LaxFriedrichs, &appctx));
73c4762a1bSJed Brown   }
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   /* Customize timestepping solver */
769566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
77c4762a1bSJed Brown   dt = 1.0 / (PetscAbsReal(appctx.a) * M);
789566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
799566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, 100));
809566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 100.0));
819566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
829566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSBEULER));
839566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   /* Evaluate initial conditions */
869566063dSJacob Faibussowitsch   PetscCall(InitialConditions(ts, U, &appctx));
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   /* For testing accuracy of TS with already known solution, e.g., '-ts_monitor_lg_error' */
899566063dSJacob Faibussowitsch   PetscCall(TSSetSolutionFunction(ts, (PetscErrorCode(*)(TS, PetscReal, Vec, void *))Solution, &appctx));
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   /* Run the timestepping solver */
929566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, U));
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   /* Free work space */
959566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
969566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
979566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
98c4762a1bSJed Brown 
999566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
100b122ec5aSJacob Faibussowitsch   return 0;
101c4762a1bSJed Brown }
102c4762a1bSJed Brown /* --------------------------------------------------------------------- */
103c4762a1bSJed Brown /*
104c4762a1bSJed Brown    InitialConditions - Computes the solution at the initial time.
105c4762a1bSJed Brown 
106c4762a1bSJed Brown    Input Parameter:
107c4762a1bSJed Brown    u - uninitialized solution vector (global)
108c4762a1bSJed Brown    appctx - user-defined application context
109c4762a1bSJed Brown 
110c4762a1bSJed Brown    Output Parameter:
111c4762a1bSJed Brown    u - vector with solution at initial time (global)
112c4762a1bSJed Brown */
113*d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialConditions(TS ts, Vec U, AppCtx *appctx)
114*d71ae5a4SJacob Faibussowitsch {
115c4762a1bSJed Brown   PetscScalar *u;
116c4762a1bSJed Brown   PetscInt     i, mstart, mend, um, M;
117c4762a1bSJed Brown   DM           da;
118c4762a1bSJed Brown   PetscReal    h;
119c4762a1bSJed Brown 
1209566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
1219566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &um, 0, 0));
1229566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
123c4762a1bSJed Brown   h    = 1.0 / M;
124c4762a1bSJed Brown   mend = mstart + um;
125c4762a1bSJed Brown   /*
126c4762a1bSJed Brown     Get a pointer to vector data.
127c4762a1bSJed Brown     - For default PETSc vectors, VecGetArray() returns a pointer to
128c4762a1bSJed Brown       the data array.  Otherwise, the routine is implementation dependent.
129c4762a1bSJed Brown     - You MUST call VecRestoreArray() when you no longer need access to
130c4762a1bSJed Brown       the array.
131c4762a1bSJed Brown     - Note that the Fortran interface to VecGetArray() differs from the
132c4762a1bSJed Brown       C version.  See the users manual for details.
133c4762a1bSJed Brown   */
1349566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, U, &u));
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   /*
137c4762a1bSJed Brown      We initialize the solution array by simply writing the solution
138c4762a1bSJed Brown      directly into the array locations.  Alternatively, we could use
139c4762a1bSJed Brown      VecSetValues() or VecSetValuesLocal().
140c4762a1bSJed Brown   */
141c4762a1bSJed Brown   for (i = mstart; i < mend; i++) u[i] = PetscSinReal(PETSC_PI * i * 6. * h) + 3. * PetscSinReal(PETSC_PI * i * 2. * h);
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   /* Restore vector */
1449566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, U, &u));
145c4762a1bSJed Brown   return 0;
146c4762a1bSJed Brown }
147c4762a1bSJed Brown /* --------------------------------------------------------------------- */
148c4762a1bSJed Brown /*
149c4762a1bSJed Brown    Solution - Computes the exact solution at a given time
150c4762a1bSJed Brown 
151c4762a1bSJed Brown    Input Parameters:
152c4762a1bSJed Brown    t - current time
153c4762a1bSJed Brown    solution - vector in which exact solution will be computed
154c4762a1bSJed Brown    appctx - user-defined application context
155c4762a1bSJed Brown 
156c4762a1bSJed Brown    Output Parameter:
157c4762a1bSJed Brown    solution - vector with the newly computed exact solution
158c4762a1bSJed Brown               u(x,t) = sin(6*PI*(x - a*t)) + 3 * sin(2*PI*(x - a*t))
159c4762a1bSJed Brown */
160*d71ae5a4SJacob Faibussowitsch PetscErrorCode Solution(TS ts, PetscReal t, Vec U, AppCtx *appctx)
161*d71ae5a4SJacob Faibussowitsch {
162c4762a1bSJed Brown   PetscScalar *u;
163c4762a1bSJed Brown   PetscReal    a = appctx->a, h, PI6, PI2;
164c4762a1bSJed Brown   PetscInt     i, mstart, mend, um, M;
165c4762a1bSJed Brown   DM           da;
166c4762a1bSJed Brown 
1679566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
1689566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &um, 0, 0));
1699566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
170c4762a1bSJed Brown   h    = 1.0 / M;
171c4762a1bSJed Brown   mend = mstart + um;
172c4762a1bSJed Brown 
173c4762a1bSJed Brown   /* Get a pointer to vector data. */
1749566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, U, &u));
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   /* u[i] = sin(6*PI*(x[i] - a*t)) + 3 * sin(2*PI*(x[i] - a*t)) */
177c4762a1bSJed Brown   PI6 = PETSC_PI * 6.;
178c4762a1bSJed Brown   PI2 = PETSC_PI * 2.;
179ad540459SPierre Jolivet   for (i = mstart; i < mend; i++) u[i] = PetscSinReal(PI6 * (i * h - a * t)) + 3. * PetscSinReal(PI2 * (i * h - a * t));
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   /* Restore vector */
1829566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, U, &u));
183c4762a1bSJed Brown   return 0;
184c4762a1bSJed Brown }
185c4762a1bSJed Brown 
186c4762a1bSJed Brown /* --------------------------------------------------------------------- */
187c4762a1bSJed Brown /*
188c4762a1bSJed Brown  Use Lax-Friedrichs method to evaluate F(u,t) = du/dt + a *  du/dx
189c4762a1bSJed Brown 
190c4762a1bSJed Brown  See https://en.wikipedia.org/wiki/Lax%E2%80%93Friedrichs_method
191c4762a1bSJed Brown  */
192*d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction_LaxFriedrichs(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx)
193*d71ae5a4SJacob Faibussowitsch {
194c4762a1bSJed Brown   AppCtx      *appctx = (AppCtx *)ctx;
195c4762a1bSJed Brown   PetscInt     mstart, mend, M, i, um;
196c4762a1bSJed Brown   DM           da;
197c4762a1bSJed Brown   Vec          Uold, localUold;
198c4762a1bSJed Brown   PetscScalar *uarray, *f, *uoldarray, h, uave, c;
199c4762a1bSJed Brown   PetscReal    dt;
200c4762a1bSJed Brown 
201c4762a1bSJed Brown   PetscFunctionBegin;
2029566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
2039566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &Uold));
204c4762a1bSJed Brown 
2059566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
2069566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
2079566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &um, 0, 0));
208c4762a1bSJed Brown   h    = 1.0 / M;
209c4762a1bSJed Brown   mend = mstart + um;
210c4762a1bSJed Brown   /* printf(" mstart %d, um %d\n",mstart,um); */
211c4762a1bSJed Brown 
2129566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localUold));
2139566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, Uold, INSERT_VALUES, localUold));
2149566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, Uold, INSERT_VALUES, localUold));
215c4762a1bSJed Brown 
216c4762a1bSJed Brown   /* Get pointers to vector data */
2179566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, U, &uarray));
2189566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, localUold, &uoldarray));
2199566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, F, &f));
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   /* advection */
222c4762a1bSJed Brown   c = appctx->a * dt / h; /* Courant-Friedrichs-Lewy number (CFL number) */
223c4762a1bSJed Brown 
224c4762a1bSJed Brown   for (i = mstart; i < mend; i++) {
225c4762a1bSJed Brown     uave = 0.5 * (uoldarray[i - 1] + uoldarray[i + 1]);
226c4762a1bSJed Brown     f[i] = uarray[i] - uave + c * 0.5 * (uoldarray[i + 1] - uoldarray[i - 1]);
227c4762a1bSJed Brown   }
228c4762a1bSJed Brown 
229c4762a1bSJed Brown   /* Restore vectors */
2309566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, U, &uarray));
2319566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, localUold, &uoldarray));
2329566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, F, &f));
2339566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localUold));
234c4762a1bSJed Brown   PetscFunctionReturn(0);
235c4762a1bSJed Brown }
236c4762a1bSJed Brown 
237c4762a1bSJed Brown /*
238c4762a1bSJed Brown  Use Lax-Wendroff method to evaluate F(u,t) = du/dt + a *  du/dx
239c4762a1bSJed Brown */
240*d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction_LaxWendroff(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx)
241*d71ae5a4SJacob Faibussowitsch {
242c4762a1bSJed Brown   AppCtx      *appctx = (AppCtx *)ctx;
243c4762a1bSJed Brown   PetscInt     mstart, mend, M, i, um;
244c4762a1bSJed Brown   DM           da;
245c4762a1bSJed Brown   Vec          Uold, localUold;
246c4762a1bSJed Brown   PetscScalar *uarray, *f, *uoldarray, h, RFlux, LFlux, lambda;
247c4762a1bSJed Brown   PetscReal    dt, a;
248c4762a1bSJed Brown 
249c4762a1bSJed Brown   PetscFunctionBegin;
2509566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
2519566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &Uold));
252c4762a1bSJed Brown 
2539566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
2549566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
2559566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &um, 0, 0));
256c4762a1bSJed Brown   h    = 1.0 / M;
257c4762a1bSJed Brown   mend = mstart + um;
258c4762a1bSJed Brown   /* printf(" mstart %d, um %d\n",mstart,um); */
259c4762a1bSJed Brown 
2609566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localUold));
2619566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, Uold, INSERT_VALUES, localUold));
2629566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, Uold, INSERT_VALUES, localUold));
263c4762a1bSJed Brown 
264c4762a1bSJed Brown   /* Get pointers to vector data */
2659566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, U, &uarray));
2669566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, localUold, &uoldarray));
2679566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, F, &f));
268c4762a1bSJed Brown 
269c4762a1bSJed Brown   /* advection -- finite volume (appctx->a < 0 -- can be relaxed?) */
270c4762a1bSJed Brown   lambda = dt / h;
271c4762a1bSJed Brown   a      = appctx->a;
272c4762a1bSJed Brown 
273c4762a1bSJed Brown   for (i = mstart; i < mend; i++) {
274c4762a1bSJed Brown     RFlux = 0.5 * a * (uoldarray[i + 1] + uoldarray[i]) - a * a * 0.5 * lambda * (uoldarray[i + 1] - uoldarray[i]);
275c4762a1bSJed Brown     LFlux = 0.5 * a * (uoldarray[i - 1] + uoldarray[i]) - a * a * 0.5 * lambda * (uoldarray[i] - uoldarray[i - 1]);
276c4762a1bSJed Brown     f[i]  = uarray[i] - uoldarray[i] + lambda * (RFlux - LFlux);
277c4762a1bSJed Brown   }
278c4762a1bSJed Brown 
279c4762a1bSJed Brown   /* Restore vectors */
2809566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, U, &uarray));
2819566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, localUold, &uoldarray));
2829566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, F, &f));
2839566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localUold));
284c4762a1bSJed Brown   PetscFunctionReturn(0);
285c4762a1bSJed Brown }
286c4762a1bSJed Brown 
287c4762a1bSJed Brown /*TEST
288c4762a1bSJed Brown 
289c4762a1bSJed Brown    test:
290c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_monitor
291c4762a1bSJed Brown 
292c4762a1bSJed Brown    test:
293c4762a1bSJed Brown       suffix: 2
294c4762a1bSJed Brown       nsize: 3
295c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_monitor
296c4762a1bSJed Brown       output_file: output/ex6_1.out
297c4762a1bSJed Brown 
298c4762a1bSJed Brown    test:
299c4762a1bSJed Brown       suffix: 3
300c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_monitor -useLaxWendroff false
301c4762a1bSJed Brown 
302c4762a1bSJed Brown    test:
303c4762a1bSJed Brown       suffix: 4
304c4762a1bSJed Brown       nsize: 3
305c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_monitor -useLaxWendroff false
306c4762a1bSJed Brown       output_file: output/ex6_3.out
307c4762a1bSJed Brown 
308c4762a1bSJed Brown TEST*/
309