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