Actual source code: trajbasic.c
petsc-3.9.4 2018-09-11
2: #include <petsc/private/tsimpl.h>
4: static PetscErrorCode OutputBIN(MPI_Comm comm,const char *filename,PetscViewer *viewer)
5: {
9: PetscViewerCreate(comm,viewer);
10: PetscViewerSetType(*viewer,PETSCVIEWERBINARY);
11: PetscViewerFileSetMode(*viewer,FILE_MODE_WRITE);
12: PetscViewerFileSetName(*viewer,filename);
13: return(0);
14: }
16: static PetscErrorCode TSTrajectorySet_Basic(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X)
17: {
18: PetscViewer viewer;
19: PetscInt ns,i;
20: Vec *Y;
21: char filename[PETSC_MAX_PATH_LEN];
22: PetscReal tprev;
24: MPI_Comm comm;
27: PetscObjectGetComm((PetscObject)ts,&comm);
28: TSGetStepNumber(ts,&stepnum);
29: if (stepnum == 0) {
30: PetscMPIInt rank;
31: MPI_Comm_rank(comm,&rank);
32: if (!rank) {
33: PetscRMTree(tj->dirname);
34: PetscMkdir(tj->dirname);
35: }
36: PetscSNPrintf(filename,sizeof(filename),tj->dirfiletemplate,stepnum);
37: OutputBIN(comm,filename,&viewer);
38: VecView(X,viewer);
39: PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL,PETSC_FALSE);
40: PetscViewerDestroy(&viewer);
41: return(0);
42: }
43: PetscSNPrintf(filename,sizeof(filename),tj->dirfiletemplate,stepnum);
44: OutputBIN(comm,filename,&viewer);
45: VecView(X,viewer);
46: PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL,PETSC_FALSE);
48: TSGetStages(ts,&ns,&Y);
49: for (i=0;i<ns;i++) {
50: VecView(Y[i],viewer);
51: }
53: TSGetPrevTime(ts,&tprev);
54: PetscViewerBinaryWrite(viewer,&tprev,1,PETSC_REAL,PETSC_FALSE);
56: PetscViewerDestroy(&viewer);
57: return(0);
58: }
60: static PetscErrorCode TSTrajectoryGet_Basic(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal *t)
61: {
62: Vec Sol,*Y;
63: PetscInt Nr,i;
64: PetscViewer viewer;
65: PetscReal timepre;
66: char filename[PETSC_MAX_PATH_LEN];
70: PetscSNPrintf(filename,sizeof filename,tj->dirfiletemplate,stepnum);
71: PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
73: TSGetSolution(ts,&Sol);
74: VecLoad(Sol,viewer);
76: PetscViewerBinaryRead(viewer,t,1,NULL,PETSC_REAL);
78: if (stepnum != 0) {
79: TSGetStages(ts,&Nr,&Y);
80: for (i=0;i<Nr ;i++) {
81: VecLoad(Y[i],viewer);
82: }
83: PetscViewerBinaryRead(viewer,&timepre,1,NULL,PETSC_REAL);
84: TSSetTimeStep(ts,-(*t)+timepre);
85: }
87: PetscViewerDestroy(&viewer);
88: return(0);
89: }
91: static PetscErrorCode TSTrajectoryDestroy_Basic(TSTrajectory tj)
92: {
94: PetscMPIInt rank;
95: MPI_Comm comm;
98: if (!tj->keepfiles) {
99: PetscObjectGetComm((PetscObject)tj,&comm);
100: MPI_Comm_rank(comm,&rank);
101: if (!rank) {
102: PetscRMTree(tj->dirname);
103: }
104: }
105: return(0);
106: }
108: /*MC
109: TSTRAJECTORYBASIC - Stores each solution of the ODE/DAE in a file
111: Saves each timestep into a seperate file named SA-data/SA-%06d.bin. The file name can be changed.
113: This version saves the solutions at all the stages
115: $PETSC_DIR/share/petsc/matlab/PetscReadBinaryTrajectory.m can read in files created with this format
117: Level: intermediate
119: .seealso: TSTrajectoryCreate(), TS, TSTrajectorySetType(), TSTrajectorySetDirname(), TSTrajectorySetFile()
121: M*/
122: PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Basic(TSTrajectory tj,TS ts)
123: {
125: tj->keepfiles = PETSC_FALSE;
126: tj->ops->set = TSTrajectorySet_Basic;
127: tj->ops->get = TSTrajectoryGet_Basic;
128: tj->ops->destroy = TSTrajectoryDestroy_Basic;
129: return(0);
130: }