1c4762a1bSJed Brown /* 2c4762a1bSJed Brown Formatted test for TS routines. 3c4762a1bSJed Brown 4c4762a1bSJed Brown Solves U_t=F(t,u) 5c4762a1bSJed Brown Where: 6c4762a1bSJed Brown 7c04c6898SStefano Zampini [2*u1+u2 ] 8c04c6898SStefano Zampini F(t,u)= [u1+2*u2+u3] 9c04c6898SStefano Zampini [ u2+2*u3] 10c4762a1bSJed Brown 118c797eb4SStefano Zampini When run in parallel, each process solves the same set of equations separately. 12c4762a1bSJed Brown */ 13c4762a1bSJed Brown 14c4762a1bSJed Brown static char help[] = "Solves a linear ODE. \n\n"; 15c4762a1bSJed Brown 16c4762a1bSJed Brown #include <petscts.h> 17c4762a1bSJed Brown #include <petscpc.h> 18c4762a1bSJed Brown 19c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS, PetscReal, Vec, Vec, void *); 20c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS, PetscReal, Vec, Mat, Mat, void *); 21c4762a1bSJed Brown extern PetscErrorCode Monitor(TS, PetscInt, PetscReal, Vec, void *); 22c4762a1bSJed Brown extern PetscErrorCode Initial(Vec, void *); 23c4762a1bSJed Brown extern PetscErrorCode MyMatMult(Mat, Vec, Vec); 24c4762a1bSJed Brown 25c4762a1bSJed Brown extern PetscReal solx(PetscReal); 26c4762a1bSJed Brown extern PetscReal soly(PetscReal); 27c4762a1bSJed Brown extern PetscReal solz(PetscReal); 28c4762a1bSJed Brown 29d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 30d71ae5a4SJacob Faibussowitsch { 31c4762a1bSJed Brown PetscInt time_steps = 100, steps; 32c4762a1bSJed Brown Vec global; 33c4762a1bSJed Brown PetscReal dt, ftime; 34c4762a1bSJed Brown TS ts; 358c797eb4SStefano Zampini Mat A, S; 368c797eb4SStefano Zampini PetscBool nest = PETSC_FALSE; 37c4762a1bSJed Brown 38327415f7SBarry Smith PetscFunctionBeginUser; 39c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-time", &time_steps, NULL)); 41eaf6e66dSStefano Zampini PetscCall(PetscOptionsGetBool(NULL, NULL, "-nest", &nest, NULL)); 42c4762a1bSJed Brown 43eaf6e66dSStefano Zampini /* create vector to hold state */ 44eaf6e66dSStefano Zampini if (nest) { 45eaf6e66dSStefano Zampini Vec g[3]; 46eaf6e66dSStefano Zampini 47eaf6e66dSStefano Zampini PetscCall(VecCreate(PETSC_COMM_WORLD, &g[0])); 488c797eb4SStefano Zampini PetscCall(VecSetSizes(g[0], 1, PETSC_DECIDE)); 49eaf6e66dSStefano Zampini PetscCall(VecSetFromOptions(g[0])); 50eaf6e66dSStefano Zampini PetscCall(VecDuplicate(g[0], &g[1])); 51eaf6e66dSStefano Zampini PetscCall(VecDuplicate(g[0], &g[2])); 52eaf6e66dSStefano Zampini PetscCall(VecCreateNest(PETSC_COMM_WORLD, 3, NULL, g, &global)); 53eaf6e66dSStefano Zampini PetscCall(VecDestroy(&g[0])); 54eaf6e66dSStefano Zampini PetscCall(VecDestroy(&g[1])); 55eaf6e66dSStefano Zampini PetscCall(VecDestroy(&g[2])); 56eaf6e66dSStefano Zampini } else { 579566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &global)); 588c797eb4SStefano Zampini PetscCall(VecSetSizes(global, 3, PETSC_DECIDE)); 599566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(global)); 60eaf6e66dSStefano Zampini } 61eaf6e66dSStefano Zampini 62eaf6e66dSStefano Zampini /* set initial conditions */ 639566063dSJacob Faibussowitsch PetscCall(Initial(global, NULL)); 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* make timestep context */ 669566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 679566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 689566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts, Monitor, NULL, NULL)); 69c4762a1bSJed Brown dt = 0.001; 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* 72c4762a1bSJed Brown The user provides the RHS and Jacobian 73c4762a1bSJed Brown */ 749566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, NULL)); 759566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 768c797eb4SStefano Zampini PetscCall(MatSetSizes(A, 3, 3, PETSC_DECIDE, PETSC_DECIDE)); 779566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 789566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 799566063dSJacob Faibussowitsch PetscCall(RHSJacobian(ts, 0.0, global, A, A, NULL)); 809566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, A, A, RHSJacobian, NULL)); 81c4762a1bSJed Brown 828c797eb4SStefano Zampini PetscCall(MatCreateShell(PETSC_COMM_WORLD, 3, 3, PETSC_DECIDE, PETSC_DECIDE, NULL, &S)); 8357d50842SBarry Smith PetscCall(MatShellSetOperation(S, MATOP_MULT, (PetscErrorCodeFn *)MyMatMult)); 849566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, S, A, RHSJacobian, NULL)); 85c4762a1bSJed Brown 869566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 879566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 88c4762a1bSJed Brown 899566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt)); 909566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, time_steps)); 919566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 1)); 929566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, global)); 93c4762a1bSJed Brown 949566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, global)); 959566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 969566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 97c4762a1bSJed Brown 988c797eb4SStefano Zampini /* free the memory */ 999566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&global)); 1019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1029566063dSJacob Faibussowitsch PetscCall(MatDestroy(&S)); 103c4762a1bSJed Brown 1049566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 105b122ec5aSJacob Faibussowitsch return 0; 106c4762a1bSJed Brown } 107c4762a1bSJed Brown 108d71ae5a4SJacob Faibussowitsch PetscErrorCode MyMatMult(Mat S, Vec x, Vec y) 109d71ae5a4SJacob Faibussowitsch { 110c4762a1bSJed Brown const PetscScalar *inptr; 111c4762a1bSJed Brown PetscScalar *outptr; 112c4762a1bSJed Brown 1137510d9b0SBarry Smith PetscFunctionBeginUser; 1149566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &inptr)); 1159566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(y, &outptr)); 116c4762a1bSJed Brown 117c4762a1bSJed Brown outptr[0] = 2.0 * inptr[0] + inptr[1]; 118c4762a1bSJed Brown outptr[1] = inptr[0] + 2.0 * inptr[1] + inptr[2]; 119c4762a1bSJed Brown outptr[2] = inptr[1] + 2.0 * inptr[2]; 120c4762a1bSJed Brown 1219566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &inptr)); 1229566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(y, &outptr)); 1233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126*2a8381b2SBarry Smith PetscErrorCode Initial(Vec global, PetscCtx ctx) 127d71ae5a4SJacob Faibussowitsch { 128c4762a1bSJed Brown PetscScalar *localptr; 129c4762a1bSJed Brown 1303ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 1319566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(global, &localptr)); 1328c797eb4SStefano Zampini localptr[0] = solx(0.0); 1338c797eb4SStefano Zampini localptr[1] = soly(0.0); 1348c797eb4SStefano Zampini localptr[2] = solz(0.0); 1359566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(global, &localptr)); 1363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 137c4762a1bSJed Brown } 138c4762a1bSJed Brown 139*2a8381b2SBarry Smith PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal time, Vec global, PetscCtx ctx) 140d71ae5a4SJacob Faibussowitsch { 141c4762a1bSJed Brown const PetscScalar *tmp; 1428c797eb4SStefano Zampini PetscScalar exact[] = {solx(time), soly(time), solz(time)}; 143c4762a1bSJed Brown 1443ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 1458c797eb4SStefano Zampini PetscCall(VecGetArrayRead(global, &tmp)); 1469371c9d4SSatish Balay PetscCall(PetscPrintf(PETSC_COMM_WORLD, "At t =%14.6e u = %14.6e %14.6e %14.6e \n", (double)time, (double)PetscRealPart(tmp[0]), (double)PetscRealPart(tmp[1]), (double)PetscRealPart(tmp[2]))); 1478c797eb4SStefano Zampini PetscCall(PetscPrintf(PETSC_COMM_WORLD, "At t =%14.6e errors = %14.6e %14.6e %14.6e \n", (double)time, (double)PetscRealPart(tmp[0] - exact[0]), (double)PetscRealPart(tmp[1] - exact[1]), (double)PetscRealPart(tmp[2] - exact[2]))); 1488c797eb4SStefano Zampini PetscCall(VecRestoreArrayRead(global, &tmp)); 1493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 150c4762a1bSJed Brown } 151c4762a1bSJed Brown 152*2a8381b2SBarry Smith PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec globalin, Vec globalout, PetscCtx ctx) 153d71ae5a4SJacob Faibussowitsch { 154c4762a1bSJed Brown PetscScalar *outptr; 155c4762a1bSJed Brown const PetscScalar *inptr; 156c4762a1bSJed Brown 1573ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 158c4762a1bSJed Brown /*Extract income array */ 1598c797eb4SStefano Zampini PetscCall(VecGetArrayRead(globalin, &inptr)); 160c4762a1bSJed Brown 161c4762a1bSJed Brown /* Extract outcome array*/ 1628c797eb4SStefano Zampini PetscCall(VecGetArrayWrite(globalout, &outptr)); 163c4762a1bSJed Brown 164c4762a1bSJed Brown outptr[0] = 2.0 * inptr[0] + inptr[1]; 165c4762a1bSJed Brown outptr[1] = inptr[0] + 2.0 * inptr[1] + inptr[2]; 166c4762a1bSJed Brown outptr[2] = inptr[1] + 2.0 * inptr[2]; 167c4762a1bSJed Brown 1688c797eb4SStefano Zampini PetscCall(VecRestoreArrayRead(globalin, &inptr)); 1698c797eb4SStefano Zampini PetscCall(VecRestoreArrayWrite(globalout, &outptr)); 1703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 171c4762a1bSJed Brown } 172c4762a1bSJed Brown 173*2a8381b2SBarry Smith PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec x, Mat A, Mat BB, PetscCtx ctx) 174d71ae5a4SJacob Faibussowitsch { 175c4762a1bSJed Brown PetscScalar v[3]; 1768c797eb4SStefano Zampini PetscInt idx[3], rst; 177c4762a1bSJed Brown 1783ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 1798c797eb4SStefano Zampini PetscCall(VecGetOwnershipRange(x, &rst, NULL)); 1808c797eb4SStefano Zampini idx[0] = 0 + rst; 1818c797eb4SStefano Zampini idx[1] = 1 + rst; 1828c797eb4SStefano Zampini idx[2] = 2 + rst; 183c4762a1bSJed Brown 1849371c9d4SSatish Balay v[0] = 2.0; 1859371c9d4SSatish Balay v[1] = 1.0; 1869371c9d4SSatish Balay v[2] = 0.0; 1878c797eb4SStefano Zampini PetscCall(MatSetValues(BB, 1, idx, 3, idx, v, INSERT_VALUES)); 188c4762a1bSJed Brown 1899371c9d4SSatish Balay v[0] = 1.0; 1909371c9d4SSatish Balay v[1] = 2.0; 1919371c9d4SSatish Balay v[2] = 1.0; 1928c797eb4SStefano Zampini PetscCall(MatSetValues(BB, 1, idx + 1, 3, idx, v, INSERT_VALUES)); 193c4762a1bSJed Brown 1949371c9d4SSatish Balay v[0] = 0.0; 1959371c9d4SSatish Balay v[1] = 1.0; 1969371c9d4SSatish Balay v[2] = 2.0; 1978c797eb4SStefano Zampini PetscCall(MatSetValues(BB, 1, idx + 2, 3, idx, v, INSERT_VALUES)); 198c4762a1bSJed Brown 1999566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(BB, MAT_FINAL_ASSEMBLY)); 2009566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(BB, MAT_FINAL_ASSEMBLY)); 201c4762a1bSJed Brown 202c4762a1bSJed Brown if (A != BB) { 2039566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 2049566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 205c4762a1bSJed Brown } 2063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 207c4762a1bSJed Brown } 208c4762a1bSJed Brown 209c4762a1bSJed Brown /* 210c4762a1bSJed Brown The exact solutions 211c4762a1bSJed Brown */ 212d71ae5a4SJacob Faibussowitsch PetscReal solx(PetscReal t) 213d71ae5a4SJacob Faibussowitsch { 2149371c9d4SSatish Balay return PetscExpReal((2.0 - PetscSqrtReal(2.0)) * t) / 2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0)) * t) / (2.0 * PetscSqrtReal(2.0)) + PetscExpReal((2.0 + PetscSqrtReal(2.0)) * t) / 2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0)) * t) / (2.0 * PetscSqrtReal(2.0)); 215c4762a1bSJed Brown } 216c4762a1bSJed Brown 217d71ae5a4SJacob Faibussowitsch PetscReal soly(PetscReal t) 218d71ae5a4SJacob Faibussowitsch { 2199371c9d4SSatish Balay return PetscExpReal((2.0 - PetscSqrtReal(2.0)) * t) / 2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0)) * t) / PetscSqrtReal(2.0) + PetscExpReal((2.0 + PetscSqrtReal(2.0)) * t) / 2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0)) * t) / PetscSqrtReal(2.0); 220c4762a1bSJed Brown } 221c4762a1bSJed Brown 222d71ae5a4SJacob Faibussowitsch PetscReal solz(PetscReal t) 223d71ae5a4SJacob Faibussowitsch { 2249371c9d4SSatish Balay return PetscExpReal((2.0 - PetscSqrtReal(2.0)) * t) / 2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0)) * t) / (2.0 * PetscSqrtReal(2.0)) + PetscExpReal((2.0 + PetscSqrtReal(2.0)) * t) / 2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0)) * t) / (2.0 * PetscSqrtReal(2.0)); 225c4762a1bSJed Brown } 226c4762a1bSJed Brown 227c4762a1bSJed Brown /*TEST 228c4762a1bSJed Brown 229c4762a1bSJed Brown test: 230c4762a1bSJed Brown suffix: euler 231eaf6e66dSStefano Zampini args: -ts_type euler -nest {{0 1}} 232c4762a1bSJed Brown requires: !single 233c4762a1bSJed Brown 234c4762a1bSJed Brown test: 235c4762a1bSJed Brown suffix: beuler 236eaf6e66dSStefano Zampini args: -ts_type beuler -nest {{0 1}} 237eaf6e66dSStefano Zampini requires: !single 238eaf6e66dSStefano Zampini 239eaf6e66dSStefano Zampini test: 240eaf6e66dSStefano Zampini suffix: rk 241eaf6e66dSStefano Zampini args: -ts_type rk -nest {{0 1}} -ts_adapt_monitor 242c4762a1bSJed Brown requires: !single 243c4762a1bSJed Brown 2448c797eb4SStefano Zampini test: 2458c797eb4SStefano Zampini diff_args: -j 2468c797eb4SStefano Zampini requires: double !complex 2478c797eb4SStefano Zampini output_file: output/ex2_be_adapt.out 2488c797eb4SStefano Zampini suffix: bdf_1_adapt 2498c797eb4SStefano Zampini args: -ts_type bdf -ts_bdf_order 1 -ts_adapt_type basic -ts_adapt_clip 0,2 2508c797eb4SStefano Zampini 2518c797eb4SStefano Zampini test: 2528c797eb4SStefano Zampini diff_args: -j 2538c797eb4SStefano Zampini requires: double !complex 2548c797eb4SStefano Zampini suffix: be_adapt 2558c797eb4SStefano Zampini args: -ts_type beuler -ts_adapt_type basic -ts_adapt_clip 0,2 2568c797eb4SStefano Zampini 257c4762a1bSJed Brown TEST*/ 258