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