xref: /petsc/src/ts/tutorials/ex74.c (revision d2567f346889bc69f462cb815202b652ec2f297c)
1*d2567f34SHong Zhang static char help[] = "Solves the constant-coefficient 1D heat equation \n\
2*d2567f34SHong Zhang with an Implicit Runge-Kutta method using MatKAIJ.                  \n\
3*d2567f34SHong Zhang                                                                     \n\
4*d2567f34SHong Zhang     du      d^2 u                                                   \n\
5*d2567f34SHong Zhang     --  = a ----- ; 0 <= x <= 1;                                    \n\
6*d2567f34SHong Zhang     dt      dx^2                                                    \n\
7*d2567f34SHong Zhang                                                                     \n\
8*d2567f34SHong Zhang   with periodic boundary conditions                                 \n\
9*d2567f34SHong Zhang                                                                     \n\
10*d2567f34SHong Zhang 2nd order central discretization in space:                          \n\
11*d2567f34SHong Zhang                                                                     \n\
12*d2567f34SHong Zhang    [ d^2 u ]     u_{i+1} - 2u_i + u_{i-1}                           \n\
13*d2567f34SHong Zhang    [ ----- ]  =  ------------------------                           \n\
14*d2567f34SHong Zhang    [ dx^2  ]i              h^2                                      \n\
15*d2567f34SHong Zhang                                                                     \n\
16*d2567f34SHong Zhang     i = grid index;    h = x_{i+1}-x_i (Uniform)                    \n\
17*d2567f34SHong Zhang     0 <= i < n         h = 1.0/n                                    \n\
18*d2567f34SHong Zhang                                                                     \n\
19*d2567f34SHong Zhang Thus,                                                               \n\
20*d2567f34SHong Zhang                                                                     \n\
21*d2567f34SHong Zhang    du                                                               \n\
22*d2567f34SHong Zhang    --  = Ju;  J = (a/h^2) tridiagonal(1,-2,1)_n                     \n\
23*d2567f34SHong Zhang    dt                                                               \n\
24*d2567f34SHong Zhang                                                                     \n\
25*d2567f34SHong Zhang This example is a TS version of the KSP ex74.c tutorial.            \n";
26*d2567f34SHong Zhang 
27*d2567f34SHong Zhang #include <petscts.h>
28*d2567f34SHong Zhang 
29*d2567f34SHong Zhang typedef enum {
30*d2567f34SHong Zhang   PHYSICS_DIFFUSION,
31*d2567f34SHong Zhang   PHYSICS_ADVECTION
32*d2567f34SHong Zhang } PhysicsType;
33*d2567f34SHong Zhang const char *const PhysicsTypes[] = {"DIFFUSION","ADVECTION","PhysicsType","PHYSICS_",NULL};
34*d2567f34SHong Zhang 
35*d2567f34SHong Zhang typedef struct Context {
36*d2567f34SHong Zhang   PetscReal     a;              /* diffusion coefficient      */
37*d2567f34SHong Zhang   PetscReal     xmin,xmax;      /* domain bounds              */
38*d2567f34SHong Zhang   PetscInt      imax;           /* number of grid points      */
39*d2567f34SHong Zhang   PhysicsType   physics_type;
40*d2567f34SHong Zhang } UserContext;
41*d2567f34SHong Zhang 
42*d2567f34SHong Zhang static PetscErrorCode ExactSolution(Vec,void*,PetscReal);
43*d2567f34SHong Zhang static PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*);
44*d2567f34SHong Zhang 
45*d2567f34SHong Zhang int main(int argc, char **argv)
46*d2567f34SHong Zhang {
47*d2567f34SHong Zhang   TS             ts;
48*d2567f34SHong Zhang   Mat            A;
49*d2567f34SHong Zhang   Vec            u,uex;
50*d2567f34SHong Zhang   UserContext    ctxt;
51*d2567f34SHong Zhang   PetscReal      err,ftime;
52*d2567f34SHong Zhang   PetscErrorCode ierr;
53*d2567f34SHong Zhang 
54*d2567f34SHong Zhang   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
55*d2567f34SHong Zhang 
56*d2567f34SHong Zhang   /* default value */
57*d2567f34SHong Zhang   ctxt.a       = 0.1;
58*d2567f34SHong Zhang   ctxt.xmin    = 0.0;
59*d2567f34SHong Zhang   ctxt.xmax    = 1.0;
60*d2567f34SHong Zhang   ctxt.imax    = 40;
61*d2567f34SHong Zhang   ctxt.physics_type = PHYSICS_DIFFUSION;
62*d2567f34SHong Zhang 
63*d2567f34SHong Zhang   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"IRK options","");CHKERRQ(ierr);
64*d2567f34SHong Zhang   ierr = PetscOptionsReal("-a","diffusion coefficient","<1.0>",ctxt.a,&ctxt.a,NULL);CHKERRQ(ierr);
65*d2567f34SHong Zhang   ierr = PetscOptionsInt ("-imax","grid size","<20>",ctxt.imax,&ctxt.imax,NULL);CHKERRQ(ierr);
66*d2567f34SHong Zhang   ierr = PetscOptionsReal("-xmin","xmin","<0.0>",ctxt.xmin,&ctxt.xmin,NULL);CHKERRQ(ierr);
67*d2567f34SHong Zhang   ierr = PetscOptionsReal("-xmax","xmax","<1.0>",ctxt.xmax,&ctxt.xmax,NULL);CHKERRQ(ierr);
68*d2567f34SHong Zhang   ierr = PetscOptionsEnum("-physics_type","Type of process to discretize","",PhysicsTypes,(PetscEnum)ctxt.physics_type,(PetscEnum*)&ctxt.physics_type,NULL);CHKERRQ(ierr);
69*d2567f34SHong Zhang   ierr = PetscOptionsEnd();CHKERRQ(ierr);
70*d2567f34SHong Zhang 
71*d2567f34SHong Zhang   /* allocate and initialize solution vector and exact solution */
72*d2567f34SHong Zhang   ierr = VecCreate(PETSC_COMM_WORLD,&u);CHKERRQ(ierr);
73*d2567f34SHong Zhang   ierr = VecSetSizes(u,PETSC_DECIDE,ctxt.imax);CHKERRQ(ierr);
74*d2567f34SHong Zhang   ierr = VecSetFromOptions(u);CHKERRQ(ierr);
75*d2567f34SHong Zhang   ierr = VecDuplicate(u,&uex);CHKERRQ(ierr);
76*d2567f34SHong Zhang   /* initial solution */
77*d2567f34SHong Zhang   ierr = ExactSolution(u,&ctxt,0.0);CHKERRQ(ierr);
78*d2567f34SHong Zhang 
79*d2567f34SHong Zhang   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
80*d2567f34SHong Zhang   ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr);
81*d2567f34SHong Zhang   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,ctxt.imax,ctxt.imax);CHKERRQ(ierr);
82*d2567f34SHong Zhang   ierr = MatSetUp(A);CHKERRQ(ierr);
83*d2567f34SHong Zhang 
84*d2567f34SHong Zhang   /* Create and set options for TS */
85*d2567f34SHong Zhang   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
86*d2567f34SHong Zhang   ierr = TSSetProblemType(ts,TS_LINEAR);CHKERRQ(ierr);
87*d2567f34SHong Zhang   ierr = TSSetTimeStep(ts,0.125);CHKERRQ(ierr);
88*d2567f34SHong Zhang   ierr = TSSetSolution(ts,u);CHKERRQ(ierr);
89*d2567f34SHong Zhang   ierr = TSSetMaxSteps(ts,10);CHKERRQ(ierr);
90*d2567f34SHong Zhang   ierr = TSSetMaxTime(ts,1.0);CHKERRQ(ierr);
91*d2567f34SHong Zhang   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
92*d2567f34SHong Zhang   ierr = TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&ctxt);CHKERRQ(ierr);
93*d2567f34SHong Zhang   ierr = RHSJacobian(ts,0,u,A,A,&ctxt);CHKERRQ(ierr);
94*d2567f34SHong Zhang   ierr = TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,&ctxt);CHKERRQ(ierr);
95*d2567f34SHong Zhang   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
96*d2567f34SHong Zhang   ierr = TSSolve(ts,u);CHKERRQ(ierr);
97*d2567f34SHong Zhang 
98*d2567f34SHong Zhang   ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
99*d2567f34SHong Zhang   /* exact   solution */
100*d2567f34SHong Zhang   ierr = ExactSolution(uex,&ctxt,ftime);CHKERRQ(ierr);
101*d2567f34SHong Zhang 
102*d2567f34SHong Zhang   /* Calculate error in final solution */
103*d2567f34SHong Zhang   ierr = VecAYPX(uex,-1.0,u);CHKERRQ(ierr);
104*d2567f34SHong Zhang   ierr = VecNorm(uex,NORM_2,&err);CHKERRQ(ierr);
105*d2567f34SHong Zhang   err  = PetscSqrtReal(err*err/((PetscReal)ctxt.imax));
106*d2567f34SHong Zhang   ierr = PetscPrintf(PETSC_COMM_WORLD,"L2 norm of the numerical error = %g (time=%g)\n",(double)err,(double)ftime);CHKERRQ(ierr);
107*d2567f34SHong Zhang 
108*d2567f34SHong Zhang   /* Free up memory */
109*d2567f34SHong Zhang   ierr = TSDestroy(&ts);CHKERRQ(ierr);
110*d2567f34SHong Zhang   ierr = MatDestroy(&A);CHKERRQ(ierr);
111*d2567f34SHong Zhang   ierr = VecDestroy(&uex);CHKERRQ(ierr);
112*d2567f34SHong Zhang   ierr = VecDestroy(&u);CHKERRQ(ierr);
113*d2567f34SHong Zhang   ierr = PetscFinalize();
114*d2567f34SHong Zhang   return ierr;
115*d2567f34SHong Zhang }
116*d2567f34SHong Zhang 
117*d2567f34SHong Zhang PetscErrorCode ExactSolution(Vec u,void *c,PetscReal t)
118*d2567f34SHong Zhang {
119*d2567f34SHong Zhang   UserContext     *ctxt = (UserContext*) c;
120*d2567f34SHong Zhang   PetscErrorCode  ierr;
121*d2567f34SHong Zhang   PetscInt        i,is,ie;
122*d2567f34SHong Zhang   PetscScalar     *uarr;
123*d2567f34SHong Zhang   PetscReal       x,dx,a=ctxt->a,pi=PETSC_PI;
124*d2567f34SHong Zhang 
125*d2567f34SHong Zhang   PetscFunctionBegin;
126*d2567f34SHong Zhang   dx = (ctxt->xmax - ctxt->xmin)/((PetscReal) ctxt->imax);
127*d2567f34SHong Zhang   ierr = VecGetOwnershipRange(u,&is,&ie);CHKERRQ(ierr);
128*d2567f34SHong Zhang   ierr = VecGetArray(u,&uarr);CHKERRQ(ierr);
129*d2567f34SHong Zhang   for (i=is; i<ie; i++) {
130*d2567f34SHong Zhang     x          = i * dx;
131*d2567f34SHong Zhang     switch (ctxt->physics_type) {
132*d2567f34SHong Zhang     case PHYSICS_DIFFUSION:
133*d2567f34SHong Zhang       uarr[i-is] = PetscExpScalar(-4.0*pi*pi*a*t)*PetscSinScalar(2*pi*x);
134*d2567f34SHong Zhang       break;
135*d2567f34SHong Zhang     case PHYSICS_ADVECTION:
136*d2567f34SHong Zhang       uarr[i-is] = PetscSinScalar(2*pi*(x - a*t));
137*d2567f34SHong Zhang       break;
138*d2567f34SHong Zhang     default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for physics type %s",PhysicsTypes[ctxt->physics_type]);
139*d2567f34SHong Zhang     }
140*d2567f34SHong Zhang   }
141*d2567f34SHong Zhang   ierr = VecRestoreArray(u,&uarr);CHKERRQ(ierr);
142*d2567f34SHong Zhang   PetscFunctionReturn(0);
143*d2567f34SHong Zhang }
144*d2567f34SHong Zhang 
145*d2567f34SHong Zhang static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat J,Mat Jpre,void *ctx)
146*d2567f34SHong Zhang {
147*d2567f34SHong Zhang   UserContext    *user = (UserContext*) ctx;
148*d2567f34SHong Zhang   PetscInt       matis,matie,i;
149*d2567f34SHong Zhang   PetscReal      dx,dx2;
150*d2567f34SHong Zhang   PetscErrorCode ierr;
151*d2567f34SHong Zhang 
152*d2567f34SHong Zhang   PetscFunctionBegin;
153*d2567f34SHong Zhang   dx = (user->xmax - user->xmin)/((PetscReal)user->imax); dx2 = dx*dx;
154*d2567f34SHong Zhang   ierr = MatGetOwnershipRange(J,&matis,&matie);CHKERRQ(ierr);
155*d2567f34SHong Zhang   for (i=matis; i<matie; i++) {
156*d2567f34SHong Zhang     PetscScalar values[3];
157*d2567f34SHong Zhang     PetscInt    col[3];
158*d2567f34SHong Zhang     switch (user->physics_type) {
159*d2567f34SHong Zhang     case PHYSICS_DIFFUSION:
160*d2567f34SHong Zhang       values[0] = user->a*1.0/dx2;
161*d2567f34SHong Zhang       values[1] = -user->a*2.0/dx2;
162*d2567f34SHong Zhang       values[2] = user->a*1.0/dx2;
163*d2567f34SHong Zhang       break;
164*d2567f34SHong Zhang     case PHYSICS_ADVECTION:
165*d2567f34SHong Zhang       values[0] = user->a*.5/dx;
166*d2567f34SHong Zhang       values[1] = 0.;
167*d2567f34SHong Zhang       values[2] = -user->a*.5/dx;
168*d2567f34SHong Zhang       break;
169*d2567f34SHong Zhang     default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for physics type %s",PhysicsTypes[user->physics_type]);
170*d2567f34SHong Zhang     }
171*d2567f34SHong Zhang     /* periodic boundaries */
172*d2567f34SHong Zhang     if (i == 0) {
173*d2567f34SHong Zhang       col[0] = user->imax-1;
174*d2567f34SHong Zhang       col[1] = i;
175*d2567f34SHong Zhang       col[2] = i+1;
176*d2567f34SHong Zhang     } else if (i == user->imax-1) {
177*d2567f34SHong Zhang       col[0] = i-1;
178*d2567f34SHong Zhang       col[1] = i;
179*d2567f34SHong Zhang       col[2] = 0;
180*d2567f34SHong Zhang     } else {
181*d2567f34SHong Zhang       col[0] = i-1;
182*d2567f34SHong Zhang       col[1] = i;
183*d2567f34SHong Zhang       col[2] = i+1;
184*d2567f34SHong Zhang     }
185*d2567f34SHong Zhang     ierr = MatSetValues(J,1,&i,3,col,values,INSERT_VALUES);CHKERRQ(ierr);
186*d2567f34SHong Zhang   }
187*d2567f34SHong Zhang   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
188*d2567f34SHong Zhang   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
189*d2567f34SHong Zhang   PetscFunctionReturn(0);
190*d2567f34SHong Zhang }
191*d2567f34SHong Zhang 
192*d2567f34SHong Zhang /*TEST
193*d2567f34SHong Zhang 
194*d2567f34SHong Zhang   test:
195*d2567f34SHong Zhang     requires: double
196*d2567f34SHong Zhang     suffix: 1
197*d2567f34SHong Zhang     nsize: {{1 2}}
198*d2567f34SHong Zhang     args: -ts_max_steps 5 -ts_monitor -ksp_monitor_short -pc_type pbjacobi -ksp_atol 1e-6 -ts_type irk -ts_irk_nstages 2
199*d2567f34SHong Zhang 
200*d2567f34SHong Zhang   test:
201*d2567f34SHong Zhang     requires: double
202*d2567f34SHong Zhang     suffix: 2
203*d2567f34SHong Zhang     args: -ts_max_steps 5 -ts_monitor -ksp_monitor_short -pc_type pbjacobi -ksp_atol 1e-6 -ts_type irk -ts_irk_nstages 3
204*d2567f34SHong Zhang 
205*d2567f34SHong Zhang   test:
206*d2567f34SHong Zhang     requires: hpddm double
207*d2567f34SHong Zhang     suffix: 3
208*d2567f34SHong Zhang     args: -ts_max_steps 5 -ts_monitor -ksp_monitor_short -pc_type pbjacobi -ksp_atol 1e-6 -ts_type irk -ts_irk_nstages 3 -ksp_view_final_residual -ksp_hpddm_type bgcrodr -ksp_type hpddm
209*d2567f34SHong Zhang TEST*/
210