xref: /petsc/src/ts/tests/ex13.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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 
9*d71ae5a4SJacob Faibussowitsch PetscScalar func(PetscInt p, PetscReal t)
10*d71ae5a4SJacob Faibussowitsch {
119371c9d4SSatish Balay   return p ? t * func(p - 1, t) : 1.0;
129371c9d4SSatish Balay }
13*d71ae5a4SJacob Faibussowitsch PetscScalar dfunc(PetscInt p, PetscReal t)
14*d71ae5a4SJacob Faibussowitsch {
159371c9d4SSatish Balay   return p > 0 ? (PetscReal)p * func(p - 1, t) : 0.0;
169371c9d4SSatish Balay }
17c4762a1bSJed Brown 
18*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
19*d71ae5a4SJacob 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;
309566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, 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));
389566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, 10));
399566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(ts));
409566063dSJacob Faibussowitsch   PetscCall(TSGetTrajectory(ts, &tj));
419566063dSJacob Faibussowitsch   PetscCall(TSTrajectorySetType(tj, ts, TSTRAJECTORYBASIC));
429566063dSJacob Faibussowitsch   PetscCall(TSTrajectorySetFromOptions(tj, ts));
439566063dSJacob Faibussowitsch   PetscCall(TSTrajectorySetSolutionOnly(tj, PETSC_TRUE));
449566063dSJacob Faibussowitsch   PetscCall(TSTrajectorySetUp(tj, ts));
459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-check", &check, NULL));
469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-p", &p, NULL));
479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetRealArray(NULL, NULL, "-interptimes", times, &Nt, NULL));
48c4762a1bSJed Brown   sort = PETSC_FALSE;
499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-sorttimes", &sort, NULL));
501baa6e33SBarry Smith   if (sort) PetscCall(PetscSortReal(10, TT));
51c4762a1bSJed Brown   sort = PETSC_FALSE;
529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-sortkeys", &sort, NULL));
531baa6e33SBarry Smith   if (sort) PetscCall(PetscSortInt(10, II));
54c4762a1bSJed Brown   p = PetscMax(p, -p);
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   /* populate trajectory */
57c4762a1bSJed Brown   for (i = 0; i < 10; i++) {
589566063dSJacob Faibussowitsch     PetscCall(VecSet(W, func(p, TT[i])));
599566063dSJacob Faibussowitsch     PetscCall(TSSetStepNumber(ts, II[i]));
609566063dSJacob Faibussowitsch     PetscCall(TSTrajectorySet(tj, ts, II[i], TT[i], W));
61c4762a1bSJed Brown   }
62c4762a1bSJed Brown   for (i = 0; i < Nt; i++) {
63c4762a1bSJed Brown     PetscReal          testtime = times[i], serr, derr;
64c4762a1bSJed Brown     const PetscScalar *aW, *aWdot;
65c4762a1bSJed Brown 
669566063dSJacob Faibussowitsch     PetscCall(TSTrajectoryGetVecs(tj, ts, PETSC_DECIDE, &testtime, W, Wdot));
679566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(W, &aW));
689566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(Wdot, &aWdot));
699566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, " f(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(func(p, testtime)), (double)PetscRealPart(aW[0])));
709566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "df(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(dfunc(p, testtime)), (double)PetscRealPart(aWdot[0])));
71c4762a1bSJed Brown     serr = PetscAbsScalar(func(p, testtime) - aW[0]);
72c4762a1bSJed Brown     derr = PetscAbsScalar(dfunc(p, testtime) - aWdot[0]);
739566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(W, &aW));
749566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(Wdot, &aWdot));
753c633725SBarry Smith     PetscCheck(!check || serr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state %g > %g", (double)serr, (double)tol);
763c633725SBarry Smith     PetscCheck(!check || derr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state derivative %g > %g", (double)derr, (double)tol);
77c4762a1bSJed Brown   }
78c4762a1bSJed Brown   for (i = Nt - 1; i >= 0; i--) {
79c4762a1bSJed Brown     PetscReal          testtime = times[i], serr;
80c4762a1bSJed Brown     const PetscScalar *aW;
81c4762a1bSJed Brown 
829566063dSJacob Faibussowitsch     PetscCall(TSTrajectoryGetVecs(tj, ts, PETSC_DECIDE, &testtime, W, NULL));
839566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(W, &aW));
849566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, " f(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(func(p, testtime)), (double)PetscRealPart(aW[0])));
85c4762a1bSJed Brown     serr = PetscAbsScalar(func(p, testtime) - aW[0]);
869566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(W, &aW));
873c633725SBarry Smith     PetscCheck(!check || serr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state %g > %g", (double)serr, (double)tol);
88c4762a1bSJed Brown   }
89c4762a1bSJed Brown   for (i = Nt - 1; i >= 0; i--) {
90c4762a1bSJed Brown     PetscReal          testtime = times[i], derr;
91c4762a1bSJed Brown     const PetscScalar *aWdot;
92c4762a1bSJed Brown 
939566063dSJacob Faibussowitsch     PetscCall(TSTrajectoryGetVecs(tj, ts, PETSC_DECIDE, &testtime, NULL, Wdot));
949566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(Wdot, &aWdot));
959566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "df(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(dfunc(p, testtime)), (double)PetscRealPart(aWdot[0])));
96c4762a1bSJed Brown     derr = PetscAbsScalar(dfunc(p, testtime) - aWdot[0]);
979566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(Wdot, &aWdot));
983c633725SBarry Smith     PetscCheck(!check || derr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state derivative %g > %g", (double)derr, (double)tol);
99c4762a1bSJed Brown   }
100c4762a1bSJed Brown   for (i = 0; i < Nt; i++) {
101c4762a1bSJed Brown     PetscReal          testtime = times[i], serr, derr;
102c4762a1bSJed Brown     const PetscScalar *aW, *aWdot;
103c4762a1bSJed Brown     Vec                hW, hWdot;
104c4762a1bSJed Brown 
1059566063dSJacob Faibussowitsch     PetscCall(TSTrajectoryGetUpdatedHistoryVecs(tj, ts, testtime, &hW, &hWdot));
1069566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(hW, &aW));
1079566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(hWdot, &aWdot));
1089566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, " f(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(func(p, testtime)), (double)PetscRealPart(aW[0])));
1099566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "df(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(dfunc(p, testtime)), (double)PetscRealPart(aWdot[0])));
110c4762a1bSJed Brown     serr = PetscAbsScalar(func(p, testtime) - aW[0]);
111c4762a1bSJed Brown     derr = PetscAbsScalar(dfunc(p, testtime) - aWdot[0]);
1129566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(hW, &aW));
1139566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(hWdot, &aWdot));
1149566063dSJacob Faibussowitsch     PetscCall(TSTrajectoryRestoreUpdatedHistoryVecs(tj, &hW, &hWdot));
1153c633725SBarry Smith     PetscCheck(!check || serr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state %g > %g", (double)serr, (double)tol);
1163c633725SBarry Smith     PetscCheck(!check || derr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state derivative %g > %g", (double)derr, (double)tol);
117c4762a1bSJed Brown   }
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   /* Test on-the-fly reconstruction */
1209566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1219566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1229566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, W2));
1239566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, 10));
1249566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(ts));
1259566063dSJacob Faibussowitsch   PetscCall(TSGetTrajectory(ts, &tj));
1269566063dSJacob Faibussowitsch   PetscCall(TSTrajectorySetType(tj, ts, TSTRAJECTORYBASIC));
1279566063dSJacob Faibussowitsch   PetscCall(TSTrajectorySetFromOptions(tj, ts));
1289566063dSJacob Faibussowitsch   PetscCall(TSTrajectorySetSolutionOnly(tj, PETSC_TRUE));
1299566063dSJacob Faibussowitsch   PetscCall(TSTrajectorySetUp(tj, ts));
130c4762a1bSJed Brown   use1 = PETSC_FALSE;
131c4762a1bSJed Brown   use2 = PETSC_TRUE;
1329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_state", &use1, NULL));
1339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_der", &use2, NULL));
1349566063dSJacob Faibussowitsch   PetscCall(PetscSortReal(10, TT));
135c4762a1bSJed Brown   for (i = 0; i < 10; i++) {
1369566063dSJacob Faibussowitsch     PetscCall(TSSetStepNumber(ts, i));
1379566063dSJacob Faibussowitsch     PetscCall(VecSet(W, func(p, TT[i])));
1389566063dSJacob Faibussowitsch     PetscCall(TSTrajectorySet(tj, ts, i, TT[i], W));
139c4762a1bSJed Brown     if (i) {
140c4762a1bSJed Brown       const PetscScalar *aW, *aWdot;
141c4762a1bSJed Brown       Vec                hW, hWdot;
142c4762a1bSJed Brown       PetscReal          testtime = TT[i], serr, derr;
143c4762a1bSJed Brown 
1449566063dSJacob Faibussowitsch       PetscCall(TSTrajectoryGetUpdatedHistoryVecs(tj, ts, testtime, use1 ? &hW : NULL, use2 ? &hWdot : NULL));
145c4762a1bSJed Brown       if (use1) {
1469566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(hW, &aW));
1479566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, " f(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(func(p, testtime)), (double)PetscRealPart(aW[0])));
148c4762a1bSJed Brown         serr = PetscAbsScalar(func(p, testtime) - aW[0]);
1499566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(hW, &aW));
1503c633725SBarry Smith         PetscCheck(!check || serr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state %g > %g", (double)serr, (double)tol);
151c4762a1bSJed Brown       }
152c4762a1bSJed Brown       if (use2) {
1539566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(hWdot, &aWdot));
1549566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "df(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(dfunc(p, testtime)), (double)PetscRealPart(aWdot[0])));
155c4762a1bSJed Brown         derr = PetscAbsScalar(dfunc(p, testtime) - aWdot[0]);
1569566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(hWdot, &aWdot));
1573c633725SBarry Smith         PetscCheck(!check || i < p || derr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state derivative %g > %g", (double)derr, (double)tol);
158c4762a1bSJed Brown       }
1599566063dSJacob Faibussowitsch       PetscCall(TSTrajectoryRestoreUpdatedHistoryVecs(tj, use1 ? &hW : NULL, use2 ? &hWdot : NULL));
160c4762a1bSJed Brown     }
161c4762a1bSJed Brown   }
1629566063dSJacob Faibussowitsch   PetscCall(TSRemoveTrajectory(ts));
1639566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&W));
1659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&W2));
1669566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Wdot));
1679566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
168b122ec5aSJacob Faibussowitsch   return 0;
169c4762a1bSJed Brown }
170c4762a1bSJed Brown 
171c4762a1bSJed Brown /*TEST
172c4762a1bSJed Brown 
173c4762a1bSJed Brown test:
174c4762a1bSJed Brown   suffix: 1
175c4762a1bSJed Brown   requires: !single
176c4762a1bSJed Brown   args: -ts_trajectory_monitor -p 3 -ts_trajectory_reconstruction_order 3 -interptimes 1,9.9,3,1.1,1.1,5.6 -check
177c4762a1bSJed Brown 
178c4762a1bSJed Brown test:
179c4762a1bSJed Brown   suffix: 2
180c4762a1bSJed Brown   requires: !single
181c4762a1bSJed 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
182c4762a1bSJed Brown 
183c4762a1bSJed Brown test:
184c4762a1bSJed Brown   suffix: 3
185c4762a1bSJed Brown   requires: !single
186c4762a1bSJed Brown   args: -ts_trajectory_monitor -p 3 -ts_trajectory_reconstruction_order 5 -interptimes 1,9.9,3,1.1,1.1,5.6 -check
187c4762a1bSJed Brown 
188c4762a1bSJed Brown TEST*/
189