xref: /petsc/src/ts/trajectory/impls/singlefile/singlefile.c (revision 972caf09b6a954926d5275ef5d098d7c49d89ccf)
11c8c567eSBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/tsimpl.h>        /*I "petscts.h"  I*/
31c8c567eSBarry Smith 
41c8c567eSBarry Smith typedef struct {
51c8c567eSBarry Smith   PetscViewer viewer;
61c8c567eSBarry Smith } TSTrajectory_Singlefile;
71c8c567eSBarry Smith 
81c8c567eSBarry Smith #undef __FUNCT__
91c8c567eSBarry Smith #define __FUNCT__ "TSTrajectorySet_Singlefile"
10*972caf09SHong Zhang PetscErrorCode TSTrajectorySet_Singlefile(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X)
111c8c567eSBarry Smith {
12*972caf09SHong Zhang   TSTrajectory_Singlefile *sf = (TSTrajectory_Singlefile*)tj->data;
131c8c567eSBarry Smith   PetscErrorCode          ierr;
141c8c567eSBarry Smith   const char              *filename;
151c8c567eSBarry Smith 
161c8c567eSBarry Smith   PetscFunctionBeginUser;
171c8c567eSBarry Smith   if (stepnum == 0) {
181c8c567eSBarry Smith     ierr = PetscViewerCreate(PETSC_COMM_WORLD, &sf->viewer);CHKERRQ(ierr);
191c8c567eSBarry Smith     ierr = PetscViewerSetType(sf->viewer, PETSCVIEWERBINARY);CHKERRQ(ierr);
201c8c567eSBarry Smith     ierr = PetscViewerFileSetMode(sf->viewer,FILE_MODE_WRITE);CHKERRQ(ierr);
21*972caf09SHong Zhang     ierr = PetscObjectGetName((PetscObject)tj,&filename);CHKERRQ(ierr);
221c8c567eSBarry Smith     ierr = PetscViewerFileSetName(sf->viewer, filename);CHKERRQ(ierr);
231c8c567eSBarry Smith   }
24c679fc15SHong Zhang   ierr = VecView(X,sf->viewer);CHKERRQ(ierr);
25c679fc15SHong Zhang   ierr = PetscViewerBinaryWrite(sf->viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
261c8c567eSBarry Smith   PetscFunctionReturn(0);
271c8c567eSBarry Smith }
281c8c567eSBarry Smith 
291c8c567eSBarry Smith #undef __FUNCT__
301c8c567eSBarry Smith #define __FUNCT__ "TSTrajectoryDestroy_Singlefile"
31*972caf09SHong Zhang PetscErrorCode TSTrajectoryDestroy_Singlefile(TSTrajectory tj)
321c8c567eSBarry Smith {
331c8c567eSBarry Smith   PetscErrorCode          ierr;
34*972caf09SHong Zhang   TSTrajectory_Singlefile *sf = (TSTrajectory_Singlefile*)tj->data;
351c8c567eSBarry Smith 
361c8c567eSBarry Smith   PetscFunctionBegin;
371c8c567eSBarry Smith   ierr = PetscViewerDestroy(&sf->viewer);CHKERRQ(ierr);
381c8c567eSBarry Smith   ierr = PetscFree(sf);CHKERRQ(ierr);
391c8c567eSBarry Smith   PetscFunctionReturn(0);
401c8c567eSBarry Smith }
411c8c567eSBarry Smith 
421c8c567eSBarry Smith /*MC
438868bf69SBarry Smith       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
441c8c567eSBarry Smith 
451c8c567eSBarry Smith   Level: intermediate
461c8c567eSBarry Smith 
471c8c567eSBarry Smith .seealso:  TSTrajectoryCreate(), TS, TSTrajectorySetType()
481c8c567eSBarry Smith 
491c8c567eSBarry Smith M*/
501c8c567eSBarry Smith #undef __FUNCT__
511c8c567eSBarry Smith #define __FUNCT__ "TSTrajectoryCreate_Singlefile"
52*972caf09SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Singlefile(TSTrajectory tj,TS ts)
531c8c567eSBarry Smith {
541c8c567eSBarry Smith   PetscErrorCode          ierr;
551c8c567eSBarry Smith   TSTrajectory_Singlefile *sf;
561c8c567eSBarry Smith 
571c8c567eSBarry Smith   PetscFunctionBegin;
581c8c567eSBarry Smith   ierr = PetscNew(&sf);CHKERRQ(ierr);
59*972caf09SHong Zhang   tj->data         = sf;
60*972caf09SHong Zhang   tj->ops->set     = TSTrajectorySet_Singlefile;
61*972caf09SHong Zhang   tj->ops->get     = NULL;
62*972caf09SHong Zhang   tj->ops->destroy = TSTrajectoryDestroy_Singlefile;
631c8c567eSBarry Smith   PetscFunctionReturn(0);
641c8c567eSBarry Smith }
65