xref: /petsc/src/ts/trajectory/impls/singlefile/singlefile.c (revision 1aebe4ae5f460046206914a7c5a819b7cd6afa42)
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 
8560360afSLisandro Dalcin static PetscErrorCode TSTrajectorySet_Singlefile(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X)
91c8c567eSBarry Smith {
10972caf09SHong Zhang   TSTrajectory_Singlefile *sf = (TSTrajectory_Singlefile*)tj->data;
111c8c567eSBarry Smith   PetscErrorCode          ierr;
121c8c567eSBarry Smith   const char              *filename;
131c8c567eSBarry Smith 
1484e977c4SHong Zhang   PetscFunctionBegin;
151c8c567eSBarry Smith   if (stepnum == 0) {
161c8c567eSBarry Smith     ierr = PetscViewerCreate(PETSC_COMM_WORLD,&sf->viewer);CHKERRQ(ierr);
171c8c567eSBarry Smith     ierr = PetscViewerSetType(sf->viewer,PETSCVIEWERBINARY);CHKERRQ(ierr);
181c8c567eSBarry Smith     ierr = PetscViewerFileSetMode(sf->viewer,FILE_MODE_WRITE);CHKERRQ(ierr);
19972caf09SHong Zhang     ierr = PetscObjectGetName((PetscObject)tj,&filename);CHKERRQ(ierr);
201c8c567eSBarry Smith     ierr = PetscViewerFileSetName(sf->viewer,filename);CHKERRQ(ierr);
211c8c567eSBarry Smith   }
22c679fc15SHong Zhang   ierr = VecView(X,sf->viewer);CHKERRQ(ierr);
23c679fc15SHong Zhang   ierr = PetscViewerBinaryWrite(sf->viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
241c8c567eSBarry Smith   PetscFunctionReturn(0);
251c8c567eSBarry Smith }
261c8c567eSBarry Smith 
27560360afSLisandro Dalcin static PetscErrorCode TSTrajectoryDestroy_Singlefile(TSTrajectory tj)
281c8c567eSBarry Smith {
291c8c567eSBarry Smith   PetscErrorCode          ierr;
30972caf09SHong Zhang   TSTrajectory_Singlefile *sf = (TSTrajectory_Singlefile*)tj->data;
311c8c567eSBarry Smith 
321c8c567eSBarry Smith   PetscFunctionBegin;
331c8c567eSBarry Smith   ierr = PetscViewerDestroy(&sf->viewer);CHKERRQ(ierr);
341c8c567eSBarry Smith   ierr = PetscFree(sf);CHKERRQ(ierr);
351c8c567eSBarry Smith   PetscFunctionReturn(0);
361c8c567eSBarry Smith }
371c8c567eSBarry Smith 
381c8c567eSBarry Smith /*MC
398868bf69SBarry 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
401c8c567eSBarry Smith 
411c8c567eSBarry Smith   Level: intermediate
421c8c567eSBarry Smith 
431c8c567eSBarry Smith .seealso:  TSTrajectoryCreate(), TS, TSTrajectorySetType()
441c8c567eSBarry Smith 
451c8c567eSBarry Smith M*/
46972caf09SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Singlefile(TSTrajectory tj,TS ts)
471c8c567eSBarry Smith {
481c8c567eSBarry Smith   PetscErrorCode          ierr;
491c8c567eSBarry Smith   TSTrajectory_Singlefile *sf;
501c8c567eSBarry Smith 
511c8c567eSBarry Smith   PetscFunctionBegin;
521c8c567eSBarry Smith   ierr = PetscNew(&sf);CHKERRQ(ierr);
53972caf09SHong Zhang   tj->data         = sf;
54972caf09SHong Zhang   tj->ops->set     = TSTrajectorySet_Singlefile;
55972caf09SHong Zhang   tj->ops->get     = NULL;
56972caf09SHong Zhang   tj->ops->destroy = TSTrajectoryDestroy_Singlefile;
57*1aebe4aeSStefano Zampini   ts->setupcalled  = PETSC_TRUE;
581c8c567eSBarry Smith   PetscFunctionReturn(0);
591c8c567eSBarry Smith }
60