xref: /petsc/src/ts/tests/ex13.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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 
9c4762a1bSJed Brown PetscScalar func(PetscInt p, PetscReal t)  { return p ? t*func(p-1,t) : 1.0; }
10c4762a1bSJed Brown PetscScalar dfunc(PetscInt p, PetscReal t)  { return p > 0 ? (PetscReal)p*func(p-1,t) : 0.0; }
11c4762a1bSJed Brown 
12c4762a1bSJed Brown int main(int argc,char **argv)
13c4762a1bSJed Brown {
14c4762a1bSJed Brown   TS             ts;
15c4762a1bSJed Brown   Vec            W,W2,Wdot;
16c4762a1bSJed Brown   TSTrajectory   tj;
17c4762a1bSJed Brown   PetscReal      times[10], tol = PETSC_SMALL;
18c4762a1bSJed Brown   PetscReal      TT[10] = { 2, 9, 1, 3, 6, 7, 5, 10, 4, 8 };
19c4762a1bSJed Brown   PetscInt       i, p = 1, Nt = 10;
20c4762a1bSJed Brown   PetscInt       II[10] = { 1, 4, 9, 2, 3, 6, 5, 8, 0, 7 };
21c4762a1bSJed Brown   PetscBool      sort,use1,use2,check = PETSC_FALSE;
22c4762a1bSJed Brown   PetscErrorCode ierr;
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&W));
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(W,1,PETSC_DECIDE));
27*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetUp(W));
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(W,&Wdot));
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(W,&W2));
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
31*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,W2));
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxSteps(ts,10));
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSaveTrajectory(ts));
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTrajectory(ts,&tj));
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSTrajectorySetType(tj,ts,TSTRAJECTORYBASIC));
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSTrajectorySetFromOptions(tj,ts));
37*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSTrajectorySetSolutionOnly(tj,PETSC_TRUE));
38*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSTrajectorySetUp(tj,ts));
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-check",&check,NULL));
40*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL));
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetRealArray(NULL,NULL,"-interptimes",times,&Nt,NULL));
42c4762a1bSJed Brown   sort = PETSC_FALSE;
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-sorttimes",&sort,NULL));
44c4762a1bSJed Brown   if (sort) {
45*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSortReal(10,TT));
46c4762a1bSJed Brown   }
47c4762a1bSJed Brown   sort = PETSC_FALSE;
48*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-sortkeys",&sort,NULL));
49c4762a1bSJed Brown   if (sort) {
50*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSortInt(10,II));
51c4762a1bSJed Brown   }
52c4762a1bSJed Brown   p = PetscMax(p,-p);
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   /* populate trajectory */
55c4762a1bSJed Brown   for (i = 0; i < 10; i++) {
56*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(W,func(p,TT[i])));
57*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetStepNumber(ts,II[i]));
58*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSTrajectorySet(tj,ts,II[i],TT[i],W));
59c4762a1bSJed Brown   }
60c4762a1bSJed Brown   for (i = 0; i < Nt; i++) {
61c4762a1bSJed Brown     PetscReal testtime = times[i], serr, derr;
62c4762a1bSJed Brown     const PetscScalar *aW,*aWdot;
63c4762a1bSJed Brown 
64*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSTrajectoryGetVecs(tj,ts,PETSC_DECIDE,&testtime,W,Wdot));
65*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArrayRead(W,&aW));
66*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArrayRead(Wdot,&aWdot));
67*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," f(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(func(p,testtime)),(double)PetscRealPart(aW[0])));
68*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"df(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(dfunc(p,testtime)),(double)PetscRealPart(aWdot[0])));
69c4762a1bSJed Brown     serr = PetscAbsScalar(func(p,testtime)-aW[0]);
70c4762a1bSJed Brown     derr = PetscAbsScalar(dfunc(p,testtime)-aWdot[0]);
71*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArrayRead(W,&aW));
72*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArrayRead(Wdot,&aWdot));
733c633725SBarry Smith     PetscCheck(!check || serr <= tol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in state %g > %g",(double)serr,(double)tol);
743c633725SBarry Smith     PetscCheck(!check || derr <= tol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in state derivative %g > %g",(double)derr,(double)tol);
75c4762a1bSJed Brown   }
76c4762a1bSJed Brown   for (i = Nt-1; i >= 0; i--) {
77c4762a1bSJed Brown     PetscReal         testtime = times[i], serr;
78c4762a1bSJed Brown     const PetscScalar *aW;
79c4762a1bSJed Brown 
80*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSTrajectoryGetVecs(tj,ts,PETSC_DECIDE,&testtime,W,NULL));
81*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArrayRead(W,&aW));
82*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," f(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(func(p,testtime)),(double)PetscRealPart(aW[0])));
83c4762a1bSJed Brown     serr = PetscAbsScalar(func(p,testtime)-aW[0]);
84*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArrayRead(W,&aW));
853c633725SBarry Smith     PetscCheck(!check || serr <= tol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in state %g > %g",(double)serr,(double)tol);
86c4762a1bSJed Brown   }
87c4762a1bSJed Brown   for (i = Nt-1; i >= 0; i--) {
88c4762a1bSJed Brown     PetscReal         testtime = times[i], derr;
89c4762a1bSJed Brown     const PetscScalar *aWdot;
90c4762a1bSJed Brown 
91*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSTrajectoryGetVecs(tj,ts,PETSC_DECIDE,&testtime,NULL,Wdot));
92*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArrayRead(Wdot,&aWdot));
93*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"df(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(dfunc(p,testtime)),(double)PetscRealPart(aWdot[0])));
94c4762a1bSJed Brown     derr = PetscAbsScalar(dfunc(p,testtime)-aWdot[0]);
95*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArrayRead(Wdot,&aWdot));
963c633725SBarry Smith     PetscCheck(!check || derr <= tol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in state derivative %g > %g",(double)derr,(double)tol);
97c4762a1bSJed Brown   }
98c4762a1bSJed Brown   for (i = 0; i < Nt; i++) {
99c4762a1bSJed Brown     PetscReal         testtime = times[i], serr, derr;
100c4762a1bSJed Brown     const PetscScalar *aW,*aWdot;
101c4762a1bSJed Brown     Vec               hW,hWdot;
102c4762a1bSJed Brown 
103*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSTrajectoryGetUpdatedHistoryVecs(tj,ts,testtime,&hW,&hWdot));
104*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArrayRead(hW,&aW));
105*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArrayRead(hWdot,&aWdot));
106*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," f(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(func(p,testtime)),(double)PetscRealPart(aW[0])));
107*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"df(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(dfunc(p,testtime)),(double)PetscRealPart(aWdot[0])));
108c4762a1bSJed Brown     serr = PetscAbsScalar(func(p,testtime)-aW[0]);
109c4762a1bSJed Brown     derr = PetscAbsScalar(dfunc(p,testtime)-aWdot[0]);
110*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArrayRead(hW,&aW));
111*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArrayRead(hWdot,&aWdot));
112*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSTrajectoryRestoreUpdatedHistoryVecs(tj,&hW,&hWdot));
1133c633725SBarry Smith     PetscCheck(!check || serr <= tol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in state %g > %g",(double)serr,(double)tol);
1143c633725SBarry Smith     PetscCheck(!check || derr <= tol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in state derivative %g > %g",(double)derr,(double)tol);
115c4762a1bSJed Brown   }
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   /* Test on-the-fly reconstruction */
118*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
119*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
120*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,W2));
121*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxSteps(ts,10));
122*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSaveTrajectory(ts));
123*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTrajectory(ts,&tj));
124*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSTrajectorySetType(tj,ts,TSTRAJECTORYBASIC));
125*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSTrajectorySetFromOptions(tj,ts));
126*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSTrajectorySetSolutionOnly(tj,PETSC_TRUE));
127*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSTrajectorySetUp(tj,ts));
128c4762a1bSJed Brown   use1 = PETSC_FALSE;
129c4762a1bSJed Brown   use2 = PETSC_TRUE;
130*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-use_state",&use1,NULL));
131*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-use_der",&use2,NULL));
132*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSortReal(10,TT));
133c4762a1bSJed Brown   for (i = 0; i < 10; i++) {
134c4762a1bSJed Brown 
135*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetStepNumber(ts,i));
136*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(W,func(p,TT[i])));
137*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSTrajectorySet(tj,ts,i,TT[i],W));
138c4762a1bSJed Brown     if (i) {
139c4762a1bSJed Brown       const PetscScalar *aW,*aWdot;
140c4762a1bSJed Brown       Vec               hW,hWdot;
141c4762a1bSJed Brown       PetscReal         testtime = TT[i], serr, derr;
142c4762a1bSJed Brown 
143*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TSTrajectoryGetUpdatedHistoryVecs(tj,ts,testtime,use1 ? &hW : NULL,use2 ? &hWdot : NULL));
144c4762a1bSJed Brown       if (use1) {
145*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecGetArrayRead(hW,&aW));
146*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," f(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(func(p,testtime)),(double)PetscRealPart(aW[0])));
147c4762a1bSJed Brown         serr = PetscAbsScalar(func(p,testtime)-aW[0]);
148*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecRestoreArrayRead(hW,&aW));
1493c633725SBarry Smith         PetscCheck(!check || serr <= tol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in state %g > %g",(double)serr,(double)tol);
150c4762a1bSJed Brown       }
151c4762a1bSJed Brown       if (use2) {
152*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecGetArrayRead(hWdot,&aWdot));
153*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"df(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(dfunc(p,testtime)),(double)PetscRealPart(aWdot[0])));
154c4762a1bSJed Brown         derr = PetscAbsScalar(dfunc(p,testtime)-aWdot[0]);
155*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecRestoreArrayRead(hWdot,&aWdot));
1563c633725SBarry Smith         PetscCheck(!check || i < p || derr <= tol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in state derivative %g > %g",(double)derr,(double)tol);
157c4762a1bSJed Brown       }
158*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TSTrajectoryRestoreUpdatedHistoryVecs(tj,use1 ? &hW : NULL,use2 ? &hWdot : NULL));
159c4762a1bSJed Brown     }
160c4762a1bSJed Brown   }
161*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRemoveTrajectory(ts));
162*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
163*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&W));
164*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&W2));
165*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Wdot));
166c4762a1bSJed Brown   ierr = PetscFinalize();
167c4762a1bSJed Brown   return ierr;
168c4762a1bSJed Brown }
169c4762a1bSJed Brown 
170c4762a1bSJed Brown /*TEST
171c4762a1bSJed Brown 
172c4762a1bSJed Brown test:
173c4762a1bSJed Brown   suffix: 1
174c4762a1bSJed Brown   requires: !single
175c4762a1bSJed Brown   args: -ts_trajectory_monitor -p 3 -ts_trajectory_reconstruction_order 3 -interptimes 1,9.9,3,1.1,1.1,5.6 -check
176c4762a1bSJed Brown 
177c4762a1bSJed Brown test:
178c4762a1bSJed Brown   suffix: 2
179c4762a1bSJed Brown   requires: !single
180c4762a1bSJed 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
181c4762a1bSJed Brown 
182c4762a1bSJed Brown test:
183c4762a1bSJed Brown   suffix: 3
184c4762a1bSJed Brown   requires: !single
185c4762a1bSJed Brown   args: -ts_trajectory_monitor -p 3 -ts_trajectory_reconstruction_order 5 -interptimes 1,9.9,3,1.1,1.1,5.6 -check
186c4762a1bSJed Brown 
187c4762a1bSJed Brown TEST*/
188