xref: /petsc/src/ts/tutorials/ex74.c (revision d0609ced746bc51b019815ca91d747429db24893)
1d2567f34SHong Zhang static char help[] = "Solves the constant-coefficient 1D heat equation \n\
2d2567f34SHong Zhang with an Implicit Runge-Kutta method using MatKAIJ.                  \n\
3d2567f34SHong Zhang                                                                     \n\
4d2567f34SHong Zhang     du      d^2 u                                                   \n\
5d2567f34SHong Zhang     --  = a ----- ; 0 <= x <= 1;                                    \n\
6d2567f34SHong Zhang     dt      dx^2                                                    \n\
7d2567f34SHong Zhang                                                                     \n\
8d2567f34SHong Zhang   with periodic boundary conditions                                 \n\
9d2567f34SHong Zhang                                                                     \n\
10d2567f34SHong Zhang 2nd order central discretization in space:                          \n\
11d2567f34SHong Zhang                                                                     \n\
12d2567f34SHong Zhang    [ d^2 u ]     u_{i+1} - 2u_i + u_{i-1}                           \n\
13d2567f34SHong Zhang    [ ----- ]  =  ------------------------                           \n\
14d2567f34SHong Zhang    [ dx^2  ]i              h^2                                      \n\
15d2567f34SHong Zhang                                                                     \n\
16d2567f34SHong Zhang     i = grid index;    h = x_{i+1}-x_i (Uniform)                    \n\
17d2567f34SHong Zhang     0 <= i < n         h = 1.0/n                                    \n\
18d2567f34SHong Zhang                                                                     \n\
19d2567f34SHong Zhang Thus,                                                               \n\
20d2567f34SHong Zhang                                                                     \n\
21d2567f34SHong Zhang    du                                                               \n\
22d2567f34SHong Zhang    --  = Ju;  J = (a/h^2) tridiagonal(1,-2,1)_n                     \n\
23d2567f34SHong Zhang    dt                                                               \n\
24d2567f34SHong Zhang                                                                     \n\
25d2567f34SHong Zhang This example is a TS version of the KSP ex74.c tutorial.            \n";
26d2567f34SHong Zhang 
27d2567f34SHong Zhang #include <petscts.h>
28d2567f34SHong Zhang 
29d2567f34SHong Zhang typedef enum {
30d2567f34SHong Zhang   PHYSICS_DIFFUSION,
31d2567f34SHong Zhang   PHYSICS_ADVECTION
32d2567f34SHong Zhang } PhysicsType;
33d2567f34SHong Zhang const char *const PhysicsTypes[] = {"DIFFUSION","ADVECTION","PhysicsType","PHYSICS_",NULL};
34d2567f34SHong Zhang 
35d2567f34SHong Zhang typedef struct Context {
36d2567f34SHong Zhang   PetscReal     a;              /* diffusion coefficient      */
37d2567f34SHong Zhang   PetscReal     xmin,xmax;      /* domain bounds              */
38d2567f34SHong Zhang   PetscInt      imax;           /* number of grid points      */
39d2567f34SHong Zhang   PhysicsType   physics_type;
40d2567f34SHong Zhang } UserContext;
41d2567f34SHong Zhang 
42d2567f34SHong Zhang static PetscErrorCode ExactSolution(Vec,void*,PetscReal);
43d2567f34SHong Zhang static PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*);
44d2567f34SHong Zhang 
45d2567f34SHong Zhang int main(int argc, char **argv)
46d2567f34SHong Zhang {
47d2567f34SHong Zhang   TS             ts;
48d2567f34SHong Zhang   Mat            A;
49d2567f34SHong Zhang   Vec            u,uex;
50d2567f34SHong Zhang   UserContext    ctxt;
51d2567f34SHong Zhang   PetscReal      err,ftime;
52d2567f34SHong Zhang 
539566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
54d2567f34SHong Zhang   /* default value */
55d2567f34SHong Zhang   ctxt.a       = 0.1;
56d2567f34SHong Zhang   ctxt.xmin    = 0.0;
57d2567f34SHong Zhang   ctxt.xmax    = 1.0;
58d2567f34SHong Zhang   ctxt.imax    = 40;
59d2567f34SHong Zhang   ctxt.physics_type = PHYSICS_DIFFUSION;
60d2567f34SHong Zhang 
61*d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"IRK options","");
629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-a","diffusion coefficient","<1.0>",ctxt.a,&ctxt.a,NULL));
639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt ("-imax","grid size","<20>",ctxt.imax,&ctxt.imax,NULL));
649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-xmin","xmin","<0.0>",ctxt.xmin,&ctxt.xmin,NULL));
659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-xmax","xmax","<1.0>",ctxt.xmax,&ctxt.xmax,NULL));
669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-physics_type","Type of process to discretize","",PhysicsTypes,(PetscEnum)ctxt.physics_type,(PetscEnum*)&ctxt.physics_type,NULL));
67*d0609cedSBarry Smith   PetscOptionsEnd();
68d2567f34SHong Zhang 
69d2567f34SHong Zhang   /* allocate and initialize solution vector and exact solution */
709566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&u));
719566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(u,PETSC_DECIDE,ctxt.imax));
729566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(u));
739566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u,&uex));
74d2567f34SHong Zhang   /* initial solution */
759566063dSJacob Faibussowitsch   PetscCall(ExactSolution(u,&ctxt,0.0));
76d2567f34SHong Zhang 
779566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
789566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATAIJ));
799566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,ctxt.imax,ctxt.imax));
809566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
81d2567f34SHong Zhang 
82d2567f34SHong Zhang   /* Create and set options for TS */
839566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
849566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts,TS_LINEAR));
859566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,0.125));
869566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts,u));
879566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts,10));
889566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,1.0));
899566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
909566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&ctxt));
919566063dSJacob Faibussowitsch   PetscCall(RHSJacobian(ts,0,u,A,A,&ctxt));
929566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,&ctxt));
939566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
949566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,u));
95d2567f34SHong Zhang 
969566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts,&ftime));
97d2567f34SHong Zhang   /* exact   solution */
989566063dSJacob Faibussowitsch   PetscCall(ExactSolution(uex,&ctxt,ftime));
99d2567f34SHong Zhang 
100d2567f34SHong Zhang   /* Calculate error in final solution */
1019566063dSJacob Faibussowitsch   PetscCall(VecAYPX(uex,-1.0,u));
1029566063dSJacob Faibussowitsch   PetscCall(VecNorm(uex,NORM_2,&err));
103d2567f34SHong Zhang   err  = PetscSqrtReal(err*err/((PetscReal)ctxt.imax));
1049566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"L2 norm of the numerical error = %g (time=%g)\n",(double)err,(double)ftime));
105d2567f34SHong Zhang 
106d2567f34SHong Zhang   /* Free up memory */
1079566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1089566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&uex));
1109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
1119566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
112b122ec5aSJacob Faibussowitsch   return 0;
113d2567f34SHong Zhang }
114d2567f34SHong Zhang 
115d2567f34SHong Zhang PetscErrorCode ExactSolution(Vec u,void *c,PetscReal t)
116d2567f34SHong Zhang {
117d2567f34SHong Zhang   UserContext     *ctxt = (UserContext*) c;
118d2567f34SHong Zhang   PetscInt        i,is,ie;
119d2567f34SHong Zhang   PetscScalar     *uarr;
120d2567f34SHong Zhang   PetscReal       x,dx,a=ctxt->a,pi=PETSC_PI;
121d2567f34SHong Zhang 
122d2567f34SHong Zhang   PetscFunctionBegin;
123d2567f34SHong Zhang   dx = (ctxt->xmax - ctxt->xmin)/((PetscReal) ctxt->imax);
1249566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(u,&is,&ie));
1259566063dSJacob Faibussowitsch   PetscCall(VecGetArray(u,&uarr));
126d2567f34SHong Zhang   for (i=is; i<ie; i++) {
127d2567f34SHong Zhang     x          = i * dx;
128d2567f34SHong Zhang     switch (ctxt->physics_type) {
129d2567f34SHong Zhang     case PHYSICS_DIFFUSION:
130d2567f34SHong Zhang       uarr[i-is] = PetscExpScalar(-4.0*pi*pi*a*t)*PetscSinScalar(2*pi*x);
131d2567f34SHong Zhang       break;
132d2567f34SHong Zhang     case PHYSICS_ADVECTION:
133d2567f34SHong Zhang       uarr[i-is] = PetscSinScalar(2*pi*(x - a*t));
134d2567f34SHong Zhang       break;
13598921bdaSJacob Faibussowitsch     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for physics type %s",PhysicsTypes[ctxt->physics_type]);
136d2567f34SHong Zhang     }
137d2567f34SHong Zhang   }
1389566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(u,&uarr));
139d2567f34SHong Zhang   PetscFunctionReturn(0);
140d2567f34SHong Zhang }
141d2567f34SHong Zhang 
142d2567f34SHong Zhang static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat J,Mat Jpre,void *ctx)
143d2567f34SHong Zhang {
144d2567f34SHong Zhang   UserContext    *user = (UserContext*) ctx;
145d2567f34SHong Zhang   PetscInt       matis,matie,i;
146d2567f34SHong Zhang   PetscReal      dx,dx2;
147d2567f34SHong Zhang 
148d2567f34SHong Zhang   PetscFunctionBegin;
149d2567f34SHong Zhang   dx = (user->xmax - user->xmin)/((PetscReal)user->imax); dx2 = dx*dx;
1509566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(J,&matis,&matie));
151d2567f34SHong Zhang   for (i=matis; i<matie; i++) {
152d2567f34SHong Zhang     PetscScalar values[3];
153d2567f34SHong Zhang     PetscInt    col[3];
154d2567f34SHong Zhang     switch (user->physics_type) {
155d2567f34SHong Zhang     case PHYSICS_DIFFUSION:
156d2567f34SHong Zhang       values[0] = user->a*1.0/dx2;
157d2567f34SHong Zhang       values[1] = -user->a*2.0/dx2;
158d2567f34SHong Zhang       values[2] = user->a*1.0/dx2;
159d2567f34SHong Zhang       break;
160d2567f34SHong Zhang     case PHYSICS_ADVECTION:
161d2567f34SHong Zhang       values[0] = user->a*.5/dx;
162d2567f34SHong Zhang       values[1] = 0.;
163d2567f34SHong Zhang       values[2] = -user->a*.5/dx;
164d2567f34SHong Zhang       break;
16598921bdaSJacob Faibussowitsch     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for physics type %s",PhysicsTypes[user->physics_type]);
166d2567f34SHong Zhang     }
167d2567f34SHong Zhang     /* periodic boundaries */
168d2567f34SHong Zhang     if (i == 0) {
169d2567f34SHong Zhang       col[0] = user->imax-1;
170d2567f34SHong Zhang       col[1] = i;
171d2567f34SHong Zhang       col[2] = i+1;
172d2567f34SHong Zhang     } else if (i == user->imax-1) {
173d2567f34SHong Zhang       col[0] = i-1;
174d2567f34SHong Zhang       col[1] = i;
175d2567f34SHong Zhang       col[2] = 0;
176d2567f34SHong Zhang     } else {
177d2567f34SHong Zhang       col[0] = i-1;
178d2567f34SHong Zhang       col[1] = i;
179d2567f34SHong Zhang       col[2] = i+1;
180d2567f34SHong Zhang     }
1819566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J,1,&i,3,col,values,INSERT_VALUES));
182d2567f34SHong Zhang   }
1839566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
1849566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
185d2567f34SHong Zhang   PetscFunctionReturn(0);
186d2567f34SHong Zhang }
187d2567f34SHong Zhang 
188d2567f34SHong Zhang /*TEST
189d2567f34SHong Zhang 
190d2567f34SHong Zhang   test:
191d2567f34SHong Zhang     requires: double
192d2567f34SHong Zhang     suffix: 1
193d2567f34SHong Zhang     nsize: {{1 2}}
194d2567f34SHong 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
195d2567f34SHong Zhang 
196d2567f34SHong Zhang   test:
197d2567f34SHong Zhang     requires: double
198d2567f34SHong Zhang     suffix: 2
199d2567f34SHong 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
200d2567f34SHong Zhang 
2012dd49c90SPierre Jolivet   testset:
2022dd49c90SPierre Jolivet     requires: hpddm
2032dd49c90SPierre Jolivet     args: -ts_max_steps 5 -ts_monitor -ksp_monitor_short -pc_type pbjacobi -ksp_atol 1e-4 -ts_type irk -ts_irk_nstages 3 -ksp_view_final_residual -ksp_hpddm_type gcrodr -ksp_type hpddm -ksp_hpddm_precision {{single double}shared output}
204d2567f34SHong Zhang     test:
205d2567f34SHong Zhang       suffix: 3
2062dd49c90SPierre Jolivet       requires: double
2072dd49c90SPierre Jolivet     test:
2082dd49c90SPierre Jolivet       suffix: 3_single
2092dd49c90SPierre Jolivet       requires: single
210d2567f34SHong Zhang TEST*/
211