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