xref: /petsc/src/ts/tutorials/ex13.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Time-dependent PDE in 2d. Simplified from ex7.c for illustrating how to use TS on a structured domain. \n";
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown    u_t = uxx + uyy
5c4762a1bSJed Brown    0 < x < 1, 0 < y < 1;
6c4762a1bSJed Brown    At t=0: u(x,y) = exp(c*r*r*r), if r=PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5)) < .125
7c4762a1bSJed Brown            u(x,y) = 0.0           if r >= .125
8c4762a1bSJed Brown 
9c4762a1bSJed Brown     mpiexec -n 2 ./ex13 -da_grid_x 40 -da_grid_y 40 -ts_max_steps 2 -snes_monitor -ksp_monitor
10c4762a1bSJed Brown     mpiexec -n 1 ./ex13 -snes_fd_color -ts_monitor_draw_solution
11c4762a1bSJed Brown     mpiexec -n 2 ./ex13 -ts_type sundials -ts_monitor
12c4762a1bSJed Brown */
13c4762a1bSJed Brown 
14c4762a1bSJed Brown #include <petscdm.h>
15c4762a1bSJed Brown #include <petscdmda.h>
16c4762a1bSJed Brown #include <petscts.h>
17c4762a1bSJed Brown 
18c4762a1bSJed Brown /*
19c4762a1bSJed Brown    User-defined data structures and routines
20c4762a1bSJed Brown */
21c4762a1bSJed Brown typedef struct {
22c4762a1bSJed Brown   PetscReal c;
23c4762a1bSJed Brown } AppCtx;
24c4762a1bSJed Brown 
25c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS, PetscReal, Vec, Vec, void *);
26c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS, PetscReal, Vec, Mat, Mat, void *);
27c4762a1bSJed Brown extern PetscErrorCode FormInitialSolution(DM, Vec, void *);
28c4762a1bSJed Brown 
29d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
30d71ae5a4SJacob Faibussowitsch {
31c4762a1bSJed Brown   TS        ts;    /* nonlinear solver */
32c4762a1bSJed Brown   Vec       u, r;  /* solution, residual vector */
33c4762a1bSJed Brown   Mat       J;     /* Jacobian matrix */
34c4762a1bSJed Brown   PetscInt  steps; /* iterations for convergence */
35c4762a1bSJed Brown   DM        da;
36c4762a1bSJed Brown   PetscReal ftime, dt;
37c4762a1bSJed Brown   AppCtx    user; /* user-defined work context */
38c4762a1bSJed Brown 
39327415f7SBarry Smith   PetscFunctionBeginUser;
409566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
41c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
42c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
43c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
449566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, 8, 8, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &da));
459566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
469566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
47c4762a1bSJed Brown 
48c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
49c4762a1bSJed Brown      Extract global vectors from DMDA;
50c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
519566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &u));
529566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &r));
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   /* Initialize user application context */
55c4762a1bSJed Brown   user.c = -30.0;
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58c4762a1bSJed Brown      Create timestepping solver context
59c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
609566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
619566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, da));
629566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSBEULER));
639566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, r, RHSFunction, &user));
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   /* Set Jacobian */
669566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(da, MATAIJ));
679566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(da, &J));
689566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts, J, J, RHSJacobian, NULL));
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   ftime = 1.0;
719566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, ftime));
729566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75c4762a1bSJed Brown      Set initial conditions
76c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
779566063dSJacob Faibussowitsch   PetscCall(FormInitialSolution(da, u, &user));
78c4762a1bSJed Brown   dt = .01;
799566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
82c4762a1bSJed Brown      Set runtime options
83c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
849566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
87c4762a1bSJed Brown      Solve nonlinear system
88c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
899566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, u));
909566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
919566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94c4762a1bSJed Brown      Free work space.
95c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
969566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
979566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
989566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
999566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1009566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
101c4762a1bSJed Brown 
1029566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
103b122ec5aSJacob Faibussowitsch   return 0;
104c4762a1bSJed Brown }
105c4762a1bSJed Brown /* ------------------------------------------------------------------- */
106c4762a1bSJed Brown /*
107c4762a1bSJed Brown    RHSFunction - Evaluates nonlinear function, F(u).
108c4762a1bSJed Brown 
109c4762a1bSJed Brown    Input Parameters:
110c4762a1bSJed Brown .  ts - the TS context
111c4762a1bSJed Brown .  U - input vector
112c4762a1bSJed Brown .  ptr - optional user-defined context, as set by TSSetFunction()
113c4762a1bSJed Brown 
114c4762a1bSJed Brown    Output Parameter:
115c4762a1bSJed Brown .  F - function vector
116c4762a1bSJed Brown  */
117d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction(TS ts, PetscReal ftime, Vec U, Vec F, void *ptr)
118d71ae5a4SJacob Faibussowitsch {
119c4762a1bSJed Brown   /* PETSC_UNUSED AppCtx *user=(AppCtx*)ptr; */
120c4762a1bSJed Brown   DM          da;
121c4762a1bSJed Brown   PetscInt    i, j, Mx, My, xs, ys, xm, ym;
122c4762a1bSJed Brown   PetscReal   two = 2.0, hx, hy, sx, sy;
123c4762a1bSJed Brown   PetscScalar u, uxx, uyy, **uarray, **f;
124c4762a1bSJed Brown   Vec         localU;
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   PetscFunctionBeginUser;
1279566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
1289566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localU));
1299566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
130c4762a1bSJed Brown 
1319371c9d4SSatish Balay   hx = 1.0 / (PetscReal)(Mx - 1);
1329371c9d4SSatish Balay   sx = 1.0 / (hx * hx);
1339371c9d4SSatish Balay   hy = 1.0 / (PetscReal)(My - 1);
1349371c9d4SSatish Balay   sy = 1.0 / (hy * hy);
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   /*
137c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
138c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
139c4762a1bSJed Brown      By placing code between these two statements, computations can be
140c4762a1bSJed Brown      done while messages are in transition.
141c4762a1bSJed Brown   */
1429566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU));
1439566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU));
144c4762a1bSJed Brown 
145c4762a1bSJed Brown   /* Get pointers to vector data */
1469566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, localU, &uarray));
1479566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, F, &f));
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   /* Get local grid boundaries */
1509566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
153c4762a1bSJed Brown   for (j = ys; j < ys + ym; j++) {
154c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
155c4762a1bSJed Brown       if (i == 0 || j == 0 || i == Mx - 1 || j == My - 1) {
156c4762a1bSJed Brown         f[j][i] = uarray[j][i];
157c4762a1bSJed Brown         continue;
158c4762a1bSJed Brown       }
159c4762a1bSJed Brown       u       = uarray[j][i];
160c4762a1bSJed Brown       uxx     = (-two * u + uarray[j][i - 1] + uarray[j][i + 1]) * sx;
161c4762a1bSJed Brown       uyy     = (-two * u + uarray[j - 1][i] + uarray[j + 1][i]) * sy;
162c4762a1bSJed Brown       f[j][i] = uxx + uyy;
163c4762a1bSJed Brown     }
164c4762a1bSJed Brown   }
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   /* Restore vectors */
1679566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, localU, &uarray));
1689566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, F, &f));
1699566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localU));
1709566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(11.0 * ym * xm));
171*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
172c4762a1bSJed Brown }
173c4762a1bSJed Brown 
174c4762a1bSJed Brown /* --------------------------------------------------------------------- */
175c4762a1bSJed Brown /*
176c4762a1bSJed Brown    RHSJacobian - User-provided routine to compute the Jacobian of
177c4762a1bSJed Brown    the nonlinear right-hand-side function of the ODE.
178c4762a1bSJed Brown 
179c4762a1bSJed Brown    Input Parameters:
180c4762a1bSJed Brown    ts - the TS context
181c4762a1bSJed Brown    t - current time
182c4762a1bSJed Brown    U - global input vector
183c4762a1bSJed Brown    dummy - optional user-defined context, as set by TSetRHSJacobian()
184c4762a1bSJed Brown 
185c4762a1bSJed Brown    Output Parameters:
186c4762a1bSJed Brown    J - Jacobian matrix
187c4762a1bSJed Brown    Jpre - optionally different preconditioning matrix
188c4762a1bSJed Brown    str - flag indicating matrix structure
189c4762a1bSJed Brown */
190d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat J, Mat Jpre, void *ctx)
191d71ae5a4SJacob Faibussowitsch {
192c4762a1bSJed Brown   DM            da;
193c4762a1bSJed Brown   DMDALocalInfo info;
194c4762a1bSJed Brown   PetscInt      i, j;
195c4762a1bSJed Brown   PetscReal     hx, hy, sx, sy;
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   PetscFunctionBeginUser;
1989566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
1999566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
2009371c9d4SSatish Balay   hx = 1.0 / (PetscReal)(info.mx - 1);
2019371c9d4SSatish Balay   sx = 1.0 / (hx * hx);
2029371c9d4SSatish Balay   hy = 1.0 / (PetscReal)(info.my - 1);
2039371c9d4SSatish Balay   sy = 1.0 / (hy * hy);
204c4762a1bSJed Brown   for (j = info.ys; j < info.ys + info.ym; j++) {
205c4762a1bSJed Brown     for (i = info.xs; i < info.xs + info.xm; i++) {
206c4762a1bSJed Brown       PetscInt    nc = 0;
207c4762a1bSJed Brown       MatStencil  row, col[5];
208c4762a1bSJed Brown       PetscScalar val[5];
2099371c9d4SSatish Balay       row.i = i;
2109371c9d4SSatish Balay       row.j = j;
211c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info.mx - 1 || j == info.my - 1) {
2129371c9d4SSatish Balay         col[nc].i = i;
2139371c9d4SSatish Balay         col[nc].j = j;
2149371c9d4SSatish Balay         val[nc++] = 1.0;
215c4762a1bSJed Brown       } else {
2169371c9d4SSatish Balay         col[nc].i = i - 1;
2179371c9d4SSatish Balay         col[nc].j = j;
2189371c9d4SSatish Balay         val[nc++] = sx;
2199371c9d4SSatish Balay         col[nc].i = i + 1;
2209371c9d4SSatish Balay         col[nc].j = j;
2219371c9d4SSatish Balay         val[nc++] = sx;
2229371c9d4SSatish Balay         col[nc].i = i;
2239371c9d4SSatish Balay         col[nc].j = j - 1;
2249371c9d4SSatish Balay         val[nc++] = sy;
2259371c9d4SSatish Balay         col[nc].i = i;
2269371c9d4SSatish Balay         col[nc].j = j + 1;
2279371c9d4SSatish Balay         val[nc++] = sy;
2289371c9d4SSatish Balay         col[nc].i = i;
2299371c9d4SSatish Balay         col[nc].j = j;
2309371c9d4SSatish Balay         val[nc++] = -2 * sx - 2 * sy;
231c4762a1bSJed Brown       }
2329566063dSJacob Faibussowitsch       PetscCall(MatSetValuesStencil(Jpre, 1, &row, nc, col, val, INSERT_VALUES));
233c4762a1bSJed Brown     }
234c4762a1bSJed Brown   }
2359566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY));
2369566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY));
237c4762a1bSJed Brown   if (J != Jpre) {
2389566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
2399566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
240c4762a1bSJed Brown   }
241*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
242c4762a1bSJed Brown }
243c4762a1bSJed Brown 
244c4762a1bSJed Brown /* ------------------------------------------------------------------- */
245d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialSolution(DM da, Vec U, void *ptr)
246d71ae5a4SJacob Faibussowitsch {
247c4762a1bSJed Brown   AppCtx       *user = (AppCtx *)ptr;
248c4762a1bSJed Brown   PetscReal     c    = user->c;
249c4762a1bSJed Brown   PetscInt      i, j, xs, ys, xm, ym, Mx, My;
250c4762a1bSJed Brown   PetscScalar **u;
251c4762a1bSJed Brown   PetscReal     hx, hy, x, y, r;
252c4762a1bSJed Brown 
253c4762a1bSJed Brown   PetscFunctionBeginUser;
2549566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
255c4762a1bSJed Brown 
256c4762a1bSJed Brown   hx = 1.0 / (PetscReal)(Mx - 1);
257c4762a1bSJed Brown   hy = 1.0 / (PetscReal)(My - 1);
258c4762a1bSJed Brown 
259c4762a1bSJed Brown   /* Get pointers to vector data */
2609566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, U, &u));
261c4762a1bSJed Brown 
262c4762a1bSJed Brown   /* Get local grid boundaries */
2639566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
264c4762a1bSJed Brown 
265c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
266c4762a1bSJed Brown   for (j = ys; j < ys + ym; j++) {
267c4762a1bSJed Brown     y = j * hy;
268c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
269c4762a1bSJed Brown       x = i * hx;
270c4762a1bSJed Brown       r = PetscSqrtReal((x - .5) * (x - .5) + (y - .5) * (y - .5));
271c4762a1bSJed Brown       if (r < .125) u[j][i] = PetscExpReal(c * r * r * r);
272c4762a1bSJed Brown       else u[j][i] = 0.0;
273c4762a1bSJed Brown     }
274c4762a1bSJed Brown   }
275c4762a1bSJed Brown 
276c4762a1bSJed Brown   /* Restore vectors */
2779566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, U, &u));
278*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
279c4762a1bSJed Brown }
280c4762a1bSJed Brown 
281c4762a1bSJed Brown /*TEST
282c4762a1bSJed Brown 
283c4762a1bSJed Brown     test:
284c4762a1bSJed Brown       args: -ts_max_steps 5 -ts_monitor
285c4762a1bSJed Brown 
286c4762a1bSJed Brown     test:
287c4762a1bSJed Brown       suffix: 2
288c4762a1bSJed Brown       args: -ts_max_steps 5 -ts_monitor
289c4762a1bSJed Brown 
290c4762a1bSJed Brown     test:
291c4762a1bSJed Brown       suffix: 3
292c4762a1bSJed Brown       args: -ts_max_steps 5 -snes_fd_color -ts_monitor
293c4762a1bSJed Brown 
294c4762a1bSJed Brown TEST*/
295