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