1 2 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "OutputBIN" 6 static PetscErrorCode OutputBIN(const char *filename,PetscViewer *viewer) 7 { 8 PetscErrorCode ierr; 9 10 PetscFunctionBegin; 11 ierr = PetscViewerCreate(PETSC_COMM_WORLD,viewer);CHKERRQ(ierr); 12 ierr = PetscViewerSetType(*viewer,PETSCVIEWERBINARY);CHKERRQ(ierr); 13 ierr = PetscViewerFileSetMode(*viewer,FILE_MODE_WRITE);CHKERRQ(ierr); 14 ierr = PetscViewerFileSetName(*viewer,filename);CHKERRQ(ierr); 15 PetscFunctionReturn(0); 16 } 17 18 19 #undef __FUNCT__ 20 #define __FUNCT__ "TSTrajectorySet_Basic" 21 PetscErrorCode TSTrajectorySet_Basic(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X) 22 { 23 PetscViewer viewer; 24 PetscInt ns,i; 25 Vec *Y; 26 char filename[PETSC_MAX_PATH_LEN]; 27 PetscReal tprev; 28 PetscErrorCode ierr; 29 30 PetscFunctionBegin; 31 ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr); 32 if (stepnum == 0) { 33 #if defined(PETSC_HAVE_POPEN) 34 PetscMPIInt rank; 35 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);CHKERRQ(ierr); 36 if (!rank) { 37 char command[PETSC_MAX_PATH_LEN]; 38 FILE *fd; 39 int err; 40 41 ierr = PetscMemzero(command,sizeof(command));CHKERRQ(ierr); 42 ierr = PetscSNPrintf(command,PETSC_MAX_PATH_LEN,"rm -fr %s","SA-data");CHKERRQ(ierr); 43 ierr = PetscPOpen(PETSC_COMM_SELF,NULL,command,"r",&fd);CHKERRQ(ierr); 44 ierr = PetscPClose(PETSC_COMM_SELF,fd,&err);CHKERRQ(ierr); 45 ierr = PetscSNPrintf(command,PETSC_MAX_PATH_LEN,"mkdir %s","SA-data");CHKERRQ(ierr); 46 ierr = PetscPOpen(PETSC_COMM_SELF,NULL,command,"r",&fd);CHKERRQ(ierr); 47 ierr = PetscPClose(PETSC_COMM_SELF,fd,&err);CHKERRQ(ierr); 48 } 49 } 50 #endif 51 PetscFunctionReturn(0); 52 } 53 ierr = PetscSNPrintf(filename,sizeof(filename),"SA-data/SA-%06d.bin",stepnum);CHKERRQ(ierr); 54 ierr = OutputBIN(filename,&viewer);CHKERRQ(ierr); 55 ierr = VecView(X,viewer);CHKERRQ(ierr); 56 ierr = PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr); 57 58 ierr = TSGetStages(ts,&ns,&Y);CHKERRQ(ierr); 59 for (i=0;i<ns;i++) { 60 ierr = VecView(Y[i],viewer);CHKERRQ(ierr); 61 } 62 63 ierr = TSGetPrevTime(ts,&tprev);CHKERRQ(ierr); 64 ierr = PetscViewerBinaryWrite(viewer,&tprev,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr); 65 66 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 67 PetscFunctionReturn(0); 68 } 69 70 #undef __FUNCT__ 71 #define __FUNCT__ "TSTrajectoryGet_Basic" 72 PetscErrorCode TSTrajectoryGet_Basic(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal *t) 73 { 74 Vec Sol,*Y; 75 PetscInt Nr,i; 76 PetscViewer viewer; 77 PetscReal timepre; 78 char filename[PETSC_MAX_PATH_LEN]; 79 PetscErrorCode ierr; 80 81 PetscFunctionBegin; 82 ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr); 83 ierr = PetscSNPrintf(filename,sizeof filename,"SA-data/SA-%06d.bin",stepnum);CHKERRQ(ierr); 84 ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr); 85 86 ierr = TSGetSolution(ts,&Sol);CHKERRQ(ierr); 87 ierr = VecLoad(Sol,viewer);CHKERRQ(ierr); 88 89 ierr = PetscViewerBinaryRead(viewer,t,1,NULL,PETSC_REAL);CHKERRQ(ierr); 90 91 ierr = TSGetStages(ts,&Nr,&Y);CHKERRQ(ierr); 92 for (i=0;i<Nr ;i++) { 93 ierr = VecLoad(Y[i],viewer);CHKERRQ(ierr); 94 } 95 96 /* ierr = PetscRealLoad(Nr,&Nr,&timepre,viewer);CHKERRQ(ierr); */ 97 ierr = PetscViewerBinaryRead(viewer,&timepre,1,NULL,PETSC_REAL);CHKERRQ(ierr); 98 99 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 100 101 ierr = TSSetTimeStep(ts,-(*t)+timepre);CHKERRQ(ierr); 102 103 PetscFunctionReturn(0); 104 } 105 106 /*MC 107 TSTRAJECTORYBASIC - Stores each solution of the ODE/ADE in a file 108 109 Level: intermediate 110 111 .seealso: TSTrajectoryCreate(), TS, TSTrajectorySetType() 112 113 M*/ 114 #undef __FUNCT__ 115 #define __FUNCT__ "TSTrajectoryCreate_Basic" 116 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Basic(TSTrajectory tj,TS ts) 117 { 118 PetscFunctionBegin; 119 tj->ops->set = TSTrajectorySet_Basic; 120 tj->ops->get = TSTrajectoryGet_Basic; 121 PetscFunctionReturn(0); 122 } 123