xref: /petsc/src/ts/tests/ex21.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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  */
34c4762a1bSJed Brown static PetscErrorCode FormIFunctionLocal(DMDALocalInfo *info,PetscReal t,PetscScalar **x,PetscScalar **xdot,PetscScalar **f,AppCtx *app)
35c4762a1bSJed Brown {
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   }
68c4762a1bSJed Brown   PetscFunctionReturn(0);
69c4762a1bSJed Brown }
70c4762a1bSJed Brown 
71c4762a1bSJed Brown /*
72c4762a1bSJed Brown    FormIJacobianLocal - Evaluates implicit Jacobian matrix on local process patch
73c4762a1bSJed Brown */
74c4762a1bSJed Brown static PetscErrorCode FormIJacobianLocal(DMDALocalInfo *info,PetscReal t,PetscScalar **x,PetscScalar **xdot,PetscScalar shift,Mat jac,Mat jacpre,AppCtx *app)
75c4762a1bSJed Brown {
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++) {
90c4762a1bSJed Brown       row.j = j; row.i = i; k = 0;
91c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
92c4762a1bSJed Brown         /* boundary points */
93c4762a1bSJed Brown         v[0] = 1.0;
945f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValuesStencil(jacpre,1,&row,1,&row,v,INSERT_VALUES));
95c4762a1bSJed Brown       } else {
96c4762a1bSJed Brown         /* interior points */
97c4762a1bSJed Brown         v[k] = -1.0/(hy*hy); col[k].j = j-1; col[k].i = i;   k++;
98c4762a1bSJed Brown         v[k] = -1.0/(hx*hx); col[k].j = j;   col[k].i = i-1; k++;
99c4762a1bSJed Brown 
100c4762a1bSJed Brown         v[k] = shift + 2.0/(hx*hx) + 2.0/(hy*hy) - lambda*PetscExpScalar(x[j][i]);
101c4762a1bSJed Brown         col[k].j = j; col[k].i = i; k++;
102c4762a1bSJed Brown 
103c4762a1bSJed Brown         v[k] = -1.0/(hx*hx); col[k].j = j;   col[k].i = i+1; k++;
104c4762a1bSJed Brown         v[k] = -1.0/(hy*hy); col[k].j = j+1; col[k].i = i;   k++;
105c4762a1bSJed Brown 
1065f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValuesStencil(jacpre,1,&row,k,col,v,INSERT_VALUES));
107c4762a1bSJed Brown       }
108c4762a1bSJed Brown     }
109c4762a1bSJed Brown   }
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   /*
112c4762a1bSJed Brown      Assemble matrix
113c4762a1bSJed Brown   */
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(jacpre,MAT_FINAL_ASSEMBLY));
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(jacpre,MAT_FINAL_ASSEMBLY));
116c4762a1bSJed Brown   PetscFunctionReturn(0);
117c4762a1bSJed Brown }
118c4762a1bSJed Brown 
119c4762a1bSJed Brown int main(int argc,char **argv)
120c4762a1bSJed Brown {
121c4762a1bSJed Brown   TS              ts;            /* ODE integrator */
122c4762a1bSJed Brown   DM              da;            /* DM context */
123c4762a1bSJed Brown   Vec             U;             /* solution vector */
124c4762a1bSJed Brown   DMBoundaryType  bt = DM_BOUNDARY_NONE;
125c4762a1bSJed Brown   DMDAStencilType st = DMDA_STENCIL_STAR;
126c4762a1bSJed Brown   PetscInt        sw = 1;
127c4762a1bSJed Brown   PetscInt        N  = 17;
128c4762a1bSJed Brown   PetscInt        n  = PETSC_DECIDE;
129c4762a1bSJed Brown   AppCtx          app;
130c4762a1bSJed Brown   PetscErrorCode  ierr;
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133c4762a1bSJed Brown      Initialize program
134c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
135*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,NULL,help));
136c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex21 options","");CHKERRQ(ierr);
137c4762a1bSJed Brown   {
138c4762a1bSJed Brown     app.lambda = 6.8; app.lambda = 6.0;
1395f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsReal("-lambda","","",app.lambda,&app.lambda,NULL));
140c4762a1bSJed Brown   }
141c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144c4762a1bSJed Brown      Create DM context
145c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,bt,bt,st,N,N,n,n,1,sw,NULL,NULL,&da));
1475f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
1485f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
1495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0,1.0));
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152c4762a1bSJed Brown      Create timestepping solver context
153c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1545f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR));
1565f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetDM(ts,da));
1575f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
158c4762a1bSJed Brown 
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
1605f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDATSSetIFunctionLocal(da,INSERT_VALUES,(DMDATSIFunctionLocal)FormIFunctionLocal,&app));
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDATSSetIJacobianLocal(da,(DMDATSIJacobianLocal)FormIJacobianLocal,&app));
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
164c4762a1bSJed Brown      Set solver options
165c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1665f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSBDF));
1675f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,1e-4));
1685f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,1.0));
1695f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
1705f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
173c4762a1bSJed Brown      Set initial conditions
174c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
1765f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&U));
1775f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(U,0.0));
1785f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,U));
179c4762a1bSJed Brown 
180c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
181c4762a1bSJed Brown      Run timestepping solver
182c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1835f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,U));
184c4762a1bSJed Brown 
185c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
186c4762a1bSJed Brown       All PETSc objects should be destroyed when they are no longer needed.
187c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1885f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&U));
1895f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
190*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
191*b122ec5aSJacob Faibussowitsch   return 0;
192c4762a1bSJed Brown }
193c4762a1bSJed Brown 
194c4762a1bSJed Brown /*TEST
195c4762a1bSJed Brown 
196c4762a1bSJed Brown     testset:
197c4762a1bSJed Brown       requires: !single !complex
198c4762a1bSJed Brown       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
199c4762a1bSJed Brown       filter: grep -v "total number of"
200c4762a1bSJed Brown       test:
201c4762a1bSJed Brown         suffix: 1_bdf_ngmres_fas_ms
202c4762a1bSJed Brown         args: -prefix_push npc_fas_levels_ -snes_type ms -snes_max_it 5 -ksp_type preonly -prefix_pop
203c4762a1bSJed Brown       test:
204c4762a1bSJed Brown         suffix: 2_bdf_ngmres_fas_ms
205c4762a1bSJed Brown         args: -prefix_push npc_fas_levels_ -snes_type ms -snes_max_it 5 -ksp_type preonly -prefix_pop
206c4762a1bSJed Brown         nsize: 2
207c4762a1bSJed Brown       test:
208c4762a1bSJed Brown         suffix: 1_bdf_ngmres_fas_ngs
209c4762a1bSJed Brown         args: -prefix_push npc_fas_levels_ -snes_type ngs -snes_max_it 5 -prefix_pop
210c4762a1bSJed Brown       test:
211c4762a1bSJed Brown         suffix: 2_bdf_ngmres_fas_ngs
212c4762a1bSJed Brown         args: -prefix_push npc_fas_levels_ -snes_type ngs -snes_max_it 5 -prefix_pop
213c4762a1bSJed Brown         nsize: 2
214c4762a1bSJed Brown 
215c4762a1bSJed Brown TEST*/
216