Actual source code: trajbasic.c
petsc-3.6.1 2015-08-06
2: #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
6: static PetscErrorCode OutputBIN(const char *filename, PetscViewer *viewer)
7: {
11: PetscViewerCreate(PETSC_COMM_WORLD, viewer);
12: PetscViewerSetType(*viewer, PETSCVIEWERBINARY);
13: PetscViewerFileSetMode(*viewer,FILE_MODE_WRITE);
14: PetscViewerFileSetName(*viewer, filename);
15: return(0);
16: }
21: PetscErrorCode TSTrajectorySet_Basic(TSTrajectory jac,TS ts,PetscInt stepnum,PetscReal time,Vec X)
22: {
23: PetscViewer viewer;
24: PetscInt ns,i;
25: Vec *Y;
26: char filename[PETSC_MAX_PATH_LEN];
27: PetscReal tprev;
31: if (stepnum == 0) {
32: #if defined(PETSC_HAVE_POPEN)
33: TSGetTotalSteps(ts,&stepnum);
34: if (stepnum == 0) {
35: PetscMPIInt rank;
36: MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);
37: if (!rank) {
38: char command[PETSC_MAX_PATH_LEN];
39: FILE *fd;
40: int err;
42: PetscMemzero(command,sizeof(command));
43: PetscSNPrintf(command,PETSC_MAX_PATH_LEN,"rm -fr %s","SA-data");
44: PetscPOpen(PETSC_COMM_SELF,NULL,command,"r",&fd);
45: PetscPClose(PETSC_COMM_SELF,fd,&err);
46: PetscSNPrintf(command,PETSC_MAX_PATH_LEN,"mkdir %s","SA-data");
47: PetscPOpen(PETSC_COMM_SELF,NULL,command,"r",&fd);
48: PetscPClose(PETSC_COMM_SELF,fd,&err);
49: }
50: }
51: #endif
52: return(0);
53: }
54: TSGetTotalSteps(ts,&stepnum);
55: PetscSNPrintf(filename,sizeof(filename),"SA-data/SA-%06d.bin",stepnum);
56: OutputBIN(filename,&viewer);
57: VecView(X,viewer);
58: PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL,PETSC_FALSE);
60: TSGetStages(ts,&ns,&Y);
61: for (i=0;i<ns;i++) {
62: VecView(Y[i],viewer);
63: }
65: TSGetPrevTime(ts,&tprev);
66: PetscViewerBinaryWrite(viewer,&tprev,1,PETSC_REAL,PETSC_FALSE);
68: PetscViewerDestroy(&viewer);
69: return(0);
70: }
74: PetscErrorCode TSTrajectoryGet_Basic(TSTrajectory jac,TS ts,PetscInt step,PetscReal *t)
75: {
76: Vec Sol,*Y;
77: PetscInt Nr,i;
78: PetscViewer viewer;
79: PetscReal timepre;
80: char filename[PETSC_MAX_PATH_LEN];
84: TSGetTotalSteps(ts,&step);
85: PetscSNPrintf(filename,sizeof filename,"SA-data/SA-%06d.bin",step);
86: PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
88: TSGetSolution(ts,&Sol);
89: VecLoad(Sol,viewer);
91: PetscViewerBinaryRead(viewer,t,1,NULL,PETSC_REAL);
93: TSGetStages(ts,&Nr,&Y);
94: for (i=0;i<Nr ;i++) {
95: VecLoad(Y[i],viewer);
96: }
98: /* PetscRealLoad(Nr,&Nr,&timepre,viewer); */
99: PetscViewerBinaryRead(viewer,&timepre,1,NULL,PETSC_REAL);
101: PetscViewerDestroy(&viewer);
103: TSSetTimeStep(ts,-(*t)+timepre);
105: return(0);
106: }
108: /*MC
109: TSTRAJECTORYBASIC - Stores each solution of the ODE/ADE in a file
111: Level: intermediate
113: .seealso: TSTrajectoryCreate(), TS, TSTrajectorySetType()
115: M*/
118: PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Basic(TSTrajectory ts)
119: {
121: ts->ops->set = TSTrajectorySet_Basic;
122: ts->ops->get = TSTrajectoryGet_Basic;
123: return(0);
124: }