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]; 12*6affb6f8SHong 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); 19f3334a03SStefano Zampini ierr = VecView(X,tjbasic->viewer);CHKERRQ(ierr); 20f3334a03SStefano Zampini ierr = PetscViewerBinaryWrite(tjbasic->viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr); 21f3334a03SStefano Zampini if (stepnum && !tj->solution_only) { 22f3334a03SStefano Zampini Vec *Y; 23f3334a03SStefano Zampini PetscReal tprev; 24bc952696SBarry Smith 25c679fc15SHong Zhang ierr = TSGetStages(ts,&ns,&Y);CHKERRQ(ierr); 26bc952696SBarry Smith for (i=0; i<ns; i++) { 27f3334a03SStefano Zampini ierr = VecView(Y[i],tjbasic->viewer);CHKERRQ(ierr); 28f3334a03SStefano Zampini } 29f3334a03SStefano Zampini ierr = TSGetPrevTime(ts,&tprev);CHKERRQ(ierr); 30f3334a03SStefano Zampini ierr = PetscViewerBinaryWrite(tjbasic->viewer,&tprev,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr); 31f3334a03SStefano Zampini } 32*6affb6f8SHong Zhang /* Tangent linear sensitivities needed by second-order adjoint */ 33*6affb6f8SHong Zhang if (ts->forward_solve) { 34*6affb6f8SHong Zhang Mat A,*S; 35*6affb6f8SHong Zhang 36*6affb6f8SHong Zhang ierr = TSForwardGetSensitivities(ts,NULL,&A);CHKERRQ(ierr); 37*6affb6f8SHong Zhang ierr = MatView(A,tjbasic->viewer);CHKERRQ(ierr); 38*6affb6f8SHong Zhang if (stepnum) { 39*6affb6f8SHong Zhang ierr = TSForwardGetStages(ts,&ns,&S);CHKERRQ(ierr); 40*6affb6f8SHong Zhang for (i=0; i<ns; i++) { 41*6affb6f8SHong Zhang ierr = MatView(S[i],tjbasic->viewer);CHKERRQ(ierr); 42*6affb6f8SHong Zhang } 43*6affb6f8SHong Zhang } 44*6affb6f8SHong Zhang } 45f3334a03SStefano Zampini PetscFunctionReturn(0); 46bc952696SBarry Smith } 47bc952696SBarry Smith 48f3334a03SStefano Zampini static PetscErrorCode TSTrajectorySetFromOptions_Basic(PetscOptionItems *PetscOptionsObject,TSTrajectory tj) 49f3334a03SStefano Zampini { 50f3334a03SStefano Zampini PetscErrorCode ierr; 51c679fc15SHong Zhang 52f3334a03SStefano Zampini PetscFunctionBegin; 53f3334a03SStefano Zampini ierr = PetscOptionsHead(PetscOptionsObject,"TS trajectory options for Basic type");CHKERRQ(ierr); 54f3334a03SStefano Zampini ierr = PetscOptionsTail();CHKERRQ(ierr); 55bc952696SBarry Smith PetscFunctionReturn(0); 56bc952696SBarry Smith } 57bc952696SBarry Smith 58560360afSLisandro Dalcin static PetscErrorCode TSTrajectoryGet_Basic(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal *t) 59bc952696SBarry Smith { 60bc952696SBarry Smith PetscViewer viewer; 619afe7f3eSBarry Smith char filename[PETSC_MAX_PATH_LEN]; 62bc952696SBarry Smith PetscErrorCode ierr; 63f3334a03SStefano Zampini Vec Sol; 64*6affb6f8SHong Zhang PetscInt ns,i; 656cd0a353SHong Zhang 6622a86ba0SHong Zhang PetscFunctionBegin; 67f3334a03SStefano Zampini ierr = PetscSNPrintf(filename,sizeof(filename),tj->dirfiletemplate,stepnum);CHKERRQ(ierr); 68bc952696SBarry Smith ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr); 69bc952696SBarry Smith ierr = TSGetSolution(ts,&Sol);CHKERRQ(ierr); 70f3334a03SStefano Zampini ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE);CHKERRQ(ierr); 71bc952696SBarry Smith ierr = VecLoad(Sol,viewer);CHKERRQ(ierr); 72c679fc15SHong Zhang ierr = PetscViewerBinaryRead(viewer,t,1,NULL,PETSC_REAL);CHKERRQ(ierr); 73ac1a7491SHong Zhang if (stepnum && !tj->solution_only) { 74f3334a03SStefano Zampini Vec *Y; 75f3334a03SStefano Zampini PetscReal timepre; 76bc952696SBarry Smith 77*6affb6f8SHong Zhang ierr = TSGetStages(ts,&ns,&Y);CHKERRQ(ierr); 78*6affb6f8SHong Zhang for (i=0; i<ns; i++) { 79bc952696SBarry Smith ierr = VecLoad(Y[i],viewer);CHKERRQ(ierr); 80bc952696SBarry Smith } 81c679fc15SHong Zhang ierr = PetscViewerBinaryRead(viewer,&timepre,1,NULL,PETSC_REAL);CHKERRQ(ierr); 82f3334a03SStefano Zampini if (tj->adjoint_solve_mode) { 83acaebbb6SHong Zhang ierr = TSSetTimeStep(ts,-(*t)+timepre);CHKERRQ(ierr); 8426656371SHong Zhang } 85f3334a03SStefano Zampini } 86*6affb6f8SHong Zhang /* Tangent linear sensitivities needed by second-order adjoint */ 87*6affb6f8SHong Zhang if (ts->forward_solve) { 88*6affb6f8SHong Zhang Mat A,*S; 89*6affb6f8SHong Zhang 90*6affb6f8SHong Zhang ierr = TSForwardGetSensitivities(ts,NULL,&A);CHKERRQ(ierr); 91*6affb6f8SHong Zhang ierr = MatLoad(A,viewer);CHKERRQ(ierr); 92*6affb6f8SHong Zhang if (stepnum) { 93*6affb6f8SHong Zhang ierr = TSForwardGetStages(ts,&ns,&S);CHKERRQ(ierr); 94*6affb6f8SHong Zhang for (i=0; i<ns; i++) { 95*6affb6f8SHong Zhang ierr = MatLoad(S[i],viewer);CHKERRQ(ierr); 96*6affb6f8SHong Zhang } 97*6affb6f8SHong Zhang } 98*6affb6f8SHong Zhang } 9926656371SHong Zhang ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 100bc952696SBarry Smith PetscFunctionReturn(0); 101bc952696SBarry Smith } 102bc952696SBarry Smith 103f3334a03SStefano Zampini static PetscErrorCode TSTrajectorySetUp_Basic(TSTrajectory tj,TS ts) 104f3334a03SStefano Zampini { 105f3334a03SStefano Zampini MPI_Comm comm; 106f3334a03SStefano Zampini PetscMPIInt rank; 107f3334a03SStefano Zampini PetscErrorCode ierr; 108f3334a03SStefano Zampini 109f3334a03SStefano Zampini PetscFunctionBegin; 110f3334a03SStefano Zampini ierr = PetscObjectGetComm((PetscObject)tj,&comm);CHKERRQ(ierr); 111f3334a03SStefano Zampini ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 112f3334a03SStefano Zampini if (!rank) { 113f3334a03SStefano Zampini const char* dir = tj->dirname; 114f3334a03SStefano Zampini PetscBool flg; 115f3334a03SStefano Zampini 116f3334a03SStefano Zampini /* I don't like running PetscRMTree on a directory */ 117f3334a03SStefano Zampini ierr = PetscTestDirectory(dir,'w',&flg);CHKERRQ(ierr); 118f3334a03SStefano Zampini if (!flg) { 119f3334a03SStefano Zampini ierr = PetscTestFile(dir,'r',&flg);CHKERRQ(ierr); 120f3334a03SStefano Zampini if (flg) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Specified path is a file - not a dir: %s",dir); 121f3334a03SStefano Zampini ierr = PetscMkdir(dir);CHKERRQ(ierr); 122ac1a7491SHong Zhang } else { 123ac1a7491SHong Zhang ierr = PetscRMTree(tj->dirname);CHKERRQ(ierr); /* avoid having to delete the folder manually */ 124ac1a7491SHong Zhang } 125f3334a03SStefano Zampini } 126f3334a03SStefano Zampini ierr = PetscBarrier((PetscObject)tj);CHKERRQ(ierr); 127f3334a03SStefano Zampini PetscFunctionReturn(0); 128f3334a03SStefano Zampini } 129f3334a03SStefano Zampini 13064fc91eeSBarry Smith static PetscErrorCode TSTrajectoryDestroy_Basic(TSTrajectory tj) 13164fc91eeSBarry Smith { 132f3334a03SStefano Zampini TSTrajectory_Basic *tjbasic = (TSTrajectory_Basic*)tj->data; 13364fc91eeSBarry Smith PetscErrorCode ierr; 13464fc91eeSBarry Smith 13564fc91eeSBarry Smith PetscFunctionBegin; 136f3334a03SStefano Zampini ierr = PetscViewerDestroy(&tjbasic->viewer);CHKERRQ(ierr); 137f3334a03SStefano Zampini ierr = PetscFree(tjbasic);CHKERRQ(ierr); 13864fc91eeSBarry Smith PetscFunctionReturn(0); 13964fc91eeSBarry Smith } 14064fc91eeSBarry Smith 141bc952696SBarry Smith /*MC 14278fbdcc8SBarry Smith TSTRAJECTORYBASIC - Stores each solution of the ODE/DAE in a file 14378fbdcc8SBarry Smith 14464e38db7SHong Zhang Saves each timestep into a seperate file named SA-data/SA-%06d.bin. The file name can be changed. 14578fbdcc8SBarry Smith 14678fbdcc8SBarry Smith This version saves the solutions at all the stages 14778fbdcc8SBarry Smith 14878fbdcc8SBarry Smith $PETSC_DIR/share/petsc/matlab/PetscReadBinaryTrajectory.m can read in files created with this format 149bc952696SBarry Smith 150bc952696SBarry Smith Level: intermediate 151bc952696SBarry Smith 15264e38db7SHong Zhang .seealso: TSTrajectoryCreate(), TS, TSTrajectorySetType(), TSTrajectorySetDirname(), TSTrajectorySetFile() 153bc952696SBarry Smith 154bc952696SBarry Smith M*/ 155c3df6c96SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Basic(TSTrajectory tj,TS ts) 156bc952696SBarry Smith { 157f3334a03SStefano Zampini TSTrajectory_Basic *tjbasic; 158f3334a03SStefano Zampini PetscErrorCode ierr; 159f3334a03SStefano Zampini 160bc952696SBarry Smith PetscFunctionBegin; 161f3334a03SStefano Zampini ierr = PetscNew(&tjbasic);CHKERRQ(ierr); 162f3334a03SStefano Zampini 163f3334a03SStefano Zampini ierr = PetscViewerCreate(PetscObjectComm((PetscObject)tj),&tjbasic->viewer);CHKERRQ(ierr); 164f3334a03SStefano Zampini ierr = PetscViewerSetType(tjbasic->viewer,PETSCVIEWERBINARY);CHKERRQ(ierr); 165f3334a03SStefano Zampini ierr = PetscViewerPushFormat(tjbasic->viewer,PETSC_VIEWER_NATIVE);CHKERRQ(ierr); 166f3334a03SStefano Zampini ierr = PetscViewerFileSetMode(tjbasic->viewer,FILE_MODE_WRITE);CHKERRQ(ierr); 167f3334a03SStefano Zampini tj->data = tjbasic; 168f3334a03SStefano Zampini 16952fbb9a2SHong Zhang tj->ops->set = TSTrajectorySet_Basic; 17052fbb9a2SHong Zhang tj->ops->get = TSTrajectoryGet_Basic; 171f3334a03SStefano Zampini tj->ops->setup = TSTrajectorySetUp_Basic; 17264fc91eeSBarry Smith tj->ops->destroy = TSTrajectoryDestroy_Basic; 173f3334a03SStefano Zampini tj->ops->setfromoptions = TSTrajectorySetFromOptions_Basic; 174bc952696SBarry Smith PetscFunctionReturn(0); 175bc952696SBarry Smith } 176