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