xref: /petsc/src/ts/trajectory/impls/singlefile/singlefile.c (revision e8e8640d1cb9a3a2f50c0c0d7b26e5c4d521e2e4)
1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h"  I*/
21c8c567eSBarry Smith 
31c8c567eSBarry Smith typedef struct {
41c8c567eSBarry Smith   PetscViewer viewer;
51c8c567eSBarry Smith } TSTrajectory_Singlefile;
61c8c567eSBarry Smith 
TSTrajectorySet_Singlefile(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X)7d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSTrajectorySet_Singlefile(TSTrajectory tj, TS ts, PetscInt stepnum, PetscReal time, Vec X)
8d71ae5a4SJacob Faibussowitsch {
9972caf09SHong Zhang   TSTrajectory_Singlefile *sf = (TSTrajectory_Singlefile *)tj->data;
101c8c567eSBarry Smith   const char              *filename;
111c8c567eSBarry Smith 
1284e977c4SHong Zhang   PetscFunctionBegin;
131c8c567eSBarry Smith   if (stepnum == 0) {
149566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)X), &sf->viewer));
159566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(sf->viewer, PETSCVIEWERBINARY));
169566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetMode(sf->viewer, FILE_MODE_WRITE));
179566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)tj, &filename));
189566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetName(sf->viewer, filename));
191c8c567eSBarry Smith   }
209566063dSJacob Faibussowitsch   PetscCall(VecView(X, sf->viewer));
219566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWrite(sf->viewer, &time, 1, PETSC_REAL));
223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
231c8c567eSBarry Smith }
241c8c567eSBarry Smith 
TSTrajectoryDestroy_Singlefile(TSTrajectory tj)25d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSTrajectoryDestroy_Singlefile(TSTrajectory tj)
26d71ae5a4SJacob Faibussowitsch {
27972caf09SHong Zhang   TSTrajectory_Singlefile *sf = (TSTrajectory_Singlefile *)tj->data;
281c8c567eSBarry Smith 
291c8c567eSBarry Smith   PetscFunctionBegin;
309566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&sf->viewer));
319566063dSJacob Faibussowitsch   PetscCall(PetscFree(sf));
323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
331c8c567eSBarry Smith }
341c8c567eSBarry Smith 
351c8c567eSBarry Smith /*MC
36bcf0153eSBarry Smith       TSTRAJECTORYSINGLEFILE - Stores all solutions of the ODE/ADE into a single file followed by each timestep.
37bcf0153eSBarry Smith       Does not save the intermediate stages in a multistage method
381c8c567eSBarry Smith 
391c8c567eSBarry Smith   Level: intermediate
401c8c567eSBarry Smith 
41*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSTrajectoryCreate()`, `TS`, `TSTrajectorySetType()`, `TSTrajectoryType`, `TSTrajectory`
421c8c567eSBarry Smith M*/
TSTrajectoryCreate_Singlefile(TSTrajectory tj,TS ts)43d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Singlefile(TSTrajectory tj, TS ts)
44d71ae5a4SJacob Faibussowitsch {
451c8c567eSBarry Smith   TSTrajectory_Singlefile *sf;
461c8c567eSBarry Smith 
471c8c567eSBarry Smith   PetscFunctionBegin;
489566063dSJacob Faibussowitsch   PetscCall(PetscNew(&sf));
49972caf09SHong Zhang   tj->data         = sf;
50972caf09SHong Zhang   tj->ops->set     = TSTrajectorySet_Singlefile;
51972caf09SHong Zhang   tj->ops->get     = NULL;
52972caf09SHong Zhang   tj->ops->destroy = TSTrajectoryDestroy_Singlefile;
531aebe4aeSStefano Zampini   ts->setupcalled  = PETSC_TRUE;
543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
551c8c567eSBarry Smith }
56