xref: /petsc/src/ts/tests/ex13.c (revision 67a3cfb084ee9b352765cdad40029d445e453889)
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;
25c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&W);CHKERRQ(ierr);
26c4762a1bSJed Brown   ierr = VecSetSizes(W,1,PETSC_DECIDE);CHKERRQ(ierr);
27c4762a1bSJed Brown   ierr = VecSetUp(W);CHKERRQ(ierr);
28c4762a1bSJed Brown   ierr = VecDuplicate(W,&Wdot);CHKERRQ(ierr);
29c4762a1bSJed Brown   ierr = VecDuplicate(W,&W2);CHKERRQ(ierr);
30c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
31c4762a1bSJed Brown   ierr = TSSetSolution(ts,W2);CHKERRQ(ierr);
32c4762a1bSJed Brown   ierr = TSSetMaxSteps(ts,10);CHKERRQ(ierr);
33c4762a1bSJed Brown   ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr);
34c4762a1bSJed Brown   ierr = TSGetTrajectory(ts,&tj);CHKERRQ(ierr);
35c4762a1bSJed Brown   ierr = TSTrajectorySetType(tj,ts,TSTRAJECTORYBASIC);CHKERRQ(ierr);
36c4762a1bSJed Brown   ierr = TSTrajectorySetFromOptions(tj,ts);CHKERRQ(ierr);
37c4762a1bSJed Brown   ierr = TSTrajectorySetSolutionOnly(tj,PETSC_TRUE);CHKERRQ(ierr);
38c4762a1bSJed Brown   ierr = TSTrajectorySetUp(tj,ts);CHKERRQ(ierr);
39c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-check",&check,NULL);CHKERRQ(ierr);
40c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL);CHKERRQ(ierr);
41c4762a1bSJed Brown   ierr = PetscOptionsGetRealArray(NULL,NULL,"-interptimes",times,&Nt,NULL);CHKERRQ(ierr);
42c4762a1bSJed Brown   sort = PETSC_FALSE;
43c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-sorttimes",&sort,NULL);CHKERRQ(ierr);
44c4762a1bSJed Brown   if (sort) {
45c4762a1bSJed Brown     ierr = PetscSortReal(10,TT);CHKERRQ(ierr);
46c4762a1bSJed Brown   }
47c4762a1bSJed Brown   sort = PETSC_FALSE;
48c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-sortkeys",&sort,NULL);CHKERRQ(ierr);
49c4762a1bSJed Brown   if (sort) {
50c4762a1bSJed Brown     ierr = PetscSortInt(10,II);CHKERRQ(ierr);
51c4762a1bSJed Brown   }
52c4762a1bSJed Brown   p = PetscMax(p,-p);
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   /* populate trajectory */
55c4762a1bSJed Brown   for (i = 0; i < 10; i++) {
56c4762a1bSJed Brown     ierr = VecSet(W,func(p,TT[i]));CHKERRQ(ierr);
57c4762a1bSJed Brown     ierr = TSSetStepNumber(ts,II[i]);CHKERRQ(ierr);
58c4762a1bSJed Brown     ierr = TSTrajectorySet(tj,ts,II[i],TT[i],W);CHKERRQ(ierr);
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 
64c4762a1bSJed Brown     ierr = TSTrajectoryGetVecs(tj,ts,PETSC_DECIDE,&testtime,W,Wdot);CHKERRQ(ierr);
65c4762a1bSJed Brown     ierr = VecGetArrayRead(W,&aW);CHKERRQ(ierr);
66c4762a1bSJed Brown     ierr = VecGetArrayRead(Wdot,&aWdot);CHKERRQ(ierr);
67c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD," f(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(func(p,testtime)),(double)PetscRealPart(aW[0]));CHKERRQ(ierr);
68c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"df(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(dfunc(p,testtime)),(double)PetscRealPart(aWdot[0]));CHKERRQ(ierr);
69c4762a1bSJed Brown     serr = PetscAbsScalar(func(p,testtime)-aW[0]);
70c4762a1bSJed Brown     derr = PetscAbsScalar(dfunc(p,testtime)-aWdot[0]);
71c4762a1bSJed Brown     ierr = VecRestoreArrayRead(W,&aW);CHKERRQ(ierr);
72c4762a1bSJed Brown     ierr = VecRestoreArrayRead(Wdot,&aWdot);CHKERRQ(ierr);
73c4762a1bSJed Brown     if (check && serr > tol) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in state %g > %g",(double)serr,(double)tol);
74c4762a1bSJed Brown     if (check && derr > tol) SETERRQ2(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 
80c4762a1bSJed Brown     ierr = TSTrajectoryGetVecs(tj,ts,PETSC_DECIDE,&testtime,W,NULL);CHKERRQ(ierr);
81c4762a1bSJed Brown     ierr = VecGetArrayRead(W,&aW);CHKERRQ(ierr);
82c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD," f(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(func(p,testtime)),(double)PetscRealPart(aW[0]));CHKERRQ(ierr);
83c4762a1bSJed Brown     serr = PetscAbsScalar(func(p,testtime)-aW[0]);
84c4762a1bSJed Brown     ierr = VecRestoreArrayRead(W,&aW);CHKERRQ(ierr);
85c4762a1bSJed Brown     if (check && serr > tol) SETERRQ2(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 
91c4762a1bSJed Brown     ierr = TSTrajectoryGetVecs(tj,ts,PETSC_DECIDE,&testtime,NULL,Wdot);CHKERRQ(ierr);
92c4762a1bSJed Brown     ierr = VecGetArrayRead(Wdot,&aWdot);CHKERRQ(ierr);
93c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"df(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(dfunc(p,testtime)),(double)PetscRealPart(aWdot[0]));CHKERRQ(ierr);
94c4762a1bSJed Brown     derr = PetscAbsScalar(dfunc(p,testtime)-aWdot[0]);
95c4762a1bSJed Brown     ierr = VecRestoreArrayRead(Wdot,&aWdot);CHKERRQ(ierr);
96c4762a1bSJed Brown     if (check && derr > tol) SETERRQ2(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 
103c4762a1bSJed Brown     ierr = TSTrajectoryGetUpdatedHistoryVecs(tj,ts,testtime,&hW,&hWdot);CHKERRQ(ierr);
104c4762a1bSJed Brown     ierr = VecGetArrayRead(hW,&aW);CHKERRQ(ierr);
105c4762a1bSJed Brown     ierr = VecGetArrayRead(hWdot,&aWdot);CHKERRQ(ierr);
106c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD," f(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(func(p,testtime)),(double)PetscRealPart(aW[0]));CHKERRQ(ierr);
107c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"df(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(dfunc(p,testtime)),(double)PetscRealPart(aWdot[0]));CHKERRQ(ierr);
108c4762a1bSJed Brown     serr = PetscAbsScalar(func(p,testtime)-aW[0]);
109c4762a1bSJed Brown     derr = PetscAbsScalar(dfunc(p,testtime)-aWdot[0]);
110c4762a1bSJed Brown     ierr = VecRestoreArrayRead(hW,&aW);CHKERRQ(ierr);
111c4762a1bSJed Brown     ierr = VecRestoreArrayRead(hWdot,&aWdot);CHKERRQ(ierr);
112c4762a1bSJed Brown     ierr = TSTrajectoryRestoreUpdatedHistoryVecs(tj,&hW,&hWdot);CHKERRQ(ierr);
113c4762a1bSJed Brown     if (check && serr > tol) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in state %g > %g",(double)serr,(double)tol);
114c4762a1bSJed Brown     if (check && derr > tol) SETERRQ2(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 */
118c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
119c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
120c4762a1bSJed Brown   ierr = TSSetSolution(ts,W2);CHKERRQ(ierr);
121c4762a1bSJed Brown   ierr = TSSetMaxSteps(ts,10);CHKERRQ(ierr);
122c4762a1bSJed Brown   ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr);
123c4762a1bSJed Brown   ierr = TSGetTrajectory(ts,&tj);CHKERRQ(ierr);
124c4762a1bSJed Brown   ierr = TSTrajectorySetType(tj,ts,TSTRAJECTORYBASIC);CHKERRQ(ierr);
125c4762a1bSJed Brown   ierr = TSTrajectorySetFromOptions(tj,ts);CHKERRQ(ierr);
126c4762a1bSJed Brown   ierr = TSTrajectorySetSolutionOnly(tj,PETSC_TRUE);CHKERRQ(ierr);
127c4762a1bSJed Brown   ierr = TSTrajectorySetUp(tj,ts);CHKERRQ(ierr);
128c4762a1bSJed Brown   use1 = PETSC_FALSE;
129c4762a1bSJed Brown   use2 = PETSC_TRUE;
130c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-use_state",&use1,NULL);CHKERRQ(ierr);
131c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-use_der",&use2,NULL);CHKERRQ(ierr);
132c4762a1bSJed Brown   ierr = PetscSortReal(10,TT);CHKERRQ(ierr);
133c4762a1bSJed Brown   for (i = 0; i < 10; i++) {
134c4762a1bSJed Brown 
135c4762a1bSJed Brown     ierr = TSSetStepNumber(ts,i);CHKERRQ(ierr);
136c4762a1bSJed Brown     ierr = VecSet(W,func(p,TT[i]));CHKERRQ(ierr);
137c4762a1bSJed Brown     ierr = TSTrajectorySet(tj,ts,i,TT[i],W);CHKERRQ(ierr);
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 
143c4762a1bSJed Brown       ierr = TSTrajectoryGetUpdatedHistoryVecs(tj,ts,testtime,use1 ? &hW : NULL,use2 ? &hWdot : NULL);CHKERRQ(ierr);
144c4762a1bSJed Brown       if (use1) {
145c4762a1bSJed Brown         ierr = VecGetArrayRead(hW,&aW);CHKERRQ(ierr);
146c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_WORLD," f(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(func(p,testtime)),(double)PetscRealPart(aW[0]));CHKERRQ(ierr);
147c4762a1bSJed Brown         serr = PetscAbsScalar(func(p,testtime)-aW[0]);
148c4762a1bSJed Brown         ierr = VecRestoreArrayRead(hW,&aW);CHKERRQ(ierr);
149c4762a1bSJed Brown         if (check && serr > tol) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in state %g > %g",(double)serr,(double)tol);
150c4762a1bSJed Brown       }
151c4762a1bSJed Brown       if (use2) {
152c4762a1bSJed Brown         ierr = VecGetArrayRead(hWdot,&aWdot);CHKERRQ(ierr);
153c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_WORLD,"df(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(dfunc(p,testtime)),(double)PetscRealPart(aWdot[0]));CHKERRQ(ierr);
154c4762a1bSJed Brown         derr = PetscAbsScalar(dfunc(p,testtime)-aWdot[0]);
155c4762a1bSJed Brown         ierr = VecRestoreArrayRead(hWdot,&aWdot);CHKERRQ(ierr);
156c4762a1bSJed Brown         if (check && i >= p && derr > tol) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in state derivative %g > %g",(double)derr,(double)tol);
157c4762a1bSJed Brown       }
158c4762a1bSJed Brown       ierr = TSTrajectoryRestoreUpdatedHistoryVecs(tj,use1 ? &hW : NULL,use2 ? &hWdot : NULL);CHKERRQ(ierr);
159c4762a1bSJed Brown     }
160c4762a1bSJed Brown   }
161*67a3cfb0SHong Zhang   ierr = TSRemoveTrajectory(ts);CHKERRQ(ierr);
162c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
163c4762a1bSJed Brown   ierr = VecDestroy(&W);CHKERRQ(ierr);
164c4762a1bSJed Brown   ierr = VecDestroy(&W2);CHKERRQ(ierr);
165c4762a1bSJed Brown   ierr = VecDestroy(&Wdot);CHKERRQ(ierr);
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