xref: /petsc/src/ts/tests/ex13.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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*9371c9d4SSatish Balay PetscScalar func(PetscInt p, PetscReal t) {
10*9371c9d4SSatish Balay   return p ? t * func(p - 1, t) : 1.0;
11*9371c9d4SSatish Balay }
12*9371c9d4SSatish Balay PetscScalar dfunc(PetscInt p, PetscReal t) {
13*9371c9d4SSatish Balay   return p > 0 ? (PetscReal)p * func(p - 1, t) : 0.0;
14*9371c9d4SSatish Balay }
15c4762a1bSJed Brown 
16*9371c9d4SSatish Balay int main(int argc, char **argv) {
17c4762a1bSJed Brown   TS           ts;
18c4762a1bSJed Brown   Vec          W, W2, Wdot;
19c4762a1bSJed Brown   TSTrajectory tj;
20c4762a1bSJed Brown   PetscReal    times[10], tol = PETSC_SMALL;
21c4762a1bSJed Brown   PetscReal    TT[10] = {2, 9, 1, 3, 6, 7, 5, 10, 4, 8};
22c4762a1bSJed Brown   PetscInt     i, p = 1, Nt = 10;
23c4762a1bSJed Brown   PetscInt     II[10] = {1, 4, 9, 2, 3, 6, 5, 8, 0, 7};
24c4762a1bSJed Brown   PetscBool    sort, use1, use2, check = PETSC_FALSE;
25c4762a1bSJed Brown 
26327415f7SBarry Smith   PetscFunctionBeginUser;
279566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
289566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &W));
299566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(W, 1, PETSC_DECIDE));
309566063dSJacob Faibussowitsch   PetscCall(VecSetUp(W));
319566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(W, &Wdot));
329566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(W, &W2));
339566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
349566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, W2));
359566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, 10));
369566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(ts));
379566063dSJacob Faibussowitsch   PetscCall(TSGetTrajectory(ts, &tj));
389566063dSJacob Faibussowitsch   PetscCall(TSTrajectorySetType(tj, ts, TSTRAJECTORYBASIC));
399566063dSJacob Faibussowitsch   PetscCall(TSTrajectorySetFromOptions(tj, ts));
409566063dSJacob Faibussowitsch   PetscCall(TSTrajectorySetSolutionOnly(tj, PETSC_TRUE));
419566063dSJacob Faibussowitsch   PetscCall(TSTrajectorySetUp(tj, ts));
429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-check", &check, NULL));
439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-p", &p, NULL));
449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetRealArray(NULL, NULL, "-interptimes", times, &Nt, NULL));
45c4762a1bSJed Brown   sort = PETSC_FALSE;
469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-sorttimes", &sort, NULL));
471baa6e33SBarry Smith   if (sort) PetscCall(PetscSortReal(10, TT));
48c4762a1bSJed Brown   sort = PETSC_FALSE;
499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-sortkeys", &sort, NULL));
501baa6e33SBarry Smith   if (sort) PetscCall(PetscSortInt(10, II));
51c4762a1bSJed Brown   p = PetscMax(p, -p);
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   /* populate trajectory */
54c4762a1bSJed Brown   for (i = 0; i < 10; i++) {
559566063dSJacob Faibussowitsch     PetscCall(VecSet(W, func(p, TT[i])));
569566063dSJacob Faibussowitsch     PetscCall(TSSetStepNumber(ts, II[i]));
579566063dSJacob Faibussowitsch     PetscCall(TSTrajectorySet(tj, ts, II[i], TT[i], W));
58c4762a1bSJed Brown   }
59c4762a1bSJed Brown   for (i = 0; i < Nt; i++) {
60c4762a1bSJed Brown     PetscReal          testtime = times[i], serr, derr;
61c4762a1bSJed Brown     const PetscScalar *aW, *aWdot;
62c4762a1bSJed Brown 
639566063dSJacob Faibussowitsch     PetscCall(TSTrajectoryGetVecs(tj, ts, PETSC_DECIDE, &testtime, W, Wdot));
649566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(W, &aW));
659566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(Wdot, &aWdot));
669566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, " f(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(func(p, testtime)), (double)PetscRealPart(aW[0])));
679566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "df(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(dfunc(p, testtime)), (double)PetscRealPart(aWdot[0])));
68c4762a1bSJed Brown     serr = PetscAbsScalar(func(p, testtime) - aW[0]);
69c4762a1bSJed Brown     derr = PetscAbsScalar(dfunc(p, testtime) - aWdot[0]);
709566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(W, &aW));
719566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(Wdot, &aWdot));
723c633725SBarry Smith     PetscCheck(!check || serr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state %g > %g", (double)serr, (double)tol);
733c633725SBarry Smith     PetscCheck(!check || derr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state derivative %g > %g", (double)derr, (double)tol);
74c4762a1bSJed Brown   }
75c4762a1bSJed Brown   for (i = Nt - 1; i >= 0; i--) {
76c4762a1bSJed Brown     PetscReal          testtime = times[i], serr;
77c4762a1bSJed Brown     const PetscScalar *aW;
78c4762a1bSJed Brown 
799566063dSJacob Faibussowitsch     PetscCall(TSTrajectoryGetVecs(tj, ts, PETSC_DECIDE, &testtime, W, NULL));
809566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(W, &aW));
819566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, " f(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(func(p, testtime)), (double)PetscRealPart(aW[0])));
82c4762a1bSJed Brown     serr = PetscAbsScalar(func(p, testtime) - aW[0]);
839566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(W, &aW));
843c633725SBarry Smith     PetscCheck(!check || serr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state %g > %g", (double)serr, (double)tol);
85c4762a1bSJed Brown   }
86c4762a1bSJed Brown   for (i = Nt - 1; i >= 0; i--) {
87c4762a1bSJed Brown     PetscReal          testtime = times[i], derr;
88c4762a1bSJed Brown     const PetscScalar *aWdot;
89c4762a1bSJed Brown 
909566063dSJacob Faibussowitsch     PetscCall(TSTrajectoryGetVecs(tj, ts, PETSC_DECIDE, &testtime, NULL, Wdot));
919566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(Wdot, &aWdot));
929566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "df(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(dfunc(p, testtime)), (double)PetscRealPart(aWdot[0])));
93c4762a1bSJed Brown     derr = PetscAbsScalar(dfunc(p, testtime) - aWdot[0]);
949566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(Wdot, &aWdot));
953c633725SBarry Smith     PetscCheck(!check || derr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state derivative %g > %g", (double)derr, (double)tol);
96c4762a1bSJed Brown   }
97c4762a1bSJed Brown   for (i = 0; i < Nt; i++) {
98c4762a1bSJed Brown     PetscReal          testtime = times[i], serr, derr;
99c4762a1bSJed Brown     const PetscScalar *aW, *aWdot;
100c4762a1bSJed Brown     Vec                hW, hWdot;
101c4762a1bSJed Brown 
1029566063dSJacob Faibussowitsch     PetscCall(TSTrajectoryGetUpdatedHistoryVecs(tj, ts, testtime, &hW, &hWdot));
1039566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(hW, &aW));
1049566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(hWdot, &aWdot));
1059566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, " f(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(func(p, testtime)), (double)PetscRealPart(aW[0])));
1069566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "df(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(dfunc(p, testtime)), (double)PetscRealPart(aWdot[0])));
107c4762a1bSJed Brown     serr = PetscAbsScalar(func(p, testtime) - aW[0]);
108c4762a1bSJed Brown     derr = PetscAbsScalar(dfunc(p, testtime) - aWdot[0]);
1099566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(hW, &aW));
1109566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(hWdot, &aWdot));
1119566063dSJacob Faibussowitsch     PetscCall(TSTrajectoryRestoreUpdatedHistoryVecs(tj, &hW, &hWdot));
1123c633725SBarry Smith     PetscCheck(!check || serr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state %g > %g", (double)serr, (double)tol);
1133c633725SBarry Smith     PetscCheck(!check || derr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state derivative %g > %g", (double)derr, (double)tol);
114c4762a1bSJed Brown   }
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   /* Test on-the-fly reconstruction */
1179566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1189566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1199566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, W2));
1209566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, 10));
1219566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(ts));
1229566063dSJacob Faibussowitsch   PetscCall(TSGetTrajectory(ts, &tj));
1239566063dSJacob Faibussowitsch   PetscCall(TSTrajectorySetType(tj, ts, TSTRAJECTORYBASIC));
1249566063dSJacob Faibussowitsch   PetscCall(TSTrajectorySetFromOptions(tj, ts));
1259566063dSJacob Faibussowitsch   PetscCall(TSTrajectorySetSolutionOnly(tj, PETSC_TRUE));
1269566063dSJacob Faibussowitsch   PetscCall(TSTrajectorySetUp(tj, ts));
127c4762a1bSJed Brown   use1 = PETSC_FALSE;
128c4762a1bSJed Brown   use2 = PETSC_TRUE;
1299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_state", &use1, NULL));
1309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_der", &use2, NULL));
1319566063dSJacob Faibussowitsch   PetscCall(PetscSortReal(10, TT));
132c4762a1bSJed Brown   for (i = 0; i < 10; i++) {
1339566063dSJacob Faibussowitsch     PetscCall(TSSetStepNumber(ts, i));
1349566063dSJacob Faibussowitsch     PetscCall(VecSet(W, func(p, TT[i])));
1359566063dSJacob Faibussowitsch     PetscCall(TSTrajectorySet(tj, ts, i, TT[i], W));
136c4762a1bSJed Brown     if (i) {
137c4762a1bSJed Brown       const PetscScalar *aW, *aWdot;
138c4762a1bSJed Brown       Vec                hW, hWdot;
139c4762a1bSJed Brown       PetscReal          testtime = TT[i], serr, derr;
140c4762a1bSJed Brown 
1419566063dSJacob Faibussowitsch       PetscCall(TSTrajectoryGetUpdatedHistoryVecs(tj, ts, testtime, use1 ? &hW : NULL, use2 ? &hWdot : NULL));
142c4762a1bSJed Brown       if (use1) {
1439566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(hW, &aW));
1449566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, " f(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(func(p, testtime)), (double)PetscRealPart(aW[0])));
145c4762a1bSJed Brown         serr = PetscAbsScalar(func(p, testtime) - aW[0]);
1469566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(hW, &aW));
1473c633725SBarry Smith         PetscCheck(!check || serr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state %g > %g", (double)serr, (double)tol);
148c4762a1bSJed Brown       }
149c4762a1bSJed Brown       if (use2) {
1509566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(hWdot, &aWdot));
1519566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "df(%g) = %g (reconstructed %g)\n", (double)testtime, (double)PetscRealPart(dfunc(p, testtime)), (double)PetscRealPart(aWdot[0])));
152c4762a1bSJed Brown         derr = PetscAbsScalar(dfunc(p, testtime) - aWdot[0]);
1539566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(hWdot, &aWdot));
1543c633725SBarry Smith         PetscCheck(!check || i < p || derr <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in state derivative %g > %g", (double)derr, (double)tol);
155c4762a1bSJed Brown       }
1569566063dSJacob Faibussowitsch       PetscCall(TSTrajectoryRestoreUpdatedHistoryVecs(tj, use1 ? &hW : NULL, use2 ? &hWdot : NULL));
157c4762a1bSJed Brown     }
158c4762a1bSJed Brown   }
1599566063dSJacob Faibussowitsch   PetscCall(TSRemoveTrajectory(ts));
1609566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1619566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&W));
1629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&W2));
1639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Wdot));
1649566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
165b122ec5aSJacob Faibussowitsch   return 0;
166c4762a1bSJed Brown }
167c4762a1bSJed Brown 
168c4762a1bSJed Brown /*TEST
169c4762a1bSJed Brown 
170c4762a1bSJed Brown test:
171c4762a1bSJed Brown   suffix: 1
172c4762a1bSJed Brown   requires: !single
173c4762a1bSJed Brown   args: -ts_trajectory_monitor -p 3 -ts_trajectory_reconstruction_order 3 -interptimes 1,9.9,3,1.1,1.1,5.6 -check
174c4762a1bSJed Brown 
175c4762a1bSJed Brown test:
176c4762a1bSJed Brown   suffix: 2
177c4762a1bSJed Brown   requires: !single
178c4762a1bSJed 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
179c4762a1bSJed Brown 
180c4762a1bSJed Brown test:
181c4762a1bSJed Brown   suffix: 3
182c4762a1bSJed Brown   requires: !single
183c4762a1bSJed Brown   args: -ts_trajectory_monitor -p 3 -ts_trajectory_reconstruction_order 5 -interptimes 1,9.9,3,1.1,1.1,5.6 -check
184c4762a1bSJed Brown 
185c4762a1bSJed Brown TEST*/
186