xref: /petsc/src/ts/tests/ex21.c (revision 2a1887a77e7b2c6e00dd0ba96d1387c839460237)
1c4762a1bSJed Brown static char help[] = "Solves a time-dependent nonlinear PDE.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /* ------------------------------------------------------------------------
4c4762a1bSJed Brown 
5c4762a1bSJed Brown    This program solves the two-dimensional time-dependent Bratu problem
6c4762a1bSJed Brown        u_t = u_xx +  u_yy + \lambda*exp(u),
7c4762a1bSJed Brown    on the domain 0 <= x,y <= 1,
8c4762a1bSJed Brown    with the boundary conditions
9c4762a1bSJed Brown        u(t,0,y) = 0, u_x(t,1,y) = 0,
10c4762a1bSJed Brown        u(t,x,0) = 0, u_x(t,x,1) = 0,
11c4762a1bSJed Brown    and the initial condition
12c4762a1bSJed Brown        u(0,x,y) = 0.
13c4762a1bSJed Brown    We discretize the right-hand side using finite differences with
14c4762a1bSJed Brown    uniform grid spacings hx,hy:
15c4762a1bSJed Brown        u_xx = (u_{i+1} - 2u_{i} + u_{i-1})/(hx^2)
16c4762a1bSJed Brown        u_yy = (u_{j+1} - 2u_{j} + u_{j-1})/(hy^2)
17c4762a1bSJed Brown 
18c4762a1bSJed Brown   ------------------------------------------------------------------------- */
19c4762a1bSJed Brown 
20c4762a1bSJed Brown #include <petscdmda.h>
21c4762a1bSJed Brown #include <petscts.h>
22c4762a1bSJed Brown 
23c4762a1bSJed Brown /*
24c4762a1bSJed Brown    User-defined application context - contains data needed by the
25c4762a1bSJed Brown    application-provided call-back routines.
26c4762a1bSJed Brown */
27c4762a1bSJed Brown typedef struct {
28c4762a1bSJed Brown   PetscReal lambda;
29c4762a1bSJed Brown } AppCtx;
30c4762a1bSJed Brown 
31c4762a1bSJed Brown /*
32c4762a1bSJed Brown    FormIFunctionLocal - Evaluates nonlinear implicit function on local process patch
33c4762a1bSJed Brown  */
FormIFunctionLocal(DMDALocalInfo * info,PetscReal t,PetscScalar ** x,PetscScalar ** xdot,PetscScalar ** f,AppCtx * app)34d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormIFunctionLocal(DMDALocalInfo *info, PetscReal t, PetscScalar **x, PetscScalar **xdot, PetscScalar **f, AppCtx *app)
35d71ae5a4SJacob Faibussowitsch {
36c4762a1bSJed Brown   PetscInt    i, j;
37c4762a1bSJed Brown   PetscReal   lambda, hx, hy;
38c4762a1bSJed Brown   PetscScalar ut, u, ue, uw, un, us, uxx, uyy;
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   PetscFunctionBeginUser;
41c4762a1bSJed Brown   lambda = app->lambda;
42c4762a1bSJed Brown   hx     = 1.0 / (PetscReal)(info->mx - 1);
43c4762a1bSJed Brown   hy     = 1.0 / (PetscReal)(info->my - 1);
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   /*
46c4762a1bSJed Brown      Compute RHS function over the locally owned part of the grid
47c4762a1bSJed Brown   */
48c4762a1bSJed Brown   for (j = info->ys; j < info->ys + info->ym; j++) {
49c4762a1bSJed Brown     for (i = info->xs; i < info->xs + info->xm; i++) {
50c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) {
51c4762a1bSJed Brown         /* boundary points */
52c4762a1bSJed Brown         f[j][i] = x[j][i] - (PetscReal)0;
53c4762a1bSJed Brown       } else {
54c4762a1bSJed Brown         /* interior points */
55c4762a1bSJed Brown         ut = xdot[j][i];
56c4762a1bSJed Brown         u  = x[j][i];
57c4762a1bSJed Brown         uw = x[j][i - 1];
58c4762a1bSJed Brown         ue = x[j][i + 1];
59c4762a1bSJed Brown         un = x[j + 1][i];
60c4762a1bSJed Brown         us = x[j - 1][i];
61c4762a1bSJed Brown 
62c4762a1bSJed Brown         uxx     = (uw - 2.0 * u + ue) / (hx * hx);
63c4762a1bSJed Brown         uyy     = (un - 2.0 * u + us) / (hy * hy);
64c4762a1bSJed Brown         f[j][i] = ut - uxx - uyy - lambda * PetscExpScalar(u);
65c4762a1bSJed Brown       }
66c4762a1bSJed Brown     }
67c4762a1bSJed Brown   }
683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
69c4762a1bSJed Brown }
70c4762a1bSJed Brown 
71c4762a1bSJed Brown /*
72c4762a1bSJed Brown    FormIJacobianLocal - Evaluates implicit Jacobian matrix on local process patch
73c4762a1bSJed Brown */
FormIJacobianLocal(DMDALocalInfo * info,PetscReal t,PetscScalar ** x,PetscScalar ** xdot,PetscScalar shift,Mat jac,Mat jacpre,AppCtx * app)74d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormIJacobianLocal(DMDALocalInfo *info, PetscReal t, PetscScalar **x, PetscScalar **xdot, PetscScalar shift, Mat jac, Mat jacpre, AppCtx *app)
75d71ae5a4SJacob Faibussowitsch {
76c4762a1bSJed Brown   PetscInt    i, j, k;
77c4762a1bSJed Brown   MatStencil  col[5], row;
78c4762a1bSJed Brown   PetscScalar v[5], lambda, hx, hy;
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   PetscFunctionBeginUser;
81c4762a1bSJed Brown   lambda = app->lambda;
82c4762a1bSJed Brown   hx     = 1.0 / (PetscReal)(info->mx - 1);
83c4762a1bSJed Brown   hy     = 1.0 / (PetscReal)(info->my - 1);
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   /*
86c4762a1bSJed Brown      Compute Jacobian entries for the locally owned part of the grid
87c4762a1bSJed Brown   */
88c4762a1bSJed Brown   for (j = info->ys; j < info->ys + info->ym; j++) {
89c4762a1bSJed Brown     for (i = info->xs; i < info->xs + info->xm; i++) {
909371c9d4SSatish Balay       row.j = j;
919371c9d4SSatish Balay       row.i = i;
929371c9d4SSatish Balay       k     = 0;
93c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) {
94c4762a1bSJed Brown         /* boundary points */
95c4762a1bSJed Brown         v[0] = 1.0;
969566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(jacpre, 1, &row, 1, &row, v, INSERT_VALUES));
97c4762a1bSJed Brown       } else {
98c4762a1bSJed Brown         /* interior points */
999371c9d4SSatish Balay         v[k]     = -1.0 / (hy * hy);
1009371c9d4SSatish Balay         col[k].j = j - 1;
1019371c9d4SSatish Balay         col[k].i = i;
1029371c9d4SSatish Balay         k++;
1039371c9d4SSatish Balay         v[k]     = -1.0 / (hx * hx);
1049371c9d4SSatish Balay         col[k].j = j;
1059371c9d4SSatish Balay         col[k].i = i - 1;
1069371c9d4SSatish Balay         k++;
107c4762a1bSJed Brown 
108c4762a1bSJed Brown         v[k]     = shift + 2.0 / (hx * hx) + 2.0 / (hy * hy) - lambda * PetscExpScalar(x[j][i]);
1099371c9d4SSatish Balay         col[k].j = j;
1109371c9d4SSatish Balay         col[k].i = i;
1119371c9d4SSatish Balay         k++;
112c4762a1bSJed Brown 
1139371c9d4SSatish Balay         v[k]     = -1.0 / (hx * hx);
1149371c9d4SSatish Balay         col[k].j = j;
1159371c9d4SSatish Balay         col[k].i = i + 1;
1169371c9d4SSatish Balay         k++;
1179371c9d4SSatish Balay         v[k]     = -1.0 / (hy * hy);
1189371c9d4SSatish Balay         col[k].j = j + 1;
1199371c9d4SSatish Balay         col[k].i = i;
1209371c9d4SSatish Balay         k++;
121c4762a1bSJed Brown 
1229566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(jacpre, 1, &row, k, col, v, INSERT_VALUES));
123c4762a1bSJed Brown       }
124c4762a1bSJed Brown     }
125c4762a1bSJed Brown   }
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   /*
128c4762a1bSJed Brown      Assemble matrix
129c4762a1bSJed Brown   */
1309566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jacpre, MAT_FINAL_ASSEMBLY));
1319566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jacpre, MAT_FINAL_ASSEMBLY));
1323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
133c4762a1bSJed Brown }
134c4762a1bSJed Brown 
main(int argc,char ** argv)135d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
136d71ae5a4SJacob Faibussowitsch {
137c4762a1bSJed Brown   TS              ts; /* ODE integrator */
138c4762a1bSJed Brown   DM              da; /* DM context */
139c4762a1bSJed Brown   Vec             U;  /* solution vector */
140c4762a1bSJed Brown   DMBoundaryType  bt = DM_BOUNDARY_NONE;
141c4762a1bSJed Brown   DMDAStencilType st = DMDA_STENCIL_STAR;
142c4762a1bSJed Brown   PetscInt        sw = 1;
143c4762a1bSJed Brown   PetscInt        N  = 17;
144c4762a1bSJed Brown   PetscInt        n  = PETSC_DECIDE;
145c4762a1bSJed Brown   AppCtx          app;
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148c4762a1bSJed Brown      Initialize program
149c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
150327415f7SBarry Smith   PetscFunctionBeginUser;
1519566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
152d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex21 options", "");
153c4762a1bSJed Brown   {
1549371c9d4SSatish Balay     app.lambda = 6.8;
1559371c9d4SSatish Balay     app.lambda = 6.0;
1569566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-lambda", "", "", app.lambda, &app.lambda, NULL));
157c4762a1bSJed Brown   }
158d0609cedSBarry Smith   PetscOptionsEnd();
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161c4762a1bSJed Brown      Create DM context
162c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1639566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, bt, bt, st, N, N, n, n, 1, sw, NULL, NULL, &da));
1649566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
1659566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
1669566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0, 1.0));
167c4762a1bSJed Brown 
168c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
169c4762a1bSJed Brown      Create timestepping solver context
170c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1719566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1729566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
1739566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, da));
1749566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
175c4762a1bSJed Brown 
1769566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
1778434afd1SBarry Smith   PetscCall(DMDATSSetIFunctionLocal(da, INSERT_VALUES, (DMDATSIFunctionLocalFn *)FormIFunctionLocal, &app));
1788434afd1SBarry Smith   PetscCall(DMDATSSetIJacobianLocal(da, (DMDATSIJacobianLocalFn *)FormIJacobianLocal, &app));
179c4762a1bSJed Brown 
180c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
181c4762a1bSJed Brown      Set solver options
182c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1839566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSBDF));
1849566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 1e-4));
1859566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 1.0));
1869566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1879566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
190c4762a1bSJed Brown      Set initial conditions
191c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1929566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
1939566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &U));
1949566063dSJacob Faibussowitsch   PetscCall(VecSet(U, 0.0));
1959566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, U));
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198c4762a1bSJed Brown      Run timestepping solver
199c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2009566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, U));
201c4762a1bSJed Brown 
202c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
203c4762a1bSJed Brown       All PETSc objects should be destroyed when they are no longer needed.
204c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
2069566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
2079566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
208b122ec5aSJacob Faibussowitsch   return 0;
209c4762a1bSJed Brown }
210c4762a1bSJed Brown 
211c4762a1bSJed Brown /*TEST
212c4762a1bSJed Brown 
213c4762a1bSJed Brown     testset:
214c4762a1bSJed Brown       requires: !single !complex
215*188af4bfSBarry Smith       args: -da_grid_x 5 -da_grid_y 5 -da_refine 2 -dm_view -ts_type bdf -ts_adapt_type none -ts_time_step 1e-3 -ts_monitor -ts_max_steps 5 -ts_view -snes_rtol 1e-6 -snes_type ngmres -npc_snes_type fas
216c4762a1bSJed Brown       filter: grep -v "total number of"
217c4762a1bSJed Brown       test:
218c4762a1bSJed Brown         suffix: 1_bdf_ngmres_fas_ms
219c4762a1bSJed Brown         args: -prefix_push npc_fas_levels_ -snes_type ms -snes_max_it 5 -ksp_type preonly -prefix_pop
220c4762a1bSJed Brown       test:
221c4762a1bSJed Brown         suffix: 2_bdf_ngmres_fas_ms
222c4762a1bSJed Brown         args: -prefix_push npc_fas_levels_ -snes_type ms -snes_max_it 5 -ksp_type preonly -prefix_pop
223c4762a1bSJed Brown         nsize: 2
224c4762a1bSJed Brown       test:
225c4762a1bSJed Brown         suffix: 1_bdf_ngmres_fas_ngs
226c4762a1bSJed Brown         args: -prefix_push npc_fas_levels_ -snes_type ngs -snes_max_it 5 -prefix_pop
227c4762a1bSJed Brown       test:
228c4762a1bSJed Brown         suffix: 2_bdf_ngmres_fas_ngs
229c4762a1bSJed Brown         args: -prefix_push npc_fas_levels_ -snes_type ngs -snes_max_it 5 -prefix_pop
230c4762a1bSJed Brown         nsize: 2
231c4762a1bSJed Brown 
232c4762a1bSJed Brown TEST*/
233