xref: /petsc/src/ts/tests/ex21.c (revision 8fb5bd83c3955fefcf33a54e3bb66920a9fa884b)
1 static char help[] ="Solves a time-dependent nonlinear PDE.\n";
2 
3 /* ------------------------------------------------------------------------
4 
5    This program solves the two-dimensional time-dependent Bratu problem
6        u_t = u_xx +  u_yy + \lambda*exp(u),
7    on the domain 0 <= x,y <= 1,
8    with the boundary conditions
9        u(t,0,y) = 0, u_x(t,1,y) = 0,
10        u(t,x,0) = 0, u_x(t,x,1) = 0,
11    and the initial condition
12        u(0,x,y) = 0.
13    We discretize the right-hand side using finite differences with
14    uniform grid spacings hx,hy:
15        u_xx = (u_{i+1} - 2u_{i} + u_{i-1})/(hx^2)
16        u_yy = (u_{j+1} - 2u_{j} + u_{j-1})/(hy^2)
17 
18   ------------------------------------------------------------------------- */
19 
20 #include <petscdmda.h>
21 #include <petscts.h>
22 
23 /*
24    User-defined application context - contains data needed by the
25    application-provided call-back routines.
26 */
27 typedef struct {
28   PetscReal lambda;
29 } AppCtx;
30 
31 /*
32    FormIFunctionLocal - Evaluates nonlinear implicit function on local process patch
33  */
34 static PetscErrorCode FormIFunctionLocal(DMDALocalInfo *info,PetscReal t,PetscScalar **x,PetscScalar **xdot,PetscScalar **f,AppCtx *app)
35 {
36   PetscInt       i,j;
37   PetscReal      lambda,hx,hy;
38   PetscScalar    ut,u,ue,uw,un,us,uxx,uyy;
39 
40   PetscFunctionBeginUser;
41   lambda = app->lambda;
42   hx     = 1.0/(PetscReal)(info->mx-1);
43   hy     = 1.0/(PetscReal)(info->my-1);
44 
45   /*
46      Compute RHS function over the locally owned part of the grid
47   */
48   for (j=info->ys; j<info->ys+info->ym; j++) {
49     for (i=info->xs; i<info->xs+info->xm; i++) {
50       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
51         /* boundary points */
52         f[j][i] = x[j][i] - (PetscReal)0;
53       } else {
54         /* interior points */
55         ut = xdot[j][i];
56         u  = x[j][i];
57         uw = x[j][i-1];
58         ue = x[j][i+1];
59         un = x[j+1][i];
60         us = x[j-1][i];
61 
62         uxx = (uw - 2.0*u + ue)/(hx*hx);
63         uyy = (un - 2.0*u + us)/(hy*hy);
64         f[j][i] = ut - uxx - uyy - lambda*PetscExpScalar(u);
65       }
66     }
67   }
68   PetscFunctionReturn(0);
69 }
70 
71 /*
72    FormIJacobianLocal - Evaluates implicit Jacobian matrix on local process patch
73 */
74 static PetscErrorCode FormIJacobianLocal(DMDALocalInfo *info,PetscReal t,PetscScalar **x,PetscScalar **xdot,PetscScalar shift,Mat jac,Mat jacpre,AppCtx *app)
75 {
76   PetscInt       i,j,k;
77   MatStencil     col[5],row;
78   PetscScalar    v[5],lambda,hx,hy;
79 
80   PetscFunctionBeginUser;
81   lambda = app->lambda;
82   hx     = 1.0/(PetscReal)(info->mx-1);
83   hy     = 1.0/(PetscReal)(info->my-1);
84 
85   /*
86      Compute Jacobian entries for the locally owned part of the grid
87   */
88   for (j=info->ys; j<info->ys+info->ym; j++) {
89     for (i=info->xs; i<info->xs+info->xm; i++) {
90       row.j = j; row.i = i; k = 0;
91       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
92         /* boundary points */
93         v[0] = 1.0;
94         PetscCall(MatSetValuesStencil(jacpre,1,&row,1,&row,v,INSERT_VALUES));
95       } else {
96         /* interior points */
97         v[k] = -1.0/(hy*hy); col[k].j = j-1; col[k].i = i;   k++;
98         v[k] = -1.0/(hx*hx); col[k].j = j;   col[k].i = i-1; k++;
99 
100         v[k] = shift + 2.0/(hx*hx) + 2.0/(hy*hy) - lambda*PetscExpScalar(x[j][i]);
101         col[k].j = j; col[k].i = i; k++;
102 
103         v[k] = -1.0/(hx*hx); col[k].j = j;   col[k].i = i+1; k++;
104         v[k] = -1.0/(hy*hy); col[k].j = j+1; col[k].i = i;   k++;
105 
106         PetscCall(MatSetValuesStencil(jacpre,1,&row,k,col,v,INSERT_VALUES));
107       }
108     }
109   }
110 
111   /*
112      Assemble matrix
113   */
114   PetscCall(MatAssemblyBegin(jacpre,MAT_FINAL_ASSEMBLY));
115   PetscCall(MatAssemblyEnd(jacpre,MAT_FINAL_ASSEMBLY));
116   PetscFunctionReturn(0);
117 }
118 
119 int main(int argc,char **argv)
120 {
121   TS              ts;            /* ODE integrator */
122   DM              da;            /* DM context */
123   Vec             U;             /* solution vector */
124   DMBoundaryType  bt = DM_BOUNDARY_NONE;
125   DMDAStencilType st = DMDA_STENCIL_STAR;
126   PetscInt        sw = 1;
127   PetscInt        N  = 17;
128   PetscInt        n  = PETSC_DECIDE;
129   AppCtx          app;
130 
131   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
132      Initialize program
133      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
134   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
135   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex21 options","");
136   {
137     app.lambda = 6.8; app.lambda = 6.0;
138     PetscCall(PetscOptionsReal("-lambda","","",app.lambda,&app.lambda,NULL));
139   }
140   PetscOptionsEnd();
141 
142   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143      Create DM context
144      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145   PetscCall(DMDACreate2d(PETSC_COMM_WORLD,bt,bt,st,N,N,n,n,1,sw,NULL,NULL,&da));
146   PetscCall(DMSetFromOptions(da));
147   PetscCall(DMSetUp(da));
148   PetscCall(DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0,1.0));
149 
150   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151      Create timestepping solver context
152      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
154   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
155   PetscCall(TSSetDM(ts,da));
156   PetscCall(DMDestroy(&da));
157 
158   PetscCall(TSGetDM(ts,&da));
159   PetscCall(DMDATSSetIFunctionLocal(da,INSERT_VALUES,(DMDATSIFunctionLocal)FormIFunctionLocal,&app));
160   PetscCall(DMDATSSetIJacobianLocal(da,(DMDATSIJacobianLocal)FormIJacobianLocal,&app));
161 
162   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
163      Set solver options
164    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
165   PetscCall(TSSetType(ts,TSBDF));
166   PetscCall(TSSetTimeStep(ts,1e-4));
167   PetscCall(TSSetMaxTime(ts,1.0));
168   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
169   PetscCall(TSSetFromOptions(ts));
170 
171   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
172      Set initial conditions
173    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
174   PetscCall(TSGetDM(ts,&da));
175   PetscCall(DMCreateGlobalVector(da,&U));
176   PetscCall(VecSet(U,0.0));
177   PetscCall(TSSetSolution(ts,U));
178 
179   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
180      Run timestepping solver
181      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
182   PetscCall(TSSolve(ts,U));
183 
184   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185       All PETSc objects should be destroyed when they are no longer needed.
186    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
187   PetscCall(VecDestroy(&U));
188   PetscCall(TSDestroy(&ts));
189   PetscCall(PetscFinalize());
190   return 0;
191 }
192 
193 /*TEST
194 
195     testset:
196       requires: !single !complex
197       args: -da_grid_x 5 -da_grid_y 5 -da_refine 2 -dm_view -ts_type bdf -ts_adapt_type none -ts_dt 1e-3 -ts_monitor -ts_max_steps 5 -ts_view -snes_rtol 1e-6 -snes_type ngmres -npc_snes_type fas
198       filter: grep -v "total number of"
199       test:
200         suffix: 1_bdf_ngmres_fas_ms
201         args: -prefix_push npc_fas_levels_ -snes_type ms -snes_max_it 5 -ksp_type preonly -prefix_pop
202       test:
203         suffix: 2_bdf_ngmres_fas_ms
204         args: -prefix_push npc_fas_levels_ -snes_type ms -snes_max_it 5 -ksp_type preonly -prefix_pop
205         nsize: 2
206       test:
207         suffix: 1_bdf_ngmres_fas_ngs
208         args: -prefix_push npc_fas_levels_ -snes_type ngs -snes_max_it 5 -prefix_pop
209       test:
210         suffix: 2_bdf_ngmres_fas_ngs
211         args: -prefix_push npc_fas_levels_ -snes_type ngs -snes_max_it 5 -prefix_pop
212         nsize: 2
213 
214 TEST*/
215