1c4762a1bSJed Brown /* 2c4762a1bSJed Brown Formatted test for TS routines. 3c4762a1bSJed Brown 4c4762a1bSJed Brown Solves U_t=F(t,u) 5c4762a1bSJed Brown Where: 6c4762a1bSJed Brown 7c4762a1bSJed Brown [2*u1+u2 8c4762a1bSJed Brown F(t,u)= [u1+2*u2+u3 9c4762a1bSJed Brown [ u2+2*u3 10c4762a1bSJed Brown We can compare the solutions from euler, beuler and SUNDIALS to 11c4762a1bSJed Brown see what is the difference. 12c4762a1bSJed Brown 13c4762a1bSJed Brown */ 14c4762a1bSJed Brown 15c4762a1bSJed Brown static char help[] = "Solves a linear ODE. \n\n"; 16c4762a1bSJed Brown 17c4762a1bSJed Brown #include <petscts.h> 18c4762a1bSJed Brown #include <petscpc.h> 19c4762a1bSJed Brown 20c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS, PetscReal, Vec, Vec, void *); 21c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS, PetscReal, Vec, Mat, Mat, void *); 22c4762a1bSJed Brown extern PetscErrorCode Monitor(TS, PetscInt, PetscReal, Vec, void *); 23c4762a1bSJed Brown extern PetscErrorCode Initial(Vec, void *); 24c4762a1bSJed Brown extern PetscErrorCode MyMatMult(Mat, Vec, Vec); 25c4762a1bSJed Brown 26c4762a1bSJed Brown extern PetscReal solx(PetscReal); 27c4762a1bSJed Brown extern PetscReal soly(PetscReal); 28c4762a1bSJed Brown extern PetscReal solz(PetscReal); 29c4762a1bSJed Brown 30*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 31*d71ae5a4SJacob Faibussowitsch { 32c4762a1bSJed Brown PetscInt time_steps = 100, steps; 33c4762a1bSJed Brown Vec global; 34c4762a1bSJed Brown PetscReal dt, ftime; 35c4762a1bSJed Brown TS ts; 36c4762a1bSJed Brown Mat A = 0, S; 37c4762a1bSJed Brown 38327415f7SBarry Smith PetscFunctionBeginUser; 399566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-time", &time_steps, NULL)); 41c4762a1bSJed Brown 42c4762a1bSJed Brown /* set initial conditions */ 439566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &global)); 449566063dSJacob Faibussowitsch PetscCall(VecSetSizes(global, PETSC_DECIDE, 3)); 459566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(global)); 469566063dSJacob Faibussowitsch PetscCall(Initial(global, NULL)); 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* make timestep context */ 499566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 509566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 519566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts, Monitor, NULL, NULL)); 52c4762a1bSJed Brown dt = 0.001; 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* 55c4762a1bSJed Brown The user provides the RHS and Jacobian 56c4762a1bSJed Brown */ 579566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, NULL)); 589566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 599566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 3, 3)); 609566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 619566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 629566063dSJacob Faibussowitsch PetscCall(RHSJacobian(ts, 0.0, global, A, A, NULL)); 639566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, A, A, RHSJacobian, NULL)); 64c4762a1bSJed Brown 659566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, 3, 3, 3, 3, NULL, &S)); 669566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(S, MATOP_MULT, (void (*)(void))MyMatMult)); 679566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, S, A, RHSJacobian, NULL)); 68c4762a1bSJed Brown 699566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 709566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 71c4762a1bSJed Brown 729566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt)); 739566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, time_steps)); 749566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 1)); 759566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, global)); 76c4762a1bSJed Brown 779566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, global)); 789566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 799566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* free the memories */ 82c4762a1bSJed Brown 839566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&global)); 859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 869566063dSJacob Faibussowitsch PetscCall(MatDestroy(&S)); 87c4762a1bSJed Brown 889566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 89b122ec5aSJacob Faibussowitsch return 0; 90c4762a1bSJed Brown } 91c4762a1bSJed Brown 92*d71ae5a4SJacob Faibussowitsch PetscErrorCode MyMatMult(Mat S, Vec x, Vec y) 93*d71ae5a4SJacob Faibussowitsch { 94c4762a1bSJed Brown const PetscScalar *inptr; 95c4762a1bSJed Brown PetscScalar *outptr; 96c4762a1bSJed Brown 977510d9b0SBarry Smith PetscFunctionBeginUser; 989566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &inptr)); 999566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(y, &outptr)); 100c4762a1bSJed Brown 101c4762a1bSJed Brown outptr[0] = 2.0 * inptr[0] + inptr[1]; 102c4762a1bSJed Brown outptr[1] = inptr[0] + 2.0 * inptr[1] + inptr[2]; 103c4762a1bSJed Brown outptr[2] = inptr[1] + 2.0 * inptr[2]; 104c4762a1bSJed Brown 1059566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &inptr)); 1069566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(y, &outptr)); 107c4762a1bSJed Brown PetscFunctionReturn(0); 108c4762a1bSJed Brown } 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* this test problem has initial values (1,1,1). */ 111*d71ae5a4SJacob Faibussowitsch PetscErrorCode Initial(Vec global, void *ctx) 112*d71ae5a4SJacob Faibussowitsch { 113c4762a1bSJed Brown PetscScalar *localptr; 114c4762a1bSJed Brown PetscInt i, mybase, myend, locsize; 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* determine starting point of each processor */ 1179566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(global, &mybase, &myend)); 1189566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(global, &locsize)); 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* Initialize the array */ 1219566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(global, &localptr)); 122c4762a1bSJed Brown for (i = 0; i < locsize; i++) localptr[i] = 1.0; 123c4762a1bSJed Brown 124c4762a1bSJed Brown if (mybase == 0) localptr[0] = 1.0; 125c4762a1bSJed Brown 1269566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(global, &localptr)); 127c4762a1bSJed Brown return 0; 128c4762a1bSJed Brown } 129c4762a1bSJed Brown 130*d71ae5a4SJacob Faibussowitsch PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal time, Vec global, void *ctx) 131*d71ae5a4SJacob Faibussowitsch { 132c4762a1bSJed Brown VecScatter scatter; 133c4762a1bSJed Brown IS from, to; 134c4762a1bSJed Brown PetscInt i, n, *idx; 135c4762a1bSJed Brown Vec tmp_vec; 136c4762a1bSJed Brown const PetscScalar *tmp; 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* Get the size of the vector */ 1399566063dSJacob Faibussowitsch PetscCall(VecGetSize(global, &n)); 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* Set the index sets */ 1429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &idx)); 143c4762a1bSJed Brown for (i = 0; i < n; i++) idx[i] = i; 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* Create local sequential vectors */ 1469566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &tmp_vec)); 147c4762a1bSJed Brown 148c4762a1bSJed Brown /* Create scatter context */ 1499566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, idx, PETSC_COPY_VALUES, &from)); 1509566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, idx, PETSC_COPY_VALUES, &to)); 1519566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(global, from, tmp_vec, to, &scatter)); 1529566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scatter, global, tmp_vec, INSERT_VALUES, SCATTER_FORWARD)); 1539566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scatter, global, tmp_vec, INSERT_VALUES, SCATTER_FORWARD)); 154c4762a1bSJed Brown 1559566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(tmp_vec, &tmp)); 1569371c9d4SSatish 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]))); 1579371c9d4SSatish Balay PetscCall(PetscPrintf(PETSC_COMM_WORLD, "At t =%14.6e errors = %14.6e %14.6e %14.6e \n", (double)time, (double)PetscRealPart(tmp[0] - solx(time)), (double)PetscRealPart(tmp[1] - soly(time)), (double)PetscRealPart(tmp[2] - solz(time)))); 1589566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(tmp_vec, &tmp)); 1599566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&scatter)); 1609566063dSJacob Faibussowitsch PetscCall(ISDestroy(&from)); 1619566063dSJacob Faibussowitsch PetscCall(ISDestroy(&to)); 1629566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 1639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tmp_vec)); 164c4762a1bSJed Brown return 0; 165c4762a1bSJed Brown } 166c4762a1bSJed Brown 167*d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec globalin, Vec globalout, void *ctx) 168*d71ae5a4SJacob Faibussowitsch { 169c4762a1bSJed Brown PetscScalar *outptr; 170c4762a1bSJed Brown const PetscScalar *inptr; 171c4762a1bSJed Brown PetscInt i, n, *idx; 172c4762a1bSJed Brown IS from, to; 173c4762a1bSJed Brown VecScatter scatter; 174c4762a1bSJed Brown Vec tmp_in, tmp_out; 175c4762a1bSJed Brown 176c4762a1bSJed Brown /* Get the length of parallel vector */ 1779566063dSJacob Faibussowitsch PetscCall(VecGetSize(globalin, &n)); 178c4762a1bSJed Brown 179c4762a1bSJed Brown /* Set the index sets */ 1809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &idx)); 181c4762a1bSJed Brown for (i = 0; i < n; i++) idx[i] = i; 182c4762a1bSJed Brown 183c4762a1bSJed Brown /* Create local sequential vectors */ 1849566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &tmp_in)); 1859566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tmp_in, &tmp_out)); 186c4762a1bSJed Brown 187c4762a1bSJed Brown /* Create scatter context */ 1889566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, idx, PETSC_COPY_VALUES, &from)); 1899566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, idx, PETSC_COPY_VALUES, &to)); 1909566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(globalin, from, tmp_in, to, &scatter)); 1919566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scatter, globalin, tmp_in, INSERT_VALUES, SCATTER_FORWARD)); 1929566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scatter, globalin, tmp_in, INSERT_VALUES, SCATTER_FORWARD)); 1939566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&scatter)); 194c4762a1bSJed Brown 195c4762a1bSJed Brown /*Extract income array */ 1969566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(tmp_in, &inptr)); 197c4762a1bSJed Brown 198c4762a1bSJed Brown /* Extract outcome array*/ 1999566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(tmp_out, &outptr)); 200c4762a1bSJed Brown 201c4762a1bSJed Brown outptr[0] = 2.0 * inptr[0] + inptr[1]; 202c4762a1bSJed Brown outptr[1] = inptr[0] + 2.0 * inptr[1] + inptr[2]; 203c4762a1bSJed Brown outptr[2] = inptr[1] + 2.0 * inptr[2]; 204c4762a1bSJed Brown 2059566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(tmp_in, &inptr)); 2069566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(tmp_out, &outptr)); 207c4762a1bSJed Brown 2089566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tmp_out, from, globalout, to, &scatter)); 2099566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scatter, tmp_out, globalout, INSERT_VALUES, SCATTER_FORWARD)); 2109566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scatter, tmp_out, globalout, INSERT_VALUES, SCATTER_FORWARD)); 211c4762a1bSJed Brown 212c4762a1bSJed Brown /* Destroy idx aand scatter */ 2139566063dSJacob Faibussowitsch PetscCall(ISDestroy(&from)); 2149566063dSJacob Faibussowitsch PetscCall(ISDestroy(&to)); 2159566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&scatter)); 2169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tmp_in)); 2179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tmp_out)); 2189566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 219c4762a1bSJed Brown return 0; 220c4762a1bSJed Brown } 221c4762a1bSJed Brown 222*d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec x, Mat A, Mat BB, void *ctx) 223*d71ae5a4SJacob Faibussowitsch { 224c4762a1bSJed Brown PetscScalar v[3]; 225c4762a1bSJed Brown const PetscScalar *tmp; 226c4762a1bSJed Brown PetscInt idx[3], i; 227c4762a1bSJed Brown 2289371c9d4SSatish Balay idx[0] = 0; 2299371c9d4SSatish Balay idx[1] = 1; 2309371c9d4SSatish Balay idx[2] = 2; 2319566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &tmp)); 232c4762a1bSJed Brown 233c4762a1bSJed Brown i = 0; 2349371c9d4SSatish Balay v[0] = 2.0; 2359371c9d4SSatish Balay v[1] = 1.0; 2369371c9d4SSatish Balay v[2] = 0.0; 2379566063dSJacob Faibussowitsch PetscCall(MatSetValues(BB, 1, &i, 3, idx, v, INSERT_VALUES)); 238c4762a1bSJed Brown 239c4762a1bSJed Brown i = 1; 2409371c9d4SSatish Balay v[0] = 1.0; 2419371c9d4SSatish Balay v[1] = 2.0; 2429371c9d4SSatish Balay v[2] = 1.0; 2439566063dSJacob Faibussowitsch PetscCall(MatSetValues(BB, 1, &i, 3, idx, v, INSERT_VALUES)); 244c4762a1bSJed Brown 245c4762a1bSJed Brown i = 2; 2469371c9d4SSatish Balay v[0] = 0.0; 2479371c9d4SSatish Balay v[1] = 1.0; 2489371c9d4SSatish Balay v[2] = 2.0; 2499566063dSJacob Faibussowitsch PetscCall(MatSetValues(BB, 1, &i, 3, idx, v, INSERT_VALUES)); 250c4762a1bSJed Brown 2519566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(BB, MAT_FINAL_ASSEMBLY)); 2529566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(BB, MAT_FINAL_ASSEMBLY)); 253c4762a1bSJed Brown 254c4762a1bSJed Brown if (A != BB) { 2559566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 2569566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 257c4762a1bSJed Brown } 2589566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &tmp)); 259c4762a1bSJed Brown 260c4762a1bSJed Brown return 0; 261c4762a1bSJed Brown } 262c4762a1bSJed Brown 263c4762a1bSJed Brown /* 264c4762a1bSJed Brown The exact solutions 265c4762a1bSJed Brown */ 266*d71ae5a4SJacob Faibussowitsch PetscReal solx(PetscReal t) 267*d71ae5a4SJacob Faibussowitsch { 2689371c9d4SSatish 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)); 269c4762a1bSJed Brown } 270c4762a1bSJed Brown 271*d71ae5a4SJacob Faibussowitsch PetscReal soly(PetscReal t) 272*d71ae5a4SJacob Faibussowitsch { 2739371c9d4SSatish 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); 274c4762a1bSJed Brown } 275c4762a1bSJed Brown 276*d71ae5a4SJacob Faibussowitsch PetscReal solz(PetscReal t) 277*d71ae5a4SJacob Faibussowitsch { 2789371c9d4SSatish 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)); 279c4762a1bSJed Brown } 280c4762a1bSJed Brown 281c4762a1bSJed Brown /*TEST 282c4762a1bSJed Brown 283c4762a1bSJed Brown test: 284c4762a1bSJed Brown suffix: euler 285c4762a1bSJed Brown args: -ts_type euler 286c4762a1bSJed Brown requires: !single 287c4762a1bSJed Brown 288c4762a1bSJed Brown test: 289c4762a1bSJed Brown suffix: beuler 290c4762a1bSJed Brown args: -ts_type beuler 291c4762a1bSJed Brown requires: !single 292c4762a1bSJed Brown 293c4762a1bSJed Brown TEST*/ 294