Actual source code: singlefile.c
1: #include <petsc/private/tsimpl.h>
3: typedef struct {
4: PetscViewer viewer;
5: } TSTrajectory_Singlefile;
7: static PetscErrorCode TSTrajectorySet_Singlefile(TSTrajectory tj, TS ts, PetscInt stepnum, PetscReal time, Vec X)
8: {
9: TSTrajectory_Singlefile *sf = (TSTrajectory_Singlefile *)tj->data;
10: const char *filename;
12: PetscFunctionBegin;
13: if (stepnum == 0) {
14: PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)X), &sf->viewer));
15: PetscCall(PetscViewerSetType(sf->viewer, PETSCVIEWERBINARY));
16: PetscCall(PetscViewerFileSetMode(sf->viewer, FILE_MODE_WRITE));
17: PetscCall(PetscObjectGetName((PetscObject)tj, &filename));
18: PetscCall(PetscViewerFileSetName(sf->viewer, filename));
19: }
20: PetscCall(VecView(X, sf->viewer));
21: PetscCall(PetscViewerBinaryWrite(sf->viewer, &time, 1, PETSC_REAL));
22: PetscFunctionReturn(PETSC_SUCCESS);
23: }
25: static PetscErrorCode TSTrajectoryDestroy_Singlefile(TSTrajectory tj)
26: {
27: TSTrajectory_Singlefile *sf = (TSTrajectory_Singlefile *)tj->data;
29: PetscFunctionBegin;
30: PetscCall(PetscViewerDestroy(&sf->viewer));
31: PetscCall(PetscFree(sf));
32: PetscFunctionReturn(PETSC_SUCCESS);
33: }
35: /*MC
36: TSTRAJECTORYSINGLEFILE - Stores all solutions of the ODE/ADE into a single file followed by each timestep.
37: Does not save the intermediate stages in a multistage method
39: Level: intermediate
41: .seealso: [](ch_ts), `TSTrajectoryCreate()`, `TS`, `TSTrajectorySetType()`, `TSTrajectoryType`, `TSTrajectory`
42: M*/
43: PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Singlefile(TSTrajectory tj, TS ts)
44: {
45: TSTrajectory_Singlefile *sf;
47: PetscFunctionBegin;
48: PetscCall(PetscNew(&sf));
49: tj->data = sf;
50: tj->ops->set = TSTrajectorySet_Singlefile;
51: tj->ops->get = NULL;
52: tj->ops->destroy = TSTrajectoryDestroy_Singlefile;
53: ts->setupcalled = PETSC_TRUE;
54: PetscFunctionReturn(PETSC_SUCCESS);
55: }