xref: /petsc/src/ts/trajectory/impls/basic/trajbasic.c (revision 57cc9fe29874d1d66254d0cd1576ebc54e08bc36)
1 
2 #include <petsc/private/tsimpl.h>        /*I "petscts.h"  I*/
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "OutputBIN"
6 static PetscErrorCode OutputBIN(const char *filename, PetscViewer *viewer)
7 {
8   PetscErrorCode ierr;
9 
10   PetscFunctionBeginUser;
11   ierr = PetscViewerCreate(PETSC_COMM_WORLD, viewer);CHKERRQ(ierr);
12   ierr = PetscViewerSetType(*viewer, PETSCVIEWERBINARY);CHKERRQ(ierr);
13   ierr = PetscViewerFileSetMode(*viewer,FILE_MODE_WRITE);CHKERRQ(ierr);
14   ierr = PetscViewerFileSetName(*viewer, filename);CHKERRQ(ierr);
15   PetscFunctionReturn(0);
16 }
17 
18 
19 #undef __FUNCT__
20 #define __FUNCT__ "TSTrajectorySet_Basic"
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;
28   PetscErrorCode ierr;
29   PetscInt       test=stepnum;
30   PetscFunctionBeginUser;
31   if (stepnum == 0) {
32 #if defined(PETSC_HAVE_POPEN)
33     ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr);
34     if (test!=stepnum) printf("test=%d,stepnum=%d\n",test,stepnum);
35     if (stepnum == 0) {
36       PetscMPIInt rank;
37       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);CHKERRQ(ierr);
38       if (!rank) {
39         char command[PETSC_MAX_PATH_LEN];
40         FILE *fd;
41         int  err;
42 
43         ierr = PetscMemzero(command,sizeof(command));CHKERRQ(ierr);
44         ierr = PetscSNPrintf(command,PETSC_MAX_PATH_LEN,"rm -fr %s","SA-data");CHKERRQ(ierr);
45         ierr = PetscPOpen(PETSC_COMM_SELF,NULL,command,"r",&fd);CHKERRQ(ierr);
46         ierr = PetscPClose(PETSC_COMM_SELF,fd,&err);CHKERRQ(ierr);
47         ierr = PetscSNPrintf(command,PETSC_MAX_PATH_LEN,"mkdir %s","SA-data");CHKERRQ(ierr);
48         ierr = PetscPOpen(PETSC_COMM_SELF,NULL,command,"r",&fd);CHKERRQ(ierr);
49         ierr = PetscPClose(PETSC_COMM_SELF,fd,&err);CHKERRQ(ierr);
50       }
51     }
52 #endif
53     ierr = PetscSNPrintf(filename,sizeof(filename),"SA-data/SA-%06d.bin",stepnum);CHKERRQ(ierr);
54     ierr = OutputBIN(filename,&viewer);CHKERRQ(ierr);
55     ierr = VecView(X,viewer);CHKERRQ(ierr);
56     ierr = PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
57     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
58     PetscFunctionReturn(0);
59   }
60   ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr);
61   if (test!=stepnum) printf("test=%d,stepnum=%d\n",test,stepnum);
62   ierr = PetscSNPrintf(filename,sizeof(filename),"SA-data/SA-%06d.bin",stepnum);CHKERRQ(ierr);
63   ierr = OutputBIN(filename,&viewer);CHKERRQ(ierr);
64   ierr = VecView(X,viewer);CHKERRQ(ierr);
65   ierr = PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
66 
67   ierr = TSGetStages(ts,&ns,&Y);CHKERRQ(ierr);
68   for (i=0;i<ns;i++) {
69     ierr = VecView(Y[i],viewer);CHKERRQ(ierr);
70   }
71 
72   ierr = TSGetPrevTime(ts,&tprev);CHKERRQ(ierr);
73   ierr = PetscViewerBinaryWrite(viewer,&tprev,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
74 
75   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
76   PetscFunctionReturn(0);
77 }
78 
79 #undef __FUNCT__
80 #define __FUNCT__ "TSTrajectoryGet_Basic"
81 PetscErrorCode TSTrajectoryGet_Basic(TSTrajectory jac,TS ts,PetscInt stepnum,PetscReal *t)
82 {
83   Vec            Sol,*Y;
84   PetscInt       Nr,i;
85   PetscViewer    viewer;
86   PetscReal      timepre;
87   char           filename[PETSC_MAX_PATH_LEN];
88   PetscErrorCode ierr;
89   PetscInt       test=stepnum;
90   PetscFunctionBeginUser;
91   ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr);
92   if (test!=stepnum) printf("test=%d,stepnum=%d\n",test,stepnum);
93   ierr = PetscSNPrintf(filename,sizeof filename,"SA-data/SA-%06d.bin",stepnum);CHKERRQ(ierr);
94   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
95 
96   ierr = TSGetSolution(ts,&Sol);CHKERRQ(ierr);
97   ierr = VecLoad(Sol,viewer);CHKERRQ(ierr);
98 
99   ierr = PetscViewerBinaryRead(viewer,t,1,NULL,PETSC_REAL);CHKERRQ(ierr);
100 
101   if (stepnum != 0) {
102     ierr = TSGetStages(ts,&Nr,&Y);CHKERRQ(ierr);
103     for (i=0;i<Nr ;i++) {
104       ierr = VecLoad(Y[i],viewer);CHKERRQ(ierr);
105     }
106     ierr = PetscViewerBinaryRead(viewer,&timepre,1,NULL,PETSC_REAL);CHKERRQ(ierr);
107     ierr = TSSetTimeStep(ts,-(*t)+timepre);CHKERRQ(ierr);
108   }
109 
110   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
111   PetscFunctionReturn(0);
112 }
113 
114 /*MC
115       TSTRAJECTORYBASIC - Stores each solution of the ODE/ADE in a file
116 
117   Level: intermediate
118 
119 .seealso:  TSTrajectoryCreate(), TS, TSTrajectorySetType()
120 
121 M*/
122 #undef __FUNCT__
123 #define __FUNCT__ "TSTrajectoryCreate_Basic"
124 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Basic(TSTrajectory ts)
125 {
126   PetscFunctionBegin;
127   ts->ops->set  = TSTrajectorySet_Basic;
128   ts->ops->get  = TSTrajectoryGet_Basic;
129   PetscFunctionReturn(0);
130 }
131