xref: /petsc/src/ts/trajectory/impls/singlefile/singlefile.c (revision 1cc06b555e92f8ec64db10330b8bbd830e5bc876)
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 
8d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSTrajectorySet_Singlefile(TSTrajectory tj, TS ts, PetscInt stepnum, PetscReal time, Vec X)
9d71ae5a4SJacob Faibussowitsch {
10972caf09SHong Zhang   TSTrajectory_Singlefile *sf = (TSTrajectory_Singlefile *)tj->data;
111c8c567eSBarry Smith   const char              *filename;
121c8c567eSBarry Smith 
1384e977c4SHong Zhang   PetscFunctionBegin;
141c8c567eSBarry Smith   if (stepnum == 0) {
159566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)X), &sf->viewer));
169566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(sf->viewer, PETSCVIEWERBINARY));
179566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetMode(sf->viewer, FILE_MODE_WRITE));
189566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)tj, &filename));
199566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetName(sf->viewer, filename));
201c8c567eSBarry Smith   }
219566063dSJacob Faibussowitsch   PetscCall(VecView(X, sf->viewer));
229566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWrite(sf->viewer, &time, 1, PETSC_REAL));
233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
241c8c567eSBarry Smith }
251c8c567eSBarry Smith 
26d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSTrajectoryDestroy_Singlefile(TSTrajectory tj)
27d71ae5a4SJacob Faibussowitsch {
28972caf09SHong Zhang   TSTrajectory_Singlefile *sf = (TSTrajectory_Singlefile *)tj->data;
291c8c567eSBarry Smith 
301c8c567eSBarry Smith   PetscFunctionBegin;
319566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&sf->viewer));
329566063dSJacob Faibussowitsch   PetscCall(PetscFree(sf));
333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
341c8c567eSBarry Smith }
351c8c567eSBarry Smith 
361c8c567eSBarry Smith /*MC
37bcf0153eSBarry Smith       TSTRAJECTORYSINGLEFILE - Stores all solutions of the ODE/ADE into a single file followed by each timestep.
38bcf0153eSBarry Smith       Does not save the intermediate stages in a multistage method
391c8c567eSBarry Smith 
401c8c567eSBarry Smith   Level: intermediate
411c8c567eSBarry Smith 
42*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSTrajectoryCreate()`, `TS`, `TSTrajectorySetType()`, `TSTrajectoryType`, `TSTrajectory`
431c8c567eSBarry Smith M*/
44d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Singlefile(TSTrajectory tj, TS ts)
45d71ae5a4SJacob Faibussowitsch {
461c8c567eSBarry Smith   TSTrajectory_Singlefile *sf;
471c8c567eSBarry Smith 
481c8c567eSBarry Smith   PetscFunctionBegin;
499566063dSJacob Faibussowitsch   PetscCall(PetscNew(&sf));
50972caf09SHong Zhang   tj->data         = sf;
51972caf09SHong Zhang   tj->ops->set     = TSTrajectorySet_Singlefile;
52972caf09SHong Zhang   tj->ops->get     = NULL;
53972caf09SHong Zhang   tj->ops->destroy = TSTrajectoryDestroy_Singlefile;
541aebe4aeSStefano Zampini   ts->setupcalled  = PETSC_TRUE;
553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
561c8c567eSBarry Smith }
57