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 23*327415f7SBarry Smith PetscFunctionBeginUser; 249566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 259566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&W)); 269566063dSJacob Faibussowitsch PetscCall(VecSetSizes(W,1,PETSC_DECIDE)); 279566063dSJacob Faibussowitsch PetscCall(VecSetUp(W)); 289566063dSJacob Faibussowitsch PetscCall(VecDuplicate(W,&Wdot)); 299566063dSJacob Faibussowitsch PetscCall(VecDuplicate(W,&W2)); 309566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 319566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts,W2)); 329566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts,10)); 339566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts)); 349566063dSJacob Faibussowitsch PetscCall(TSGetTrajectory(ts,&tj)); 359566063dSJacob Faibussowitsch PetscCall(TSTrajectorySetType(tj,ts,TSTRAJECTORYBASIC)); 369566063dSJacob Faibussowitsch PetscCall(TSTrajectorySetFromOptions(tj,ts)); 379566063dSJacob Faibussowitsch PetscCall(TSTrajectorySetSolutionOnly(tj,PETSC_TRUE)); 389566063dSJacob Faibussowitsch PetscCall(TSTrajectorySetUp(tj,ts)); 399566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-check",&check,NULL)); 409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL)); 419566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetRealArray(NULL,NULL,"-interptimes",times,&Nt,NULL)); 42c4762a1bSJed Brown sort = PETSC_FALSE; 439566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-sorttimes",&sort,NULL)); 441baa6e33SBarry Smith if (sort) PetscCall(PetscSortReal(10,TT)); 45c4762a1bSJed Brown sort = PETSC_FALSE; 469566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-sortkeys",&sort,NULL)); 471baa6e33SBarry Smith if (sort) PetscCall(PetscSortInt(10,II)); 48c4762a1bSJed Brown p = PetscMax(p,-p); 49c4762a1bSJed Brown 50c4762a1bSJed Brown /* populate trajectory */ 51c4762a1bSJed Brown for (i = 0; i < 10; i++) { 529566063dSJacob Faibussowitsch PetscCall(VecSet(W,func(p,TT[i]))); 539566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts,II[i])); 549566063dSJacob Faibussowitsch PetscCall(TSTrajectorySet(tj,ts,II[i],TT[i],W)); 55c4762a1bSJed Brown } 56c4762a1bSJed Brown for (i = 0; i < Nt; i++) { 57c4762a1bSJed Brown PetscReal testtime = times[i], serr, derr; 58c4762a1bSJed Brown const PetscScalar *aW,*aWdot; 59c4762a1bSJed Brown 609566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGetVecs(tj,ts,PETSC_DECIDE,&testtime,W,Wdot)); 619566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(W,&aW)); 629566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Wdot,&aWdot)); 639566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," f(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(func(p,testtime)),(double)PetscRealPart(aW[0]))); 649566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"df(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(dfunc(p,testtime)),(double)PetscRealPart(aWdot[0]))); 65c4762a1bSJed Brown serr = PetscAbsScalar(func(p,testtime)-aW[0]); 66c4762a1bSJed Brown derr = PetscAbsScalar(dfunc(p,testtime)-aWdot[0]); 679566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(W,&aW)); 689566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Wdot,&aWdot)); 693c633725SBarry Smith PetscCheck(!check || serr <= tol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in state %g > %g",(double)serr,(double)tol); 703c633725SBarry Smith PetscCheck(!check || derr <= tol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in state derivative %g > %g",(double)derr,(double)tol); 71c4762a1bSJed Brown } 72c4762a1bSJed Brown for (i = Nt-1; i >= 0; i--) { 73c4762a1bSJed Brown PetscReal testtime = times[i], serr; 74c4762a1bSJed Brown const PetscScalar *aW; 75c4762a1bSJed Brown 769566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGetVecs(tj,ts,PETSC_DECIDE,&testtime,W,NULL)); 779566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(W,&aW)); 789566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," f(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(func(p,testtime)),(double)PetscRealPart(aW[0]))); 79c4762a1bSJed Brown serr = PetscAbsScalar(func(p,testtime)-aW[0]); 809566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(W,&aW)); 813c633725SBarry Smith PetscCheck(!check || serr <= tol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in state %g > %g",(double)serr,(double)tol); 82c4762a1bSJed Brown } 83c4762a1bSJed Brown for (i = Nt-1; i >= 0; i--) { 84c4762a1bSJed Brown PetscReal testtime = times[i], derr; 85c4762a1bSJed Brown const PetscScalar *aWdot; 86c4762a1bSJed Brown 879566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGetVecs(tj,ts,PETSC_DECIDE,&testtime,NULL,Wdot)); 889566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Wdot,&aWdot)); 899566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"df(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(dfunc(p,testtime)),(double)PetscRealPart(aWdot[0]))); 90c4762a1bSJed Brown derr = PetscAbsScalar(dfunc(p,testtime)-aWdot[0]); 919566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Wdot,&aWdot)); 923c633725SBarry Smith PetscCheck(!check || derr <= tol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in state derivative %g > %g",(double)derr,(double)tol); 93c4762a1bSJed Brown } 94c4762a1bSJed Brown for (i = 0; i < Nt; i++) { 95c4762a1bSJed Brown PetscReal testtime = times[i], serr, derr; 96c4762a1bSJed Brown const PetscScalar *aW,*aWdot; 97c4762a1bSJed Brown Vec hW,hWdot; 98c4762a1bSJed Brown 999566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGetUpdatedHistoryVecs(tj,ts,testtime,&hW,&hWdot)); 1009566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(hW,&aW)); 1019566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(hWdot,&aWdot)); 1029566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," f(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(func(p,testtime)),(double)PetscRealPart(aW[0]))); 1039566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"df(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(dfunc(p,testtime)),(double)PetscRealPart(aWdot[0]))); 104c4762a1bSJed Brown serr = PetscAbsScalar(func(p,testtime)-aW[0]); 105c4762a1bSJed Brown derr = PetscAbsScalar(dfunc(p,testtime)-aWdot[0]); 1069566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(hW,&aW)); 1079566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(hWdot,&aWdot)); 1089566063dSJacob Faibussowitsch PetscCall(TSTrajectoryRestoreUpdatedHistoryVecs(tj,&hW,&hWdot)); 1093c633725SBarry Smith PetscCheck(!check || serr <= tol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in state %g > %g",(double)serr,(double)tol); 1103c633725SBarry Smith PetscCheck(!check || derr <= tol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in state derivative %g > %g",(double)derr,(double)tol); 111c4762a1bSJed Brown } 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* Test on-the-fly reconstruction */ 1149566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1159566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 1169566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts,W2)); 1179566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts,10)); 1189566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts)); 1199566063dSJacob Faibussowitsch PetscCall(TSGetTrajectory(ts,&tj)); 1209566063dSJacob Faibussowitsch PetscCall(TSTrajectorySetType(tj,ts,TSTRAJECTORYBASIC)); 1219566063dSJacob Faibussowitsch PetscCall(TSTrajectorySetFromOptions(tj,ts)); 1229566063dSJacob Faibussowitsch PetscCall(TSTrajectorySetSolutionOnly(tj,PETSC_TRUE)); 1239566063dSJacob Faibussowitsch PetscCall(TSTrajectorySetUp(tj,ts)); 124c4762a1bSJed Brown use1 = PETSC_FALSE; 125c4762a1bSJed Brown use2 = PETSC_TRUE; 1269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-use_state",&use1,NULL)); 1279566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-use_der",&use2,NULL)); 1289566063dSJacob Faibussowitsch PetscCall(PetscSortReal(10,TT)); 129c4762a1bSJed Brown for (i = 0; i < 10; i++) { 130c4762a1bSJed Brown 1319566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts,i)); 1329566063dSJacob Faibussowitsch PetscCall(VecSet(W,func(p,TT[i]))); 1339566063dSJacob Faibussowitsch PetscCall(TSTrajectorySet(tj,ts,i,TT[i],W)); 134c4762a1bSJed Brown if (i) { 135c4762a1bSJed Brown const PetscScalar *aW,*aWdot; 136c4762a1bSJed Brown Vec hW,hWdot; 137c4762a1bSJed Brown PetscReal testtime = TT[i], serr, derr; 138c4762a1bSJed Brown 1399566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGetUpdatedHistoryVecs(tj,ts,testtime,use1 ? &hW : NULL,use2 ? &hWdot : NULL)); 140c4762a1bSJed Brown if (use1) { 1419566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(hW,&aW)); 1429566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," f(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(func(p,testtime)),(double)PetscRealPart(aW[0]))); 143c4762a1bSJed Brown serr = PetscAbsScalar(func(p,testtime)-aW[0]); 1449566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(hW,&aW)); 1453c633725SBarry Smith PetscCheck(!check || serr <= tol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in state %g > %g",(double)serr,(double)tol); 146c4762a1bSJed Brown } 147c4762a1bSJed Brown if (use2) { 1489566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(hWdot,&aWdot)); 1499566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"df(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(dfunc(p,testtime)),(double)PetscRealPart(aWdot[0]))); 150c4762a1bSJed Brown derr = PetscAbsScalar(dfunc(p,testtime)-aWdot[0]); 1519566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(hWdot,&aWdot)); 1523c633725SBarry Smith PetscCheck(!check || i < p || derr <= tol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in state derivative %g > %g",(double)derr,(double)tol); 153c4762a1bSJed Brown } 1549566063dSJacob Faibussowitsch PetscCall(TSTrajectoryRestoreUpdatedHistoryVecs(tj,use1 ? &hW : NULL,use2 ? &hWdot : NULL)); 155c4762a1bSJed Brown } 156c4762a1bSJed Brown } 1579566063dSJacob Faibussowitsch PetscCall(TSRemoveTrajectory(ts)); 1589566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&W)); 1609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&W2)); 1619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Wdot)); 1629566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 163b122ec5aSJacob Faibussowitsch return 0; 164c4762a1bSJed Brown } 165c4762a1bSJed Brown 166c4762a1bSJed Brown /*TEST 167c4762a1bSJed Brown 168c4762a1bSJed Brown test: 169c4762a1bSJed Brown suffix: 1 170c4762a1bSJed Brown requires: !single 171c4762a1bSJed Brown args: -ts_trajectory_monitor -p 3 -ts_trajectory_reconstruction_order 3 -interptimes 1,9.9,3,1.1,1.1,5.6 -check 172c4762a1bSJed Brown 173c4762a1bSJed Brown test: 174c4762a1bSJed Brown suffix: 2 175c4762a1bSJed Brown requires: !single 176c4762a1bSJed 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 177c4762a1bSJed Brown 178c4762a1bSJed Brown test: 179c4762a1bSJed Brown suffix: 3 180c4762a1bSJed Brown requires: !single 181c4762a1bSJed Brown args: -ts_trajectory_monitor -p 3 -ts_trajectory_reconstruction_order 5 -interptimes 1,9.9,3,1.1,1.1,5.6 -check 182c4762a1bSJed Brown 183c4762a1bSJed Brown TEST*/ 184