Actual source code: singlefile.c
petsc-3.7.3 2016-08-01
2: #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
4: typedef struct {
5: PetscViewer viewer;
6: } TSTrajectory_Singlefile;
10: static PetscErrorCode TSTrajectorySet_Singlefile(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X)
11: {
12: TSTrajectory_Singlefile *sf = (TSTrajectory_Singlefile*)tj->data;
13: PetscErrorCode ierr;
14: const char *filename;
17: if (stepnum == 0) {
18: PetscViewerCreate(PETSC_COMM_WORLD,&sf->viewer);
19: PetscViewerSetType(sf->viewer,PETSCVIEWERBINARY);
20: PetscViewerFileSetMode(sf->viewer,FILE_MODE_WRITE);
21: PetscObjectGetName((PetscObject)tj,&filename);
22: PetscViewerFileSetName(sf->viewer,filename);
23: }
24: VecView(X,sf->viewer);
25: PetscViewerBinaryWrite(sf->viewer,&time,1,PETSC_REAL,PETSC_FALSE);
26: return(0);
27: }
31: static PetscErrorCode TSTrajectoryDestroy_Singlefile(TSTrajectory tj)
32: {
33: PetscErrorCode ierr;
34: TSTrajectory_Singlefile *sf = (TSTrajectory_Singlefile*)tj->data;
37: PetscViewerDestroy(&sf->viewer);
38: PetscFree(sf);
39: return(0);
40: }
42: /*MC
43: 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
45: Level: intermediate
47: .seealso: TSTrajectoryCreate(), TS, TSTrajectorySetType()
49: M*/
52: PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Singlefile(TSTrajectory tj,TS ts)
53: {
54: PetscErrorCode ierr;
55: TSTrajectory_Singlefile *sf;
58: PetscNew(&sf);
59: tj->data = sf;
60: tj->ops->set = TSTrajectorySet_Singlefile;
61: tj->ops->get = NULL;
62: tj->ops->destroy = TSTrajectoryDestroy_Singlefile;
63: return(0);
64: }