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 #endif 50 PetscFunctionReturn(0); 51 } 52 ierr = PetscSNPrintf(filename,sizeof(filename),"SA-data/SA-%06d.bin",stepnum);CHKERRQ(ierr); 53 ierr = OutputBIN(filename,&viewer);CHKERRQ(ierr); 54 ierr = VecView(X,viewer);CHKERRQ(ierr); 55 ierr = PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr); 56 57 ierr = TSGetStages(ts,&ns,&Y);CHKERRQ(ierr); 58 for (i=0;i<ns;i++) { 59 ierr = VecView(Y[i],viewer);CHKERRQ(ierr); 60 } 61 62 ierr = TSGetPrevTime(ts,&tprev);CHKERRQ(ierr); 63 ierr = PetscViewerBinaryWrite(viewer,&tprev,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr); 64 65 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 66 PetscFunctionReturn(0); 67 } 68 69 #undef __FUNCT__ 70 #define __FUNCT__ "TSTrajectoryGet_Basic" 71 PetscErrorCode TSTrajectoryGet_Basic(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal *t) 72 { 73 Vec Sol,*Y; 74 PetscInt Nr,i; 75 PetscViewer viewer; 76 PetscReal timepre; 77 char filename[PETSC_MAX_PATH_LEN]; 78 PetscErrorCode ierr; 79 80 PetscFunctionBegin; 81 ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr); 82 ierr = PetscSNPrintf(filename,sizeof filename,"SA-data/SA-%06d.bin",stepnum);CHKERRQ(ierr); 83 ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr); 84 85 ierr = TSGetSolution(ts,&Sol);CHKERRQ(ierr); 86 ierr = VecLoad(Sol,viewer);CHKERRQ(ierr); 87 88 ierr = PetscViewerBinaryRead(viewer,t,1,NULL,PETSC_REAL);CHKERRQ(ierr); 89 90 ierr = TSGetStages(ts,&Nr,&Y);CHKERRQ(ierr); 91 for (i=0;i<Nr ;i++) { 92 ierr = VecLoad(Y[i],viewer);CHKERRQ(ierr); 93 } 94 95 /* ierr = PetscRealLoad(Nr,&Nr,&timepre,viewer);CHKERRQ(ierr); */ 96 ierr = PetscViewerBinaryRead(viewer,&timepre,1,NULL,PETSC_REAL);CHKERRQ(ierr); 97 98 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 99 100 ierr = TSSetTimeStep(ts,-(*t)+timepre);CHKERRQ(ierr); 101 102 PetscFunctionReturn(0); 103 } 104 105 /*MC 106 TSTRAJECTORYBASIC - Stores each solution of the ODE/ADE in a file 107 108 Level: intermediate 109 110 .seealso: TSTrajectoryCreate(), TS, TSTrajectorySetType() 111 112 M*/ 113 #undef __FUNCT__ 114 #define __FUNCT__ "TSTrajectoryCreate_Basic" 115 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Basic(TSTrajectory tj,TS ts) 116 { 117 PetscFunctionBegin; 118 tj->ops->set = TSTrajectorySet_Basic; 119 tj->ops->get = TSTrajectoryGet_Basic; 120 PetscFunctionReturn(0); 121 } 122