1 2 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 3 4 typedef struct { 5 PetscViewer viewer; 6 } TSTrajectory_Singlefile; 7 8 #undef __FUNCT__ 9 #define __FUNCT__ "TSTrajectorySet_Singlefile" 10 PetscErrorCode TSTrajectorySet_Singlefile(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X) 11 { 12 TSTrajectory_Singlefile *sf = (TSTrajectory_Singlefile*)tj->data; 13 PetscErrorCode ierr; 14 const char *filename; 15 16 PetscFunctionBegin; 17 if (stepnum == 0) { 18 ierr = PetscViewerCreate(PETSC_COMM_WORLD,&sf->viewer);CHKERRQ(ierr); 19 ierr = PetscViewerSetType(sf->viewer,PETSCVIEWERBINARY);CHKERRQ(ierr); 20 ierr = PetscViewerFileSetMode(sf->viewer,FILE_MODE_WRITE);CHKERRQ(ierr); 21 ierr = PetscObjectGetName((PetscObject)tj,&filename);CHKERRQ(ierr); 22 ierr = PetscViewerFileSetName(sf->viewer,filename);CHKERRQ(ierr); 23 } 24 ierr = VecView(X,sf->viewer);CHKERRQ(ierr); 25 ierr = PetscViewerBinaryWrite(sf->viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr); 26 PetscFunctionReturn(0); 27 } 28 29 #undef __FUNCT__ 30 #define __FUNCT__ "TSTrajectoryDestroy_Singlefile" 31 PetscErrorCode TSTrajectoryDestroy_Singlefile(TSTrajectory tj) 32 { 33 PetscErrorCode ierr; 34 TSTrajectory_Singlefile *sf = (TSTrajectory_Singlefile*)tj->data; 35 36 PetscFunctionBegin; 37 ierr = PetscViewerDestroy(&sf->viewer);CHKERRQ(ierr); 38 ierr = PetscFree(sf);CHKERRQ(ierr); 39 PetscFunctionReturn(0); 40 } 41 42 /*MC 43 TSTRAJECTORYSINGLEFILE - Stores all solutions of the ODE/ADE into a single file followed by each timestep. Does not save the intermediate stages in a multistage method 44 45 Level: intermediate 46 47 .seealso: TSTrajectoryCreate(), TS, TSTrajectorySetType() 48 49 M*/ 50 #undef __FUNCT__ 51 #define __FUNCT__ "TSTrajectoryCreate_Singlefile" 52 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Singlefile(TSTrajectory tj,TS ts) 53 { 54 PetscErrorCode ierr; 55 TSTrajectory_Singlefile *sf; 56 57 PetscFunctionBegin; 58 ierr = PetscNew(&sf);CHKERRQ(ierr); 59 tj->data = sf; 60 tj->ops->set = TSTrajectorySet_Singlefile; 61 tj->ops->get = NULL; 62 tj->ops->destroy = TSTrajectoryDestroy_Singlefile; 63 PetscFunctionReturn(0); 64 } 65