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