Actual source code: singlefile.c


  2: #include <petsc/private/tsimpl.h>

  4: typedef struct {
  5:   PetscViewer viewer;
  6: } TSTrajectory_Singlefile;

  8: static PetscErrorCode TSTrajectorySet_Singlefile(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X)
  9: {
 10:   TSTrajectory_Singlefile *sf = (TSTrajectory_Singlefile*)tj->data;
 11:   const char              *filename;

 13:   if (stepnum == 0) {
 14:     PetscViewerCreate(PetscObjectComm((PetscObject)X),&sf->viewer);
 15:     PetscViewerSetType(sf->viewer,PETSCVIEWERBINARY);
 16:     PetscViewerFileSetMode(sf->viewer,FILE_MODE_WRITE);
 17:     PetscObjectGetName((PetscObject)tj,&filename);
 18:     PetscViewerFileSetName(sf->viewer,filename);
 19:   }
 20:   VecView(X,sf->viewer);
 21:   PetscViewerBinaryWrite(sf->viewer,&time,1,PETSC_REAL);
 22:   return 0;
 23: }

 25: static PetscErrorCode TSTrajectoryDestroy_Singlefile(TSTrajectory tj)
 26: {
 27:   TSTrajectory_Singlefile *sf = (TSTrajectory_Singlefile*)tj->data;

 29:   PetscViewerDestroy(&sf->viewer);
 30:   PetscFree(sf);
 31:   return 0;
 32: }

 34: /*MC
 35:       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

 37:   Level: intermediate

 39: .seealso:  TSTrajectoryCreate(), TS, TSTrajectorySetType()

 41: M*/
 42: PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Singlefile(TSTrajectory tj,TS ts)
 43: {
 44:   TSTrajectory_Singlefile *sf;

 46:   PetscNew(&sf);
 47:   tj->data         = sf;
 48:   tj->ops->set     = TSTrajectorySet_Singlefile;
 49:   tj->ops->get     = NULL;
 50:   tj->ops->destroy = TSTrajectoryDestroy_Singlefile;
 51:   ts->setupcalled  = PETSC_TRUE;
 52:   return 0;
 53: }