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 18d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx) 19d71ae5a4SJacob Faibussowitsch { 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] = -PetscSinScalar(ctx->a * t); 27c4762a1bSJed Brown f[1] = ctx->b * PetscCosScalar(ctx->b * t) * u[0] - PetscSinScalar(ctx->b * t) * PetscSinScalar(ctx->a * t); 289566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 299566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 31c4762a1bSJed Brown } 32c4762a1bSJed Brown 33d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSFunctionslow(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx) 34d71ae5a4SJacob Faibussowitsch { 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] = -PetscSinScalar(ctx->a * t); 429566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 439566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 45c4762a1bSJed Brown } 46c4762a1bSJed Brown 47d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSFunctionfast(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx) 48d71ae5a4SJacob Faibussowitsch { 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] = ctx->b * PetscCosScalar(ctx->b * t) * u[0] - PetscSinScalar(ctx->b * t) * PetscSinScalar(ctx->a * t); 569566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 579566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 59c4762a1bSJed Brown } 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* 62c4762a1bSJed Brown Define the analytic solution for check method easily 63c4762a1bSJed Brown */ 64d71ae5a4SJacob Faibussowitsch static PetscErrorCode sol_true(PetscReal t, Vec U, AppCtx *ctx) 65d71ae5a4SJacob Faibussowitsch { 66c4762a1bSJed Brown PetscScalar *u; 67c4762a1bSJed Brown 68c4762a1bSJed Brown PetscFunctionBegin; 699566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 70c4762a1bSJed Brown u[0] = PetscCosScalar(ctx->a * t) / ctx->a; 71c4762a1bSJed Brown u[1] = PetscSinScalar(ctx->b * t) * PetscCosScalar(ctx->a * t) / ctx->a; 729566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 74c4762a1bSJed Brown } 75c4762a1bSJed Brown 76d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 77d71ae5a4SJacob Faibussowitsch { 78c4762a1bSJed Brown TS ts; /* ODE integrator */ 79c4762a1bSJed Brown Vec U; /* solution will be stored here */ 80c4762a1bSJed Brown Vec Utrue; 81c4762a1bSJed Brown PetscMPIInt size; 82c4762a1bSJed Brown AppCtx ctx; 83c4762a1bSJed Brown PetscScalar *u; 84c4762a1bSJed Brown IS iss; 85c4762a1bSJed Brown IS isf; 86c4762a1bSJed Brown PetscInt *indicess; 87c4762a1bSJed Brown PetscInt *indicesf; 88c4762a1bSJed Brown PetscInt n = 2; 89c4762a1bSJed Brown PetscScalar error, tt; 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 92c4762a1bSJed Brown Initialize program 93c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 94327415f7SBarry Smith PetscFunctionBeginUser; 95c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 973c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs"); 98c4762a1bSJed Brown 99c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 100c4762a1bSJed Brown Create index for slow part and fast part 101c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1, &indicess)); 103c4762a1bSJed Brown indicess[0] = 0; 1049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1, &indicesf)); 105c4762a1bSJed Brown indicesf[0] = 1; 1069566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 1, indicess, PETSC_COPY_VALUES, &iss)); 1079566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 1, indicesf, PETSC_COPY_VALUES, &isf)); 108c4762a1bSJed Brown 109c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 1106aad120cSJose E. Roman Create necessary vector 111c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1129566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &U)); 1139566063dSJacob Faibussowitsch PetscCall(VecSetSizes(U, n, PETSC_DETERMINE)); 1149566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(U)); 1159566063dSJacob Faibussowitsch PetscCall(VecDuplicate(U, &Utrue)); 1169566063dSJacob Faibussowitsch PetscCall(VecCopy(U, Utrue)); 117c4762a1bSJed Brown 118c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 119c4762a1bSJed Brown Set runtime options 120c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 121d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ODE options", ""); 122c4762a1bSJed Brown { 123c4762a1bSJed Brown ctx.a = 1.0; 124c4762a1bSJed Brown ctx.b = 25.0; 1259566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-a", "", "", ctx.a, &ctx.a, NULL)); 1269566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-b", "", "", ctx.b, &ctx.b, NULL)); 127c4762a1bSJed Brown ctx.Tf = 5.0; 128c4762a1bSJed Brown ctx.dt = 0.01; 1299566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Tf", "", "", ctx.Tf, &ctx.Tf, NULL)); 1309566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-dt", "", "", ctx.dt, &ctx.dt, NULL)); 131c4762a1bSJed Brown } 132d0609cedSBarry Smith PetscOptionsEnd(); 133c4762a1bSJed Brown 134c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 135c4762a1bSJed Brown Initialize the solution 136c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1379566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 138c4762a1bSJed Brown u[0] = 1.0 / ctx.a; 139c4762a1bSJed Brown u[1] = 0.0; 1409566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 141c4762a1bSJed Brown 142c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 143c4762a1bSJed Brown Create timestepping solver context 144c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1459566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1469566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSMPRK)); 147c4762a1bSJed Brown 1488434afd1SBarry Smith PetscCall(TSSetRHSFunction(ts, NULL, (TSRHSFunctionFn *)RHSFunction, &ctx)); 1499566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts, "slow", iss)); 1509566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts, "fast", isf)); 1518434afd1SBarry Smith PetscCall(TSRHSSplitSetRHSFunction(ts, "slow", NULL, (TSRHSFunctionFn *)RHSFunctionslow, &ctx)); 1528434afd1SBarry Smith PetscCall(TSRHSSplitSetRHSFunction(ts, "fast", NULL, (TSRHSFunctionFn *)RHSFunctionfast, &ctx)); 153c4762a1bSJed Brown 154c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 155c4762a1bSJed Brown Set initial conditions 156c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1579566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, U)); 158c4762a1bSJed Brown 159c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 160c4762a1bSJed Brown Set solver options 161c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1629566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, ctx.Tf)); 1639566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, ctx.dt)); 1649566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 1659566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 166c4762a1bSJed Brown 167c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 168c4762a1bSJed Brown Solve linear system 169c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1709566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 1719566063dSJacob Faibussowitsch PetscCall(VecView(U, PETSC_VIEWER_STDOUT_WORLD)); 172c4762a1bSJed Brown 173c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 174*f0b74427SPierre Jolivet Check the error of the PETSc solution 175c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1769566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &tt)); 1779566063dSJacob Faibussowitsch PetscCall(sol_true(tt, Utrue, &ctx)); 1789566063dSJacob Faibussowitsch PetscCall(VecAXPY(Utrue, -1.0, U)); 1799566063dSJacob Faibussowitsch PetscCall(VecNorm(Utrue, NORM_2, &error)); 180c4762a1bSJed Brown 181c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 182c4762a1bSJed Brown Print norm2 error 183c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 18463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "l2 error norm: %g\n", (double)PetscAbsScalar(error))); 185c4762a1bSJed Brown 186c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 187c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 188c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 1909566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Utrue)); 1929566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iss)); 1939566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isf)); 1949566063dSJacob Faibussowitsch PetscCall(PetscFree(indicess)); 1959566063dSJacob Faibussowitsch PetscCall(PetscFree(indicesf)); 1969566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 197b122ec5aSJacob Faibussowitsch return 0; 198c4762a1bSJed Brown } 199109dc152SHong Zhang 200109dc152SHong Zhang /*TEST 201109dc152SHong Zhang build: 202f56ea12dSJed Brown requires: !complex 203109dc152SHong Zhang 204109dc152SHong Zhang test: 205109dc152SHong Zhang 206109dc152SHong Zhang TEST*/ 207