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 PetscErrorCode ierr; 53d2567f34SHong Zhang 54*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 55d2567f34SHong Zhang 56d2567f34SHong Zhang /* default value */ 57d2567f34SHong Zhang ctxt.a = 0.1; 58d2567f34SHong Zhang ctxt.xmin = 0.0; 59d2567f34SHong Zhang ctxt.xmax = 1.0; 60d2567f34SHong Zhang ctxt.imax = 40; 61d2567f34SHong Zhang ctxt.physics_type = PHYSICS_DIFFUSION; 62d2567f34SHong Zhang 63d2567f34SHong Zhang ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"IRK options","");CHKERRQ(ierr); 645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-a","diffusion coefficient","<1.0>",ctxt.a,&ctxt.a,NULL)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt ("-imax","grid size","<20>",ctxt.imax,&ctxt.imax,NULL)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-xmin","xmin","<0.0>",ctxt.xmin,&ctxt.xmin,NULL)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-xmax","xmax","<1.0>",ctxt.xmax,&ctxt.xmax,NULL)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEnum("-physics_type","Type of process to discretize","",PhysicsTypes,(PetscEnum)ctxt.physics_type,(PetscEnum*)&ctxt.physics_type,NULL)); 69d2567f34SHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 70d2567f34SHong Zhang 71d2567f34SHong Zhang /* allocate and initialize solution vector and exact solution */ 725f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&u)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(u,PETSC_DECIDE,ctxt.imax)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(u)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&uex)); 76d2567f34SHong Zhang /* initial solution */ 775f80ce2aSJacob Faibussowitsch CHKERRQ(ExactSolution(u,&ctxt,0.0)); 78d2567f34SHong Zhang 795f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATAIJ)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,ctxt.imax,ctxt.imax)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 83d2567f34SHong Zhang 84d2567f34SHong Zhang /* Create and set options for TS */ 855f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetProblemType(ts,TS_LINEAR)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,0.125)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSolution(ts,u)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxSteps(ts,10)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,1.0)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 925f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&ctxt)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(RHSJacobian(ts,0,u,A,A,&ctxt)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,&ctxt)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,u)); 97d2567f34SHong Zhang 985f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSolveTime(ts,&ftime)); 99d2567f34SHong Zhang /* exact solution */ 1005f80ce2aSJacob Faibussowitsch CHKERRQ(ExactSolution(uex,&ctxt,ftime)); 101d2567f34SHong Zhang 102d2567f34SHong Zhang /* Calculate error in final solution */ 1035f80ce2aSJacob Faibussowitsch CHKERRQ(VecAYPX(uex,-1.0,u)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(uex,NORM_2,&err)); 105d2567f34SHong Zhang err = PetscSqrtReal(err*err/((PetscReal)ctxt.imax)); 1065f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"L2 norm of the numerical error = %g (time=%g)\n",(double)err,(double)ftime)); 107d2567f34SHong Zhang 108d2567f34SHong Zhang /* Free up memory */ 1095f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&uex)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 113*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 114*b122ec5aSJacob Faibussowitsch return 0; 115d2567f34SHong Zhang } 116d2567f34SHong Zhang 117d2567f34SHong Zhang PetscErrorCode ExactSolution(Vec u,void *c,PetscReal t) 118d2567f34SHong Zhang { 119d2567f34SHong Zhang UserContext *ctxt = (UserContext*) c; 120d2567f34SHong Zhang PetscInt i,is,ie; 121d2567f34SHong Zhang PetscScalar *uarr; 122d2567f34SHong Zhang PetscReal x,dx,a=ctxt->a,pi=PETSC_PI; 123d2567f34SHong Zhang 124d2567f34SHong Zhang PetscFunctionBegin; 125d2567f34SHong Zhang dx = (ctxt->xmax - ctxt->xmin)/((PetscReal) ctxt->imax); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(u,&is,&ie)); 1275f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(u,&uarr)); 128d2567f34SHong Zhang for (i=is; i<ie; i++) { 129d2567f34SHong Zhang x = i * dx; 130d2567f34SHong Zhang switch (ctxt->physics_type) { 131d2567f34SHong Zhang case PHYSICS_DIFFUSION: 132d2567f34SHong Zhang uarr[i-is] = PetscExpScalar(-4.0*pi*pi*a*t)*PetscSinScalar(2*pi*x); 133d2567f34SHong Zhang break; 134d2567f34SHong Zhang case PHYSICS_ADVECTION: 135d2567f34SHong Zhang uarr[i-is] = PetscSinScalar(2*pi*(x - a*t)); 136d2567f34SHong Zhang break; 13798921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for physics type %s",PhysicsTypes[ctxt->physics_type]); 138d2567f34SHong Zhang } 139d2567f34SHong Zhang } 1405f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(u,&uarr)); 141d2567f34SHong Zhang PetscFunctionReturn(0); 142d2567f34SHong Zhang } 143d2567f34SHong Zhang 144d2567f34SHong Zhang static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat J,Mat Jpre,void *ctx) 145d2567f34SHong Zhang { 146d2567f34SHong Zhang UserContext *user = (UserContext*) ctx; 147d2567f34SHong Zhang PetscInt matis,matie,i; 148d2567f34SHong Zhang PetscReal dx,dx2; 149d2567f34SHong Zhang 150d2567f34SHong Zhang PetscFunctionBegin; 151d2567f34SHong Zhang dx = (user->xmax - user->xmin)/((PetscReal)user->imax); dx2 = dx*dx; 1525f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(J,&matis,&matie)); 153d2567f34SHong Zhang for (i=matis; i<matie; i++) { 154d2567f34SHong Zhang PetscScalar values[3]; 155d2567f34SHong Zhang PetscInt col[3]; 156d2567f34SHong Zhang switch (user->physics_type) { 157d2567f34SHong Zhang case PHYSICS_DIFFUSION: 158d2567f34SHong Zhang values[0] = user->a*1.0/dx2; 159d2567f34SHong Zhang values[1] = -user->a*2.0/dx2; 160d2567f34SHong Zhang values[2] = user->a*1.0/dx2; 161d2567f34SHong Zhang break; 162d2567f34SHong Zhang case PHYSICS_ADVECTION: 163d2567f34SHong Zhang values[0] = user->a*.5/dx; 164d2567f34SHong Zhang values[1] = 0.; 165d2567f34SHong Zhang values[2] = -user->a*.5/dx; 166d2567f34SHong Zhang break; 16798921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for physics type %s",PhysicsTypes[user->physics_type]); 168d2567f34SHong Zhang } 169d2567f34SHong Zhang /* periodic boundaries */ 170d2567f34SHong Zhang if (i == 0) { 171d2567f34SHong Zhang col[0] = user->imax-1; 172d2567f34SHong Zhang col[1] = i; 173d2567f34SHong Zhang col[2] = i+1; 174d2567f34SHong Zhang } else if (i == user->imax-1) { 175d2567f34SHong Zhang col[0] = i-1; 176d2567f34SHong Zhang col[1] = i; 177d2567f34SHong Zhang col[2] = 0; 178d2567f34SHong Zhang } else { 179d2567f34SHong Zhang col[0] = i-1; 180d2567f34SHong Zhang col[1] = i; 181d2567f34SHong Zhang col[2] = i+1; 182d2567f34SHong Zhang } 1835f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(J,1,&i,3,col,values,INSERT_VALUES)); 184d2567f34SHong Zhang } 1855f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 1865f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 187d2567f34SHong Zhang PetscFunctionReturn(0); 188d2567f34SHong Zhang } 189d2567f34SHong Zhang 190d2567f34SHong Zhang /*TEST 191d2567f34SHong Zhang 192d2567f34SHong Zhang test: 193d2567f34SHong Zhang requires: double 194d2567f34SHong Zhang suffix: 1 195d2567f34SHong Zhang nsize: {{1 2}} 196d2567f34SHong 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 197d2567f34SHong Zhang 198d2567f34SHong Zhang test: 199d2567f34SHong Zhang requires: double 200d2567f34SHong Zhang suffix: 2 201d2567f34SHong 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 202d2567f34SHong Zhang 2032dd49c90SPierre Jolivet testset: 2042dd49c90SPierre Jolivet requires: hpddm 2052dd49c90SPierre 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} 206d2567f34SHong Zhang test: 207d2567f34SHong Zhang suffix: 3 2082dd49c90SPierre Jolivet requires: double 2092dd49c90SPierre Jolivet test: 2102dd49c90SPierre Jolivet suffix: 3_single 2112dd49c90SPierre Jolivet requires: single 212d2567f34SHong Zhang TEST*/ 213