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