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' = -sin(a*t)\\ 7c4762a1bSJed Brown yf' = bcos(b*t)ys-sin(b*t)sin(a*t)\\ 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 a,b,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 PetscErrorCode ierr; 21c4762a1bSJed Brown const PetscScalar *u; 22c4762a1bSJed Brown PetscScalar *f; 23c4762a1bSJed Brown 24c4762a1bSJed Brown PetscFunctionBegin; 25c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 26c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 27c4762a1bSJed Brown f[0] = -PetscSinScalar(ctx->a*t); 28c4762a1bSJed Brown f[1] = ctx->b*PetscCosScalar(ctx->b*t)*u[0]-PetscSinScalar(ctx->b*t)*PetscSinScalar(ctx->a*t); 29c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 30c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 31c4762a1bSJed Brown PetscFunctionReturn(0); 32c4762a1bSJed Brown } 33c4762a1bSJed Brown 34c4762a1bSJed Brown static PetscErrorCode RHSFunctionslow(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx) 35c4762a1bSJed Brown { 36c4762a1bSJed Brown PetscErrorCode ierr; 37c4762a1bSJed Brown const PetscScalar *u; 38c4762a1bSJed Brown PetscScalar *f; 39c4762a1bSJed Brown 40c4762a1bSJed Brown PetscFunctionBegin; 41c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 42c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 43c4762a1bSJed Brown f[0] = -PetscSinScalar(ctx->a*t); 44c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 45c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 46c4762a1bSJed Brown PetscFunctionReturn(0); 47c4762a1bSJed Brown } 48c4762a1bSJed Brown 49c4762a1bSJed Brown static PetscErrorCode RHSFunctionfast(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx) 50c4762a1bSJed Brown { 51c4762a1bSJed Brown PetscErrorCode ierr; 52c4762a1bSJed Brown const PetscScalar *u; 53c4762a1bSJed Brown PetscScalar *f; 54c4762a1bSJed Brown 55c4762a1bSJed Brown PetscFunctionBegin; 56c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 57c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 58c4762a1bSJed Brown f[0] = ctx->b*PetscCosScalar(ctx->b*t)*u[0]-PetscSinScalar(ctx->b*t)*PetscSinScalar(ctx->a*t); 59c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 60c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 61c4762a1bSJed Brown PetscFunctionReturn(0); 62c4762a1bSJed Brown } 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* 65c4762a1bSJed Brown Define the analytic solution for check method easily 66c4762a1bSJed Brown */ 67c4762a1bSJed Brown static PetscErrorCode sol_true(PetscReal t,Vec U,AppCtx *ctx) 68c4762a1bSJed Brown { 69c4762a1bSJed Brown PetscErrorCode ierr; 70c4762a1bSJed Brown PetscScalar *u; 71c4762a1bSJed Brown 72c4762a1bSJed Brown PetscFunctionBegin; 73c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 74c4762a1bSJed Brown u[0] = PetscCosScalar(ctx->a*t)/ctx->a; 75c4762a1bSJed Brown u[1] = PetscSinScalar(ctx->b*t)*PetscCosScalar(ctx->a*t)/ctx->a; 76c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 77c4762a1bSJed Brown PetscFunctionReturn(0); 78c4762a1bSJed Brown } 79c4762a1bSJed Brown 80c4762a1bSJed Brown int main(int argc,char **argv) 81c4762a1bSJed Brown { 82c4762a1bSJed Brown TS ts; /* ODE integrator */ 83c4762a1bSJed Brown Vec U; /* solution will be stored here */ 84c4762a1bSJed Brown Vec Utrue; 85c4762a1bSJed Brown PetscErrorCode ierr; 86c4762a1bSJed Brown PetscMPIInt size; 87c4762a1bSJed Brown AppCtx ctx; 88c4762a1bSJed Brown PetscScalar *u; 89c4762a1bSJed Brown IS iss; 90c4762a1bSJed Brown IS isf; 91c4762a1bSJed Brown PetscInt *indicess; 92c4762a1bSJed Brown PetscInt *indicesf; 93c4762a1bSJed Brown PetscInt n=2; 94c4762a1bSJed Brown PetscScalar error,tt; 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 97c4762a1bSJed Brown Initialize program 98c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 99c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 100c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 101c4762a1bSJed Brown if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs"); 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 104c4762a1bSJed Brown Create index for slow part and fast part 105c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 106c4762a1bSJed Brown ierr = PetscMalloc1(1,&indicess);CHKERRQ(ierr); 107c4762a1bSJed Brown indicess[0]=0; 108c4762a1bSJed Brown ierr = PetscMalloc1(1,&indicesf);CHKERRQ(ierr); 109c4762a1bSJed Brown indicesf[0]=1; 110c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,1,indicess,PETSC_COPY_VALUES,&iss);CHKERRQ(ierr); 111c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,1,indicesf,PETSC_COPY_VALUES,&isf);CHKERRQ(ierr); 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 114c4762a1bSJed Brown Create necesary vector 115c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 116c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&U);CHKERRQ(ierr); 117c4762a1bSJed Brown ierr = VecSetSizes(U,n,PETSC_DETERMINE);CHKERRQ(ierr); 118c4762a1bSJed Brown ierr = VecSetFromOptions(U);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = VecDuplicate(U,&Utrue);CHKERRQ(ierr); 120c4762a1bSJed Brown ierr = VecCopy(U,Utrue);CHKERRQ(ierr); 121c4762a1bSJed Brown 122c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 123c4762a1bSJed Brown Set runtime options 124c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 125c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ODE options","");CHKERRQ(ierr); 126c4762a1bSJed Brown { 127c4762a1bSJed Brown ctx.a = 1.0; 128c4762a1bSJed Brown ctx.b = 25.0; 129c4762a1bSJed Brown ierr = PetscOptionsScalar("-a","","",ctx.a,&ctx.a,NULL);CHKERRQ(ierr); 130c4762a1bSJed Brown ierr = PetscOptionsScalar("-b","","",ctx.b,&ctx.b,NULL);CHKERRQ(ierr); 131c4762a1bSJed Brown ctx.Tf = 5.0; 132c4762a1bSJed Brown ctx.dt = 0.01; 133c4762a1bSJed Brown ierr = PetscOptionsScalar("-Tf","","",ctx.Tf,&ctx.Tf,NULL);CHKERRQ(ierr); 134c4762a1bSJed Brown ierr = PetscOptionsScalar("-dt","","",ctx.dt,&ctx.dt,NULL);CHKERRQ(ierr); 135c4762a1bSJed Brown } 136c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 139c4762a1bSJed Brown Initialize the solution 140c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 141c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 142c4762a1bSJed Brown u[0] = 1.0/ctx.a; 143c4762a1bSJed Brown u[1] = 0.0; 144c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 145c4762a1bSJed Brown 146c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 147c4762a1bSJed Brown Create timestepping solver context 148c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 149c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 150c4762a1bSJed Brown ierr = TSSetType(ts,TSMPRK);CHKERRQ(ierr); 151c4762a1bSJed Brown 152c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,(TSRHSFunction)RHSFunction,&ctx);CHKERRQ(ierr); 153c4762a1bSJed Brown ierr = TSRHSSplitSetIS(ts,"slow",iss);CHKERRQ(ierr); 154c4762a1bSJed Brown ierr = TSRHSSplitSetIS(ts,"fast",isf);CHKERRQ(ierr); 155109dc152SHong Zhang ierr = TSRHSSplitSetRHSFunction(ts,"slow",NULL,(TSRHSFunction)RHSFunctionslow,&ctx);CHKERRQ(ierr); 156109dc152SHong Zhang ierr = TSRHSSplitSetRHSFunction(ts,"fast",NULL,(TSRHSFunction)RHSFunctionfast,&ctx);CHKERRQ(ierr); 157c4762a1bSJed Brown 158c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 159c4762a1bSJed Brown Set initial conditions 160c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 161c4762a1bSJed Brown ierr = TSSetSolution(ts,U);CHKERRQ(ierr); 162c4762a1bSJed Brown 163c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 164c4762a1bSJed Brown Set solver options 165c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 166c4762a1bSJed Brown ierr = TSSetMaxTime(ts,ctx.Tf);CHKERRQ(ierr); 167c4762a1bSJed Brown ierr = TSSetTimeStep(ts,ctx.dt);CHKERRQ(ierr); 168c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 169c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 170c4762a1bSJed Brown 171c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 172c4762a1bSJed Brown Solve linear system 173c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 174c4762a1bSJed Brown ierr = TSSolve(ts,U);CHKERRQ(ierr); 175c4762a1bSJed Brown ierr = VecView(U,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 176c4762a1bSJed Brown 177c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 178c4762a1bSJed Brown Check the error of the Petsc solution 179c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 180c4762a1bSJed Brown ierr = TSGetTime(ts,&tt);CHKERRQ(ierr); 181c4762a1bSJed Brown ierr = sol_true(tt,Utrue,&ctx);CHKERRQ(ierr); 182c4762a1bSJed Brown ierr = VecAXPY(Utrue,-1.0,U);CHKERRQ(ierr); 183c4762a1bSJed Brown ierr = VecNorm(Utrue,NORM_2,&error); 184c4762a1bSJed Brown 185c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 186c4762a1bSJed Brown Print norm2 error 187c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 188c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"l2 error norm: %g\n", error);CHKERRQ(ierr); 189c4762a1bSJed Brown 190c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 191c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 192c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 193c4762a1bSJed Brown ierr = VecDestroy(&U);CHKERRQ(ierr); 194c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 195c4762a1bSJed Brown ierr = VecDestroy(&Utrue);CHKERRQ(ierr); 196c4762a1bSJed Brown ierr = ISDestroy(&iss);CHKERRQ(ierr); 197c4762a1bSJed Brown ierr = ISDestroy(&isf);CHKERRQ(ierr); 198c4762a1bSJed Brown ierr = PetscFree(indicess);CHKERRQ(ierr); 199c4762a1bSJed Brown ierr = PetscFree(indicesf);CHKERRQ(ierr); 200c4762a1bSJed Brown ierr = PetscFinalize(); 201c4762a1bSJed Brown return ierr; 202c4762a1bSJed Brown } 203109dc152SHong Zhang 204109dc152SHong Zhang /*TEST 205109dc152SHong Zhang build: 206*f56ea12dSJed Brown requires: !complex 207109dc152SHong Zhang 208109dc152SHong Zhang test: 209109dc152SHong Zhang 210109dc152SHong Zhang TEST*/ 211