1bc952696SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 3bc952696SBarry Smith 4f3334a03SStefano Zampini typedef struct { 5f3334a03SStefano Zampini PetscViewer viewer; 6f3334a03SStefano Zampini } TSTrajectory_Basic; 7bc952696SBarry Smith 8560360afSLisandro Dalcin static PetscErrorCode TSTrajectorySet_Basic(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X) 9bc952696SBarry Smith { 10f3334a03SStefano Zampini TSTrajectory_Basic *tjbasic = (TSTrajectory_Basic*)tj->data; 119afe7f3eSBarry Smith char filename[PETSC_MAX_PATH_LEN]; 126affb6f8SHong Zhang PetscInt ns,i; 13bc952696SBarry Smith PetscErrorCode ierr; 146cd0a353SHong Zhang 1522a86ba0SHong Zhang PetscFunctionBegin; 169afe7f3eSBarry Smith ierr = PetscSNPrintf(filename,sizeof(filename),tj->dirfiletemplate,stepnum);CHKERRQ(ierr); 17ac1a7491SHong Zhang ierr = PetscViewerFileSetName(tjbasic->viewer,filename);CHKERRQ(ierr); /* this triggers PetscViewer to be set up again */ 18f3334a03SStefano Zampini ierr = PetscViewerSetUp(tjbasic->viewer);CHKERRQ(ierr); 19*cc4f23bcSHong Zhang if (!ts->stifflyaccurate) { 20f3334a03SStefano Zampini ierr = VecView(X,tjbasic->viewer);CHKERRQ(ierr); 21*cc4f23bcSHong Zhang } 22f253e43cSLisandro Dalcin ierr = PetscViewerBinaryWrite(tjbasic->viewer,&time,1,PETSC_REAL);CHKERRQ(ierr); 23f3334a03SStefano Zampini if (stepnum && !tj->solution_only) { 24f3334a03SStefano Zampini Vec *Y; 25f3334a03SStefano Zampini PetscReal tprev; 26bc952696SBarry Smith 27c679fc15SHong Zhang ierr = TSGetStages(ts,&ns,&Y);CHKERRQ(ierr); 28bc952696SBarry Smith for (i=0; i<ns; i++) { 29f3334a03SStefano Zampini ierr = VecView(Y[i],tjbasic->viewer);CHKERRQ(ierr); 30f3334a03SStefano Zampini } 31f3334a03SStefano Zampini ierr = TSGetPrevTime(ts,&tprev);CHKERRQ(ierr); 32f253e43cSLisandro Dalcin ierr = PetscViewerBinaryWrite(tjbasic->viewer,&tprev,1,PETSC_REAL);CHKERRQ(ierr); 33f3334a03SStefano Zampini } 346affb6f8SHong Zhang /* Tangent linear sensitivities needed by second-order adjoint */ 356affb6f8SHong Zhang if (ts->forward_solve) { 366affb6f8SHong Zhang Mat A,*S; 376affb6f8SHong Zhang 386affb6f8SHong Zhang ierr = TSForwardGetSensitivities(ts,NULL,&A);CHKERRQ(ierr); 396affb6f8SHong Zhang ierr = MatView(A,tjbasic->viewer);CHKERRQ(ierr); 406affb6f8SHong Zhang if (stepnum) { 416affb6f8SHong Zhang ierr = TSForwardGetStages(ts,&ns,&S);CHKERRQ(ierr); 426affb6f8SHong Zhang for (i=0; i<ns; i++) { 436affb6f8SHong Zhang ierr = MatView(S[i],tjbasic->viewer);CHKERRQ(ierr); 446affb6f8SHong Zhang } 456affb6f8SHong Zhang } 466affb6f8SHong Zhang } 47f3334a03SStefano Zampini PetscFunctionReturn(0); 48bc952696SBarry Smith } 49bc952696SBarry Smith 50f3334a03SStefano Zampini static PetscErrorCode TSTrajectorySetFromOptions_Basic(PetscOptionItems *PetscOptionsObject,TSTrajectory tj) 51f3334a03SStefano Zampini { 52f3334a03SStefano Zampini PetscErrorCode ierr; 53c679fc15SHong Zhang 54f3334a03SStefano Zampini PetscFunctionBegin; 55f3334a03SStefano Zampini ierr = PetscOptionsHead(PetscOptionsObject,"TS trajectory options for Basic type");CHKERRQ(ierr); 56f3334a03SStefano Zampini ierr = PetscOptionsTail();CHKERRQ(ierr); 57bc952696SBarry Smith PetscFunctionReturn(0); 58bc952696SBarry Smith } 59bc952696SBarry Smith 60560360afSLisandro Dalcin static PetscErrorCode TSTrajectoryGet_Basic(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal *t) 61bc952696SBarry Smith { 62bc952696SBarry Smith PetscViewer viewer; 639afe7f3eSBarry Smith char filename[PETSC_MAX_PATH_LEN]; 64bc952696SBarry Smith PetscErrorCode ierr; 65f3334a03SStefano Zampini Vec Sol; 666affb6f8SHong Zhang PetscInt ns,i; 676cd0a353SHong Zhang 6822a86ba0SHong Zhang PetscFunctionBegin; 69f3334a03SStefano Zampini ierr = PetscSNPrintf(filename,sizeof(filename),tj->dirfiletemplate,stepnum);CHKERRQ(ierr); 708a10d460SHong Zhang ierr = PetscViewerBinaryOpen(PetscObjectComm((PetscObject)tj),filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr); 71f3334a03SStefano Zampini ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE);CHKERRQ(ierr); 72*cc4f23bcSHong Zhang if (!ts->stifflyaccurate) { 73*cc4f23bcSHong Zhang ierr = TSGetSolution(ts,&Sol);CHKERRQ(ierr); 74bc952696SBarry Smith ierr = VecLoad(Sol,viewer);CHKERRQ(ierr); 75*cc4f23bcSHong Zhang } 76c679fc15SHong Zhang ierr = PetscViewerBinaryRead(viewer,t,1,NULL,PETSC_REAL);CHKERRQ(ierr); 77ac1a7491SHong Zhang if (stepnum && !tj->solution_only) { 78f3334a03SStefano Zampini Vec *Y; 79f3334a03SStefano Zampini PetscReal timepre; 80bc952696SBarry Smith 816affb6f8SHong Zhang ierr = TSGetStages(ts,&ns,&Y);CHKERRQ(ierr); 826affb6f8SHong Zhang for (i=0; i<ns; i++) { 83bc952696SBarry Smith ierr = VecLoad(Y[i],viewer);CHKERRQ(ierr); 84bc952696SBarry Smith } 85c679fc15SHong Zhang ierr = PetscViewerBinaryRead(viewer,&timepre,1,NULL,PETSC_REAL);CHKERRQ(ierr); 86f3334a03SStefano Zampini if (tj->adjoint_solve_mode) { 87acaebbb6SHong Zhang ierr = TSSetTimeStep(ts,-(*t)+timepre);CHKERRQ(ierr); 8826656371SHong Zhang } 89f3334a03SStefano Zampini } 906affb6f8SHong Zhang /* Tangent linear sensitivities needed by second-order adjoint */ 916affb6f8SHong Zhang if (ts->forward_solve) { 926affb6f8SHong Zhang 93*cc4f23bcSHong Zhang if (!ts->stifflyaccurate) { 94*cc4f23bcSHong Zhang Mat A; 956affb6f8SHong Zhang ierr = TSForwardGetSensitivities(ts,NULL,&A);CHKERRQ(ierr); 966affb6f8SHong Zhang ierr = MatLoad(A,viewer);CHKERRQ(ierr); 97*cc4f23bcSHong Zhang } 986affb6f8SHong Zhang if (stepnum) { 99*cc4f23bcSHong Zhang Mat *S; 1006affb6f8SHong Zhang ierr = TSForwardGetStages(ts,&ns,&S);CHKERRQ(ierr); 1016affb6f8SHong Zhang for (i=0; i<ns; i++) { 1026affb6f8SHong Zhang ierr = MatLoad(S[i],viewer);CHKERRQ(ierr); 1036affb6f8SHong Zhang } 1046affb6f8SHong Zhang } 1056affb6f8SHong Zhang } 10626656371SHong Zhang ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 107bc952696SBarry Smith PetscFunctionReturn(0); 108bc952696SBarry Smith } 109bc952696SBarry Smith 1108a10d460SHong Zhang PetscErrorCode TSTrajectorySetUp_Basic(TSTrajectory tj,TS ts) 111f3334a03SStefano Zampini { 112f3334a03SStefano Zampini MPI_Comm comm; 113f3334a03SStefano Zampini PetscMPIInt rank; 114f3334a03SStefano Zampini PetscErrorCode ierr; 115f3334a03SStefano Zampini 116f3334a03SStefano Zampini PetscFunctionBegin; 117f3334a03SStefano Zampini ierr = PetscObjectGetComm((PetscObject)tj,&comm);CHKERRQ(ierr); 118ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 119f3334a03SStefano Zampini if (!rank) { 1208a10d460SHong Zhang char *dir = tj->dirname; 121f3334a03SStefano Zampini PetscBool flg; 122f3334a03SStefano Zampini 1238a10d460SHong Zhang if (!dir) { 1248a10d460SHong Zhang char dtempname[16] = "TS-data-XXXXXX"; 1258a10d460SHong Zhang ierr = PetscMkdtemp(dtempname);CHKERRQ(ierr); 1268a10d460SHong Zhang ierr = PetscStrallocpy(dtempname,&tj->dirname);CHKERRQ(ierr); 1278a10d460SHong Zhang } else { 128f3334a03SStefano Zampini ierr = PetscTestDirectory(dir,'w',&flg);CHKERRQ(ierr); 129f3334a03SStefano Zampini if (!flg) { 130f3334a03SStefano Zampini ierr = PetscTestFile(dir,'r',&flg);CHKERRQ(ierr); 131f3334a03SStefano Zampini if (flg) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Specified path is a file - not a dir: %s",dir); 132f3334a03SStefano Zampini ierr = PetscMkdir(dir);CHKERRQ(ierr); 1338946c185SHong Zhang } else SETERRQ1(comm,PETSC_ERR_SUP,"Directory %s not empty",tj->dirname); 134f3334a03SStefano Zampini } 1358a10d460SHong Zhang } 136f3334a03SStefano Zampini ierr = PetscBarrier((PetscObject)tj);CHKERRQ(ierr); 137f3334a03SStefano Zampini PetscFunctionReturn(0); 138f3334a03SStefano Zampini } 139f3334a03SStefano Zampini 14064fc91eeSBarry Smith static PetscErrorCode TSTrajectoryDestroy_Basic(TSTrajectory tj) 14164fc91eeSBarry Smith { 142f3334a03SStefano Zampini TSTrajectory_Basic *tjbasic = (TSTrajectory_Basic*)tj->data; 14364fc91eeSBarry Smith PetscErrorCode ierr; 14464fc91eeSBarry Smith 14564fc91eeSBarry Smith PetscFunctionBegin; 146f3334a03SStefano Zampini ierr = PetscViewerDestroy(&tjbasic->viewer);CHKERRQ(ierr); 147f3334a03SStefano Zampini ierr = PetscFree(tjbasic);CHKERRQ(ierr); 14864fc91eeSBarry Smith PetscFunctionReturn(0); 14964fc91eeSBarry Smith } 15064fc91eeSBarry Smith 151bc952696SBarry Smith /*MC 15278fbdcc8SBarry Smith TSTRAJECTORYBASIC - Stores each solution of the ODE/DAE in a file 15378fbdcc8SBarry Smith 1548e5aa403SBarry Smith Saves each timestep into a separate file named TS-data-XXXXXX/TS-%06d.bin. The file name can be changed. 15578fbdcc8SBarry Smith 15678fbdcc8SBarry Smith This version saves the solutions at all the stages 15778fbdcc8SBarry Smith 15878fbdcc8SBarry Smith $PETSC_DIR/share/petsc/matlab/PetscReadBinaryTrajectory.m can read in files created with this format 159bc952696SBarry Smith 160bc952696SBarry Smith Level: intermediate 161bc952696SBarry Smith 16264e38db7SHong Zhang .seealso: TSTrajectoryCreate(), TS, TSTrajectorySetType(), TSTrajectorySetDirname(), TSTrajectorySetFile() 163bc952696SBarry Smith 164bc952696SBarry Smith M*/ 165c3df6c96SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Basic(TSTrajectory tj,TS ts) 166bc952696SBarry Smith { 167f3334a03SStefano Zampini TSTrajectory_Basic *tjbasic; 168f3334a03SStefano Zampini PetscErrorCode ierr; 169f3334a03SStefano Zampini 170bc952696SBarry Smith PetscFunctionBegin; 171f3334a03SStefano Zampini ierr = PetscNew(&tjbasic);CHKERRQ(ierr); 172f3334a03SStefano Zampini 173f3334a03SStefano Zampini ierr = PetscViewerCreate(PetscObjectComm((PetscObject)tj),&tjbasic->viewer);CHKERRQ(ierr); 174f3334a03SStefano Zampini ierr = PetscViewerSetType(tjbasic->viewer,PETSCVIEWERBINARY);CHKERRQ(ierr); 175f3334a03SStefano Zampini ierr = PetscViewerPushFormat(tjbasic->viewer,PETSC_VIEWER_NATIVE);CHKERRQ(ierr); 176f3334a03SStefano Zampini ierr = PetscViewerFileSetMode(tjbasic->viewer,FILE_MODE_WRITE);CHKERRQ(ierr); 177f3334a03SStefano Zampini tj->data = tjbasic; 178f3334a03SStefano Zampini 17952fbb9a2SHong Zhang tj->ops->set = TSTrajectorySet_Basic; 18052fbb9a2SHong Zhang tj->ops->get = TSTrajectoryGet_Basic; 181f3334a03SStefano Zampini tj->ops->setup = TSTrajectorySetUp_Basic; 18264fc91eeSBarry Smith tj->ops->destroy = TSTrajectoryDestroy_Basic; 183f3334a03SStefano Zampini tj->ops->setfromoptions = TSTrajectorySetFromOptions_Basic; 184bc952696SBarry Smith PetscFunctionReturn(0); 185bc952696SBarry Smith } 186