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*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 245f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&W)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(W,1,PETSC_DECIDE)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetUp(W)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(W,&Wdot)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(W,&W2)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSolution(ts,W2)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxSteps(ts,10)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSaveTrajectory(ts)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetTrajectory(ts,&tj)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(TSTrajectorySetType(tj,ts,TSTRAJECTORYBASIC)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(TSTrajectorySetFromOptions(tj,ts)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(TSTrajectorySetSolutionOnly(tj,PETSC_TRUE)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(TSTrajectorySetUp(tj,ts)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-check",&check,NULL)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetRealArray(NULL,NULL,"-interptimes",times,&Nt,NULL)); 41c4762a1bSJed Brown sort = PETSC_FALSE; 425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-sorttimes",&sort,NULL)); 43c4762a1bSJed Brown if (sort) { 445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSortReal(10,TT)); 45c4762a1bSJed Brown } 46c4762a1bSJed Brown sort = PETSC_FALSE; 475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-sortkeys",&sort,NULL)); 48c4762a1bSJed Brown if (sort) { 495f80ce2aSJacob Faibussowitsch CHKERRQ(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++) { 555f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(W,func(p,TT[i]))); 565f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetStepNumber(ts,II[i])); 575f80ce2aSJacob Faibussowitsch CHKERRQ(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 635f80ce2aSJacob Faibussowitsch CHKERRQ(TSTrajectoryGetVecs(tj,ts,PETSC_DECIDE,&testtime,W,Wdot)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(W,&aW)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(Wdot,&aWdot)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," f(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(func(p,testtime)),(double)PetscRealPart(aW[0]))); 675f80ce2aSJacob Faibussowitsch CHKERRQ(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]); 705f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(W,&aW)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(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 795f80ce2aSJacob Faibussowitsch CHKERRQ(TSTrajectoryGetVecs(tj,ts,PETSC_DECIDE,&testtime,W,NULL)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(W,&aW)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(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]); 835f80ce2aSJacob Faibussowitsch CHKERRQ(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 905f80ce2aSJacob Faibussowitsch CHKERRQ(TSTrajectoryGetVecs(tj,ts,PETSC_DECIDE,&testtime,NULL,Wdot)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(Wdot,&aWdot)); 925f80ce2aSJacob Faibussowitsch CHKERRQ(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]); 945f80ce2aSJacob Faibussowitsch CHKERRQ(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 1025f80ce2aSJacob Faibussowitsch CHKERRQ(TSTrajectoryGetUpdatedHistoryVecs(tj,ts,testtime,&hW,&hWdot)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(hW,&aW)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(hWdot,&aWdot)); 1055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," f(%g) = %g (reconstructed %g)\n",(double)testtime,(double)PetscRealPart(func(p,testtime)),(double)PetscRealPart(aW[0]))); 1065f80ce2aSJacob Faibussowitsch CHKERRQ(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]); 1095f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(hW,&aW)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(hWdot,&aWdot)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 1175f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 1195f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSolution(ts,W2)); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxSteps(ts,10)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSaveTrajectory(ts)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetTrajectory(ts,&tj)); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(TSTrajectorySetType(tj,ts,TSTRAJECTORYBASIC)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(TSTrajectorySetFromOptions(tj,ts)); 1255f80ce2aSJacob Faibussowitsch CHKERRQ(TSTrajectorySetSolutionOnly(tj,PETSC_TRUE)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(TSTrajectorySetUp(tj,ts)); 127c4762a1bSJed Brown use1 = PETSC_FALSE; 128c4762a1bSJed Brown use2 = PETSC_TRUE; 1295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-use_state",&use1,NULL)); 1305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-use_der",&use2,NULL)); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSortReal(10,TT)); 132c4762a1bSJed Brown for (i = 0; i < 10; i++) { 133c4762a1bSJed Brown 1345f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetStepNumber(ts,i)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(W,func(p,TT[i]))); 1365f80ce2aSJacob Faibussowitsch CHKERRQ(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 1425f80ce2aSJacob Faibussowitsch CHKERRQ(TSTrajectoryGetUpdatedHistoryVecs(tj,ts,testtime,use1 ? &hW : NULL,use2 ? &hWdot : NULL)); 143c4762a1bSJed Brown if (use1) { 1445f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(hW,&aW)); 1455f80ce2aSJacob Faibussowitsch CHKERRQ(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]); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(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) { 1515f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(hWdot,&aWdot)); 1525f80ce2aSJacob Faibussowitsch CHKERRQ(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]); 1545f80ce2aSJacob Faibussowitsch CHKERRQ(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 } 1575f80ce2aSJacob Faibussowitsch CHKERRQ(TSTrajectoryRestoreUpdatedHistoryVecs(tj,use1 ? &hW : NULL,use2 ? &hWdot : NULL)); 158c4762a1bSJed Brown } 159c4762a1bSJed Brown } 1605f80ce2aSJacob Faibussowitsch CHKERRQ(TSRemoveTrajectory(ts)); 1615f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 1625f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&W)); 1635f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&W2)); 1645f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Wdot)); 165*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 166*b122ec5aSJacob 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