xref: /petsc/src/ts/trajectory/impls/basic/trajbasic.c (revision acaebbb693a27a7efb4b151c67e018b05ed9fbaa) !
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 
10bc952696SBarry Smith   PetscFunctionBeginUser;
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 
19bc952696SBarry Smith #undef __FUNCT__
20bc952696SBarry Smith #define __FUNCT__ "TSTrajectorySet_Basic"
21bc952696SBarry Smith PetscErrorCode TSTrajectorySet_Basic(TSTrajectory jac,TS ts,PetscInt stepnum,PetscReal time,Vec X)
22bc952696SBarry Smith {
23bc952696SBarry Smith   PetscViewer    viewer;
24bc952696SBarry Smith   PetscInt       ns,i;
25bc952696SBarry Smith   Vec            *Y;
26bc952696SBarry Smith   char           filename[PETSC_MAX_PATH_LEN];
27bc952696SBarry Smith   PetscReal      tprev;
28bc952696SBarry Smith   PetscErrorCode ierr;
29c186a807SHong Zhang   PetscInt       test=stepnum;
30bc952696SBarry Smith   PetscFunctionBeginUser;
31bc952696SBarry Smith   if (stepnum == 0) {
32bc952696SBarry Smith #if defined(PETSC_HAVE_POPEN)
33bc952696SBarry Smith     ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr);
34c186a807SHong Zhang     if (test!=stepnum) printf("test=%d,stepnum=%d\n",test,stepnum);
35bc952696SBarry Smith     if (stepnum == 0) {
36bc952696SBarry Smith       PetscMPIInt rank;
37bc952696SBarry Smith       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);CHKERRQ(ierr);
38bc952696SBarry Smith       if (!rank) {
39bc952696SBarry Smith         char command[PETSC_MAX_PATH_LEN];
40bc952696SBarry Smith         FILE *fd;
41bc952696SBarry Smith         int  err;
42bc952696SBarry Smith 
43bc952696SBarry Smith         ierr = PetscMemzero(command,sizeof(command));CHKERRQ(ierr);
44bc952696SBarry Smith         ierr = PetscSNPrintf(command,PETSC_MAX_PATH_LEN,"rm -fr %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         ierr = PetscSNPrintf(command,PETSC_MAX_PATH_LEN,"mkdir %s","SA-data");CHKERRQ(ierr);
48bc952696SBarry Smith         ierr = PetscPOpen(PETSC_COMM_SELF,NULL,command,"r",&fd);CHKERRQ(ierr);
49bc952696SBarry Smith         ierr = PetscPClose(PETSC_COMM_SELF,fd,&err);CHKERRQ(ierr);
50bc952696SBarry Smith       }
51bc952696SBarry Smith     }
52bc952696SBarry Smith #endif
5326656371SHong Zhang     ierr = PetscSNPrintf(filename,sizeof(filename),"SA-data/SA-%06d.bin",stepnum);CHKERRQ(ierr);
5426656371SHong Zhang     ierr = OutputBIN(filename,&viewer);CHKERRQ(ierr);
5526656371SHong Zhang     ierr = VecView(X,viewer);CHKERRQ(ierr);
5626656371SHong Zhang     ierr = PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
57c186a807SHong Zhang     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
58bc952696SBarry Smith     PetscFunctionReturn(0);
59bc952696SBarry Smith   }
60bc952696SBarry Smith   ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr);
61c186a807SHong Zhang   if (test!=stepnum) printf("test=%d,stepnum=%d\n",test,stepnum);
62bc952696SBarry Smith   ierr = PetscSNPrintf(filename,sizeof(filename),"SA-data/SA-%06d.bin",stepnum);CHKERRQ(ierr);
63bc952696SBarry Smith   ierr = OutputBIN(filename,&viewer);CHKERRQ(ierr);
64bc952696SBarry Smith   ierr = VecView(X,viewer);CHKERRQ(ierr);
65c679fc15SHong Zhang   ierr = PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
66bc952696SBarry Smith 
67c679fc15SHong Zhang   ierr = TSGetStages(ts,&ns,&Y);CHKERRQ(ierr);
68bc952696SBarry Smith   for (i=0;i<ns;i++) {
69bc952696SBarry Smith     ierr = VecView(Y[i],viewer);CHKERRQ(ierr);
70bc952696SBarry Smith   }
71bc952696SBarry Smith 
72c679fc15SHong Zhang   ierr = TSGetPrevTime(ts,&tprev);CHKERRQ(ierr);
73c679fc15SHong Zhang   ierr = PetscViewerBinaryWrite(viewer,&tprev,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
74c679fc15SHong Zhang 
75bc952696SBarry Smith   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
76bc952696SBarry Smith   PetscFunctionReturn(0);
77bc952696SBarry Smith }
78bc952696SBarry Smith 
79bc952696SBarry Smith #undef __FUNCT__
80bc952696SBarry Smith #define __FUNCT__ "TSTrajectoryGet_Basic"
8126656371SHong Zhang PetscErrorCode TSTrajectoryGet_Basic(TSTrajectory jac,TS ts,PetscInt stepnum,PetscReal *t)
82bc952696SBarry Smith {
83bc952696SBarry Smith   Vec            Sol,*Y;
84060da220SMatthew G. Knepley   PetscInt       Nr,i;
85bc952696SBarry Smith   PetscViewer    viewer;
86bc952696SBarry Smith   PetscReal      timepre;
87bc952696SBarry Smith   char           filename[PETSC_MAX_PATH_LEN];
88bc952696SBarry Smith   PetscErrorCode ierr;
89c186a807SHong Zhang   PetscInt       test=stepnum;
90bc952696SBarry Smith   PetscFunctionBeginUser;
9126656371SHong Zhang   ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr);
92c186a807SHong Zhang   if (test!=stepnum) printf("test=%d,stepnum=%d\n",test,stepnum);
9326656371SHong Zhang   ierr = PetscSNPrintf(filename,sizeof filename,"SA-data/SA-%06d.bin",stepnum);CHKERRQ(ierr);
94bc952696SBarry Smith   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
95bc952696SBarry Smith 
96bc952696SBarry Smith   ierr = TSGetSolution(ts,&Sol);CHKERRQ(ierr);
97bc952696SBarry Smith   ierr = VecLoad(Sol,viewer);CHKERRQ(ierr);
98bc952696SBarry Smith 
99c679fc15SHong Zhang   ierr = PetscViewerBinaryRead(viewer,t,1,NULL,PETSC_REAL);CHKERRQ(ierr);
100bc952696SBarry Smith 
10126656371SHong Zhang   if (stepnum != 0) {
102bc952696SBarry Smith     ierr = TSGetStages(ts,&Nr,&Y);CHKERRQ(ierr);
103bc952696SBarry Smith     for (i=0;i<Nr ;i++) {
104bc952696SBarry Smith       ierr = VecLoad(Y[i],viewer);CHKERRQ(ierr);
105bc952696SBarry Smith     }
106c679fc15SHong Zhang     ierr = PetscViewerBinaryRead(viewer,&timepre,1,NULL,PETSC_REAL);CHKERRQ(ierr);
107*acaebbb6SHong Zhang     ierr = TSSetTimeStep(ts,-(*t)+timepre);CHKERRQ(ierr);
10826656371SHong Zhang   }
109bc952696SBarry Smith 
11026656371SHong Zhang   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
111bc952696SBarry Smith   PetscFunctionReturn(0);
112bc952696SBarry Smith }
113bc952696SBarry Smith 
114bc952696SBarry Smith /*MC
1151c8c567eSBarry Smith       TSTRAJECTORYBASIC - Stores each solution of the ODE/ADE in a file
116bc952696SBarry Smith 
117bc952696SBarry Smith   Level: intermediate
118bc952696SBarry Smith 
119bc952696SBarry Smith .seealso:  TSTrajectoryCreate(), TS, TSTrajectorySetType()
120bc952696SBarry Smith 
121bc952696SBarry Smith M*/
122bc952696SBarry Smith #undef __FUNCT__
123bc952696SBarry Smith #define __FUNCT__ "TSTrajectoryCreate_Basic"
124bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Basic(TSTrajectory ts)
125bc952696SBarry Smith {
126bc952696SBarry Smith   PetscFunctionBegin;
127bc952696SBarry Smith   ts->ops->set  = TSTrajectorySet_Basic;
128bc952696SBarry Smith   ts->ops->get  = TSTrajectoryGet_Basic;
129bc952696SBarry Smith   PetscFunctionReturn(0);
130bc952696SBarry Smith }
131