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; 249566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 259566063dSJacob Faibussowitsch PetscCall(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]); 289566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 299566063dSJacob Faibussowitsch PetscCall(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; 399566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 409566063dSJacob Faibussowitsch PetscCall(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]); 429566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 439566063dSJacob Faibussowitsch PetscCall(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; 539566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 549566063dSJacob Faibussowitsch PetscCall(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]); 569566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 579566063dSJacob Faibussowitsch PetscCall(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; 669566063dSJacob Faibussowitsch PetscCall(VecGetArray(U,&u)); 67c4762a1bSJed Brown u[0] = PetscSqrtScalar(1.0+PetscCosScalar(t)); 68c4762a1bSJed Brown u[1] = PetscSqrtScalar(2.0+PetscCosScalar(5.0*t)); 699566063dSJacob Faibussowitsch PetscCall(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 PetscMPIInt size; 79c4762a1bSJed Brown AppCtx ctx; 80c4762a1bSJed Brown PetscScalar *u; 81c4762a1bSJed Brown IS iss; 82c4762a1bSJed Brown IS isf; 83c4762a1bSJed Brown PetscInt *indicess; 84c4762a1bSJed Brown PetscInt *indicesf; 85c4762a1bSJed Brown PetscInt n=2; 86c4762a1bSJed Brown PetscReal error,tt; 87c4762a1bSJed Brown 88c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 89c4762a1bSJed Brown Initialize program 90c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 919566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 929566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 933c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs"); 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 96c4762a1bSJed Brown Create index for slow part and fast part 97c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1,&indicess)); 99c4762a1bSJed Brown indicess[0]=0; 1009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1,&indicesf)); 101c4762a1bSJed Brown indicesf[0]=1; 1029566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,1,indicess,PETSC_COPY_VALUES,&iss)); 1039566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,1,indicesf,PETSC_COPY_VALUES,&isf)); 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 106c4762a1bSJed Brown Create necesary vector 107c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1089566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&U)); 1099566063dSJacob Faibussowitsch PetscCall(VecSetSizes(U,n,PETSC_DETERMINE)); 1109566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(U)); 1119566063dSJacob Faibussowitsch PetscCall(VecDuplicate(U,&Utrue)); 1129566063dSJacob Faibussowitsch PetscCall(VecCopy(U,Utrue)); 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 115c4762a1bSJed Brown Set initial condition 116c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1179566063dSJacob Faibussowitsch PetscCall(VecGetArray(U,&u)); 118c4762a1bSJed Brown u[0] = PetscSqrtScalar(2.0); 119c4762a1bSJed Brown u[1] = PetscSqrtScalar(3.0); 1209566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U,&u)); 121c4762a1bSJed Brown 122c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 123c4762a1bSJed Brown Create timestepping solver context 124c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1259566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 1269566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSMPRK)); 127c4762a1bSJed Brown 1289566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,NULL,(TSRHSFunction)RHSFunction,&ctx)); 1299566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts,"slow",iss)); 1309566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts,"fast",isf)); 1319566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts,"slow",NULL,(TSRHSFunction)RHSFunctionslow,&ctx)); 1329566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts,"fast",NULL,(TSRHSFunction)RHSFunctionfast,&ctx)); 133c4762a1bSJed Brown 134c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 135c4762a1bSJed Brown Set initial conditions 136c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1379566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts,U)); 138c4762a1bSJed Brown 139c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 140c4762a1bSJed Brown Set solver options 141c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 142*d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ODE options",""); 143c4762a1bSJed Brown { 144c4762a1bSJed Brown ctx.Tf = 0.3; 145c4762a1bSJed Brown ctx.dt = 0.01; 1469566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Tf","","",ctx.Tf,&ctx.Tf,NULL)); 1479566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-dt","","",ctx.dt,&ctx.dt,NULL)); 148c4762a1bSJed Brown } 149*d0609cedSBarry Smith PetscOptionsEnd(); 1509566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,ctx.Tf)); 1519566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,ctx.dt)); 1529566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 1539566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 154c4762a1bSJed Brown 155c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 156c4762a1bSJed Brown Solve linear system 157c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1589566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,U)); 1599566063dSJacob Faibussowitsch PetscCall(VecView(U,PETSC_VIEWER_STDOUT_WORLD)); 160c4762a1bSJed Brown 161c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 162c4762a1bSJed Brown Check the error of the Petsc solution 163c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1649566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts,&tt)); 1659566063dSJacob Faibussowitsch PetscCall(sol_true(tt,Utrue)); 1669566063dSJacob Faibussowitsch PetscCall(VecAXPY(Utrue,-1.0,U)); 1679566063dSJacob Faibussowitsch PetscCall(VecNorm(Utrue,NORM_2,&error)); 168c4762a1bSJed Brown 169c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 170c4762a1bSJed Brown Print norm2 error 171c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1729566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"l2 error norm: %g\n", error)); 173c4762a1bSJed Brown 174c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 175c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 176c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 1789566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Utrue)); 1809566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iss)); 1819566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isf)); 1829566063dSJacob Faibussowitsch PetscCall(PetscFree(indicess)); 1839566063dSJacob Faibussowitsch PetscCall(PetscFree(indicesf)); 1849566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 185b122ec5aSJacob Faibussowitsch return 0; 186c4762a1bSJed Brown } 187109dc152SHong Zhang 188109dc152SHong Zhang /*TEST 189109dc152SHong Zhang build: 190f56ea12dSJed Brown requires: !complex 191109dc152SHong Zhang 192109dc152SHong Zhang test: 193109dc152SHong Zhang 194109dc152SHong Zhang TEST*/ 195