1c4762a1bSJed Brown static char help[] = "Basic problem for multi-rate method.\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /*F 4c4762a1bSJed Brown 5c4762a1bSJed Brown \begin{eqnarray} 6c4762a1bSJed Brown ys' = -2.0*\frac{-1.0+ys^2.0-\cos(t)}{2.0*ys}+0.05*\frac{-2.0+yf^2-\cos(5.0*t)}{2.0*yf}-\frac{\sin(t)}{2.0*ys}\\ 7c4762a1bSJed Brown yf' = 0.05*\frac{-1.0+ys^2-\cos(t)}{2.0*ys}-\frac{-2.0+yf^2-\cos(5.0*t)}{2.0*yf}-5.0*\frac{\sin(5.0*t)}{2.0*yf}\\ 8c4762a1bSJed Brown \end{eqnarray} 9c4762a1bSJed Brown 10c4762a1bSJed Brown F*/ 11c4762a1bSJed Brown 12c4762a1bSJed Brown #include <petscts.h> 13c4762a1bSJed Brown 14c4762a1bSJed Brown typedef struct { 15c4762a1bSJed Brown PetscReal Tf,dt; 16c4762a1bSJed Brown } AppCtx; 17c4762a1bSJed Brown 18c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx) 19c4762a1bSJed Brown { 20c4762a1bSJed Brown const PetscScalar *u; 21c4762a1bSJed Brown PetscScalar *f; 22c4762a1bSJed Brown 23c4762a1bSJed Brown PetscFunctionBegin; 245f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(U,&u)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(F,&f)); 26c4762a1bSJed Brown f[0] = -2.0*(-1.0+u[0]*u[0]-PetscCosScalar(t))/(2.0*u[0])+0.05*(-2.0+u[1]*u[1]-PetscCosScalar(5.0*t))/(2.0*u[1])-PetscSinScalar(t)/(2.0*u[0]); 27c4762a1bSJed Brown f[1] = 0.05*(-1.0+u[0]*u[0]-PetscCosScalar(t))/(2.0*u[0])-(-2.0+u[1]*u[1]-PetscCosScalar(5.0*t))/(2.0*u[1])-5.0*PetscSinScalar(5.0*t)/(2.0*u[1]); 285f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(U,&u)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(F,&f)); 30c4762a1bSJed Brown PetscFunctionReturn(0); 31c4762a1bSJed Brown } 32c4762a1bSJed Brown 33c4762a1bSJed Brown static PetscErrorCode RHSFunctionslow(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx) 34c4762a1bSJed Brown { 35c4762a1bSJed Brown const PetscScalar *u; 36c4762a1bSJed Brown PetscScalar *f; 37c4762a1bSJed Brown 38c4762a1bSJed Brown PetscFunctionBegin; 395f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(U,&u)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(F,&f)); 41c4762a1bSJed Brown f[0] = -2.0*(-1.0+u[0]*u[0]-PetscCosScalar(t))/(2.0*u[0])+0.05*(-2.0+u[1]*u[1]-PetscCosScalar(5.0*t))/(2.0*u[1])-PetscSinScalar(t)/(2.0*u[0]); 425f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(U,&u)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(F,&f)); 44c4762a1bSJed Brown PetscFunctionReturn(0); 45c4762a1bSJed Brown } 46c4762a1bSJed Brown 47c4762a1bSJed Brown static PetscErrorCode RHSFunctionfast(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx) 48c4762a1bSJed Brown { 49c4762a1bSJed Brown const PetscScalar *u; 50c4762a1bSJed Brown PetscScalar *f; 51c4762a1bSJed Brown 52c4762a1bSJed Brown PetscFunctionBegin; 535f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(U,&u)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(F,&f)); 55c4762a1bSJed Brown f[0] = 0.05*(-1.0+u[0]*u[0]-PetscCosScalar(t))/(2.0*u[0])-(-2.0+u[1]*u[1]-PetscCosScalar(5.0*t))/(2.0*u[1])-5.0*PetscSinScalar(5.0*t)/(2.0*u[1]); 565f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(U,&u)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(F,&f)); 58c4762a1bSJed Brown PetscFunctionReturn(0); 59c4762a1bSJed Brown } 60c4762a1bSJed Brown 61c4762a1bSJed Brown static PetscErrorCode sol_true(PetscReal t,Vec U) 62c4762a1bSJed Brown { 63c4762a1bSJed Brown PetscScalar *u; 64c4762a1bSJed Brown 65c4762a1bSJed Brown PetscFunctionBegin; 665f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(U,&u)); 67c4762a1bSJed Brown u[0] = PetscSqrtScalar(1.0+PetscCosScalar(t)); 68c4762a1bSJed Brown u[1] = PetscSqrtScalar(2.0+PetscCosScalar(5.0*t)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(U,&u)); 70c4762a1bSJed Brown PetscFunctionReturn(0); 71c4762a1bSJed Brown } 72c4762a1bSJed Brown 73c4762a1bSJed Brown int main(int argc,char **argv) 74c4762a1bSJed Brown { 75c4762a1bSJed Brown TS ts; /* ODE integrator */ 76c4762a1bSJed Brown Vec U; /* solution will be stored here */ 77c4762a1bSJed Brown Vec Utrue; 78c4762a1bSJed Brown PetscErrorCode ierr; 79c4762a1bSJed Brown PetscMPIInt size; 80c4762a1bSJed Brown AppCtx ctx; 81c4762a1bSJed Brown PetscScalar *u; 82c4762a1bSJed Brown IS iss; 83c4762a1bSJed Brown IS isf; 84c4762a1bSJed Brown PetscInt *indicess; 85c4762a1bSJed Brown PetscInt *indicesf; 86c4762a1bSJed Brown PetscInt n=2; 87c4762a1bSJed Brown PetscReal error,tt; 88c4762a1bSJed Brown 89c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 90c4762a1bSJed Brown Initialize program 91c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 92*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 935f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 943c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs"); 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 97c4762a1bSJed Brown Create index for slow part and fast part 98c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(1,&indicess)); 100c4762a1bSJed Brown indicess[0]=0; 1015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(1,&indicesf)); 102c4762a1bSJed Brown indicesf[0]=1; 1035f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,1,indicess,PETSC_COPY_VALUES,&iss)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,1,indicesf,PETSC_COPY_VALUES,&isf)); 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 107c4762a1bSJed Brown Create necesary vector 108c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1095f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&U)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(U,n,PETSC_DETERMINE)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(U)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(U,&Utrue)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(U,Utrue)); 114c4762a1bSJed Brown 115c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 116c4762a1bSJed Brown Set initial condition 117c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1185f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(U,&u)); 119c4762a1bSJed Brown u[0] = PetscSqrtScalar(2.0); 120c4762a1bSJed Brown u[1] = PetscSqrtScalar(3.0); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(U,&u)); 122c4762a1bSJed Brown 123c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 124c4762a1bSJed Brown Create timestepping solver context 125c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1265f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 1275f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts,TSMPRK)); 128c4762a1bSJed Brown 1295f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts,NULL,(TSRHSFunction)RHSFunction,&ctx)); 1305f80ce2aSJacob Faibussowitsch CHKERRQ(TSRHSSplitSetIS(ts,"slow",iss)); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(TSRHSSplitSetIS(ts,"fast",isf)); 1325f80ce2aSJacob Faibussowitsch CHKERRQ(TSRHSSplitSetRHSFunction(ts,"slow",NULL,(TSRHSFunction)RHSFunctionslow,&ctx)); 1335f80ce2aSJacob Faibussowitsch CHKERRQ(TSRHSSplitSetRHSFunction(ts,"fast",NULL,(TSRHSFunction)RHSFunctionfast,&ctx)); 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 136c4762a1bSJed Brown Set initial conditions 137c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1385f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSolution(ts,U)); 139c4762a1bSJed Brown 140c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 141c4762a1bSJed Brown Set solver options 142c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 143c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ODE options","");CHKERRQ(ierr); 144c4762a1bSJed Brown { 145c4762a1bSJed Brown ctx.Tf = 0.3; 146c4762a1bSJed Brown ctx.dt = 0.01; 1475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsScalar("-Tf","","",ctx.Tf,&ctx.Tf,NULL)); 1485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsScalar("-dt","","",ctx.dt,&ctx.dt,NULL)); 149c4762a1bSJed Brown } 150c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 1515f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,ctx.Tf)); 1525f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,ctx.dt)); 1535f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 1545f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 155c4762a1bSJed Brown 156c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 157c4762a1bSJed Brown Solve linear system 158c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1595f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,U)); 1605f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(U,PETSC_VIEWER_STDOUT_WORLD)); 161c4762a1bSJed Brown 162c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 163c4762a1bSJed Brown Check the error of the Petsc solution 164c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1655f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetTime(ts,&tt)); 1665f80ce2aSJacob Faibussowitsch CHKERRQ(sol_true(tt,Utrue)); 1675f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(Utrue,-1.0,U)); 1685f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(Utrue,NORM_2,&error)); 169c4762a1bSJed Brown 170c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 171c4762a1bSJed Brown Print norm2 error 172c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"l2 error norm: %g\n", error)); 174c4762a1bSJed Brown 175c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 176c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 177c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1785f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&U)); 1795f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 1805f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Utrue)); 1815f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iss)); 1825f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isf)); 1835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(indicess)); 1845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(indicesf)); 185*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 186*b122ec5aSJacob Faibussowitsch return 0; 187c4762a1bSJed Brown } 188109dc152SHong Zhang 189109dc152SHong Zhang /*TEST 190109dc152SHong Zhang build: 191f56ea12dSJed Brown requires: !complex 192109dc152SHong Zhang 193109dc152SHong Zhang test: 194109dc152SHong Zhang 195109dc152SHong Zhang TEST*/ 196