xref: /petsc/src/ts/trajectory/impls/singlefile/singlefile.c (revision bef158480efac06de457f7a665168877ab3c2fd7)
1 
2 #include <petsc/private/tsimpl.h>        /*I "petscts.h"  I*/
3 
4 typedef struct {
5   PetscViewer viewer;
6 } TSTrajectory_Singlefile;
7 
8 static PetscErrorCode TSTrajectorySet_Singlefile(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X)
9 {
10   TSTrajectory_Singlefile *sf = (TSTrajectory_Singlefile*)tj->data;
11   PetscErrorCode          ierr;
12   const char              *filename;
13 
14   PetscFunctionBegin;
15   if (stepnum == 0) {
16     ierr = PetscViewerCreate(PetscObjectComm((PetscObject)X),&sf->viewer);CHKERRQ(ierr);
17     ierr = PetscViewerSetType(sf->viewer,PETSCVIEWERBINARY);CHKERRQ(ierr);
18     ierr = PetscViewerFileSetMode(sf->viewer,FILE_MODE_WRITE);CHKERRQ(ierr);
19     ierr = PetscObjectGetName((PetscObject)tj,&filename);CHKERRQ(ierr);
20     ierr = PetscViewerFileSetName(sf->viewer,filename);CHKERRQ(ierr);
21   }
22   ierr = VecView(X,sf->viewer);CHKERRQ(ierr);
23   ierr = PetscViewerBinaryWrite(sf->viewer,&time,1,PETSC_REAL);CHKERRQ(ierr);
24   PetscFunctionReturn(0);
25 }
26 
27 static PetscErrorCode TSTrajectoryDestroy_Singlefile(TSTrajectory tj)
28 {
29   PetscErrorCode          ierr;
30   TSTrajectory_Singlefile *sf = (TSTrajectory_Singlefile*)tj->data;
31 
32   PetscFunctionBegin;
33   ierr = PetscViewerDestroy(&sf->viewer);CHKERRQ(ierr);
34   ierr = PetscFree(sf);CHKERRQ(ierr);
35   PetscFunctionReturn(0);
36 }
37 
38 /*MC
39       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
40 
41   Level: intermediate
42 
43 .seealso:  TSTrajectoryCreate(), TS, TSTrajectorySetType()
44 
45 M*/
46 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Singlefile(TSTrajectory tj,TS ts)
47 {
48   PetscErrorCode          ierr;
49   TSTrajectory_Singlefile *sf;
50 
51   PetscFunctionBegin;
52   ierr = PetscNew(&sf);CHKERRQ(ierr);
53   tj->data         = sf;
54   tj->ops->set     = TSTrajectorySet_Singlefile;
55   tj->ops->get     = NULL;
56   tj->ops->destroy = TSTrajectoryDestroy_Singlefile;
57   ts->setupcalled  = PETSC_TRUE;
58   PetscFunctionReturn(0);
59 }
60