1c4762a1bSJed Brown static char help[] = "Tests TSTrajectoryGetVecs. \n\n";
2c4762a1bSJed Brown /*
3c4762a1bSJed Brown This example tests TSTrajectory and the ability of TSTrajectoryGetVecs
4c4762a1bSJed Brown to reconstructs states and derivatives via interpolation (if necessary).
5c4762a1bSJed Brown It also tests TSTrajectory{Get|Restore}UpdatedHistoryVecs
6c4762a1bSJed Brown */
7c4762a1bSJed Brown #include <petscts.h>
8c4762a1bSJed Brown
func(PetscInt p,PetscReal t)9d71ae5a4SJacob Faibussowitsch PetscScalar func(PetscInt p, PetscReal t)
10d71ae5a4SJacob Faibussowitsch {
119371c9d4SSatish Balay return p ? t * func(p - 1, t) : 1.0;
129371c9d4SSatish Balay }
dfunc(PetscInt p,PetscReal t)13d71ae5a4SJacob Faibussowitsch PetscScalar dfunc(PetscInt p, PetscReal t)
14d71ae5a4SJacob Faibussowitsch {
159371c9d4SSatish Balay return p > 0 ? (PetscReal)p * func(p - 1, t) : 0.0;
169371c9d4SSatish Balay }
17c4762a1bSJed Brown
main(int argc,char ** argv)18d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
19d71ae5a4SJacob Faibussowitsch {
20c4762a1bSJed Brown TS ts;
21c4762a1bSJed Brown Vec W, W2, Wdot;
22c4762a1bSJed Brown TSTrajectory tj;
23c4762a1bSJed Brown PetscReal times[10], tol = PETSC_SMALL;
24c4762a1bSJed Brown PetscReal TT[10] = {2, 9, 1, 3, 6, 7, 5, 10, 4, 8};
25c4762a1bSJed Brown PetscInt i, p = 1, Nt = 10;
26c4762a1bSJed Brown PetscInt II[10] = {1, 4, 9, 2, 3, 6, 5, 8, 0, 7};
27c4762a1bSJed Brown PetscBool sort, use1, use2, check = PETSC_FALSE;
28c4762a1bSJed Brown
29327415f7SBarry Smith PetscFunctionBeginUser;
30*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
319566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &W));
329566063dSJacob Faibussowitsch PetscCall(VecSetSizes(W, 1, PETSC_DECIDE));
339566063dSJacob Faibussowitsch PetscCall(VecSetUp(W));
349566063dSJacob Faibussowitsch PetscCall(VecDuplicate(W, &Wdot));
359566063dSJacob Faibussowitsch PetscCall(VecDuplicate(W, &W2));
369566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
379566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, W2));
388337b394SStefano Zampini PetscCall(TSSetMaxTime(ts, 10.));
399566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, 10));
409566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts));
419566063dSJacob Faibussowitsch PetscCall(TSGetTrajectory(ts, &tj));
429566063dSJacob Faibussowitsch PetscCall(TSTrajectorySetType(tj, ts, TSTRAJECTORYBASIC));
439566063dSJacob Faibussowitsch PetscCall(TSTrajectorySetFromOptions(tj, ts));
449566063dSJacob Faibussowitsch PetscCall(TSTrajectorySetSolutionOnly(tj, PETSC_TRUE));
459566063dSJacob Faibussowitsch PetscCall(TSTrajectorySetUp(tj, ts));
469566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-check", &check, NULL));
479566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-p", &p, NULL));
489566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetRealArray(NULL, NULL, "-interptimes", times, &Nt, NULL));
49c4762a1bSJed Brown sort = PETSC_FALSE;
509566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-sorttimes", &sort, NULL));
511baa6e33SBarry Smith if (sort) PetscCall(PetscSortReal(10, TT));
52c4762a1bSJed Brown sort = PETSC_FALSE;
539566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-sortkeys", &sort, NULL));
541baa6e33SBarry Smith if (sort) PetscCall(PetscSortInt(10, II));
55c4762a1bSJed Brown p = PetscMax(p, -p);
56c4762a1bSJed Brown
57c4762a1bSJed Brown /* populate trajectory */
58c4762a1bSJed Brown for (i = 0; i < 10; i++) {
599566063dSJacob Faibussowitsch PetscCall(VecSet(W, func(p, TT[i])));
609566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, II[i]));
619566063dSJacob Faibussowitsch PetscCall(TSTrajectorySet(tj, ts, II[i], TT[i], W));
62c4762a1bSJed Brown }
63c4762a1bSJed Brown for (i = 0; i < Nt; i++) {
64c4762a1bSJed Brown PetscReal testtime = times[i], serr, derr;
65c4762a1bSJed Brown const PetscScalar *aW, *aWdot;
66c4762a1bSJed Brown
679566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGetVecs(tj, ts, PETSC_DECIDE, &testtime, W, Wdot));
689566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(W, &aW));
699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Wdot, &aWdot));
709566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " f(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(func(p, testtime)), (double)PetscRealPart(aW[0])));
719566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "df(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(dfunc(p, testtime)), (double)PetscRealPart(aWdot[0])));
72c4762a1bSJed Brown serr = PetscAbsScalar(func(p, testtime) - aW[0]);
73c4762a1bSJed Brown derr = PetscAbsScalar(dfunc(p, testtime) - aWdot[0]);
749566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(W, &aW));
759566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Wdot, &aWdot));
763c633725SBarry Smith PetscCheck(!check || serr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state %g > %g", (double)serr, (double)tol);
773c633725SBarry Smith PetscCheck(!check || derr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state derivative %g > %g", (double)derr, (double)tol);
78c4762a1bSJed Brown }
79c4762a1bSJed Brown for (i = Nt - 1; i >= 0; i--) {
80c4762a1bSJed Brown PetscReal testtime = times[i], serr;
81c4762a1bSJed Brown const PetscScalar *aW;
82c4762a1bSJed Brown
839566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGetVecs(tj, ts, PETSC_DECIDE, &testtime, W, NULL));
849566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(W, &aW));
859566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " f(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(func(p, testtime)), (double)PetscRealPart(aW[0])));
86c4762a1bSJed Brown serr = PetscAbsScalar(func(p, testtime) - aW[0]);
879566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(W, &aW));
883c633725SBarry Smith PetscCheck(!check || serr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state %g > %g", (double)serr, (double)tol);
89c4762a1bSJed Brown }
90c4762a1bSJed Brown for (i = Nt - 1; i >= 0; i--) {
91c4762a1bSJed Brown PetscReal testtime = times[i], derr;
92c4762a1bSJed Brown const PetscScalar *aWdot;
93c4762a1bSJed Brown
949566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGetVecs(tj, ts, PETSC_DECIDE, &testtime, NULL, Wdot));
959566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Wdot, &aWdot));
969566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "df(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(dfunc(p, testtime)), (double)PetscRealPart(aWdot[0])));
97c4762a1bSJed Brown derr = PetscAbsScalar(dfunc(p, testtime) - aWdot[0]);
989566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Wdot, &aWdot));
993c633725SBarry Smith PetscCheck(!check || derr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state derivative %g > %g", (double)derr, (double)tol);
100c4762a1bSJed Brown }
101c4762a1bSJed Brown for (i = 0; i < Nt; i++) {
102c4762a1bSJed Brown PetscReal testtime = times[i], serr, derr;
103c4762a1bSJed Brown const PetscScalar *aW, *aWdot;
104c4762a1bSJed Brown Vec hW, hWdot;
105c4762a1bSJed Brown
1069566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGetUpdatedHistoryVecs(tj, ts, testtime, &hW, &hWdot));
1079566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(hW, &aW));
1089566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(hWdot, &aWdot));
1099566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " f(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(func(p, testtime)), (double)PetscRealPart(aW[0])));
1109566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "df(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(dfunc(p, testtime)), (double)PetscRealPart(aWdot[0])));
111c4762a1bSJed Brown serr = PetscAbsScalar(func(p, testtime) - aW[0]);
112c4762a1bSJed Brown derr = PetscAbsScalar(dfunc(p, testtime) - aWdot[0]);
1139566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(hW, &aW));
1149566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(hWdot, &aWdot));
1159566063dSJacob Faibussowitsch PetscCall(TSTrajectoryRestoreUpdatedHistoryVecs(tj, &hW, &hWdot));
1163c633725SBarry Smith PetscCheck(!check || serr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state %g > %g", (double)serr, (double)tol);
1173c633725SBarry Smith PetscCheck(!check || derr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state derivative %g > %g", (double)derr, (double)tol);
118c4762a1bSJed Brown }
119c4762a1bSJed Brown
120c4762a1bSJed Brown /* Test on-the-fly reconstruction */
1219566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
1229566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1239566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, W2));
1248337b394SStefano Zampini PetscCall(TSSetMaxTime(ts, 10.));
1259566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, 10));
1269566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts));
1279566063dSJacob Faibussowitsch PetscCall(TSGetTrajectory(ts, &tj));
1289566063dSJacob Faibussowitsch PetscCall(TSTrajectorySetType(tj, ts, TSTRAJECTORYBASIC));
1299566063dSJacob Faibussowitsch PetscCall(TSTrajectorySetFromOptions(tj, ts));
1309566063dSJacob Faibussowitsch PetscCall(TSTrajectorySetSolutionOnly(tj, PETSC_TRUE));
1319566063dSJacob Faibussowitsch PetscCall(TSTrajectorySetUp(tj, ts));
132c4762a1bSJed Brown use1 = PETSC_FALSE;
133c4762a1bSJed Brown use2 = PETSC_TRUE;
1349566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_state", &use1, NULL));
1359566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_der", &use2, NULL));
1369566063dSJacob Faibussowitsch PetscCall(PetscSortReal(10, TT));
137c4762a1bSJed Brown for (i = 0; i < 10; i++) {
1389566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, i));
1399566063dSJacob Faibussowitsch PetscCall(VecSet(W, func(p, TT[i])));
1409566063dSJacob Faibussowitsch PetscCall(TSTrajectorySet(tj, ts, i, TT[i], W));
141c4762a1bSJed Brown if (i) {
142c4762a1bSJed Brown const PetscScalar *aW, *aWdot;
143c4762a1bSJed Brown Vec hW, hWdot;
144c4762a1bSJed Brown PetscReal testtime = TT[i], serr, derr;
145c4762a1bSJed Brown
1469566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGetUpdatedHistoryVecs(tj, ts, testtime, use1 ? &hW : NULL, use2 ? &hWdot : NULL));
147c4762a1bSJed Brown if (use1) {
1489566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(hW, &aW));
1499566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " f(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(func(p, testtime)), (double)PetscRealPart(aW[0])));
150c4762a1bSJed Brown serr = PetscAbsScalar(func(p, testtime) - aW[0]);
1519566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(hW, &aW));
1523c633725SBarry Smith PetscCheck(!check || serr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state %g > %g", (double)serr, (double)tol);
153c4762a1bSJed Brown }
154c4762a1bSJed Brown if (use2) {
1559566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(hWdot, &aWdot));
1569566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "df(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(dfunc(p, testtime)), (double)PetscRealPart(aWdot[0])));
157c4762a1bSJed Brown derr = PetscAbsScalar(dfunc(p, testtime) - aWdot[0]);
1589566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(hWdot, &aWdot));
1593c633725SBarry Smith PetscCheck(!check || i < p || derr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state derivative %g > %g", (double)derr, (double)tol);
160c4762a1bSJed Brown }
1619566063dSJacob Faibussowitsch PetscCall(TSTrajectoryRestoreUpdatedHistoryVecs(tj, use1 ? &hW : NULL, use2 ? &hWdot : NULL));
162c4762a1bSJed Brown }
163c4762a1bSJed Brown }
1649566063dSJacob Faibussowitsch PetscCall(TSRemoveTrajectory(ts));
1659566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
1669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&W));
1679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&W2));
1689566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Wdot));
1699566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
170b122ec5aSJacob Faibussowitsch return 0;
171c4762a1bSJed Brown }
172c4762a1bSJed Brown
173c4762a1bSJed Brown /*TEST
174c4762a1bSJed Brown
175c4762a1bSJed Brown test:
176c4762a1bSJed Brown suffix: 1
177c4762a1bSJed Brown requires: !single
178c4762a1bSJed Brown args: -ts_trajectory_monitor -p 3 -ts_trajectory_reconstruction_order 3 -interptimes 1,9.9,3,1.1,1.1,5.6 -check
179c4762a1bSJed Brown
180c4762a1bSJed Brown test:
181c4762a1bSJed Brown suffix: 2
182c4762a1bSJed Brown requires: !single
183c4762a1bSJed Brown args: -sortkeys -ts_trajectory_monitor -ts_trajectory_type memory -p 3 -ts_trajectory_reconstruction_order 3 -ts_trajectory_adjointmode 0 -interptimes 1,9.9,3,1.1,1.1,5.6 -check
184c4762a1bSJed Brown
185c4762a1bSJed Brown test:
186c4762a1bSJed Brown suffix: 3
187c4762a1bSJed Brown requires: !single
188c4762a1bSJed Brown args: -ts_trajectory_monitor -p 3 -ts_trajectory_reconstruction_order 5 -interptimes 1,9.9,3,1.1,1.1,5.6 -check
189c4762a1bSJed Brown
190c4762a1bSJed Brown TEST*/
191