xref: /petsc/src/ts/trajectory/impls/basic/trajbasic.c (revision c3df6c969aaac6746af0a3b9590f7f2f4e39d2c2)
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 
1022a86ba0SHong Zhang   PetscFunctionBegin;
11bc952696SBarry Smith   ierr = PetscViewerCreate(PETSC_COMM_WORLD,viewer);CHKERRQ(ierr);
120bce9e18SHong Zhang   ierr = PetscViewerCreate(PETSC_COMM_WORLD,viewer);CHKERRQ(ierr);
13bc952696SBarry Smith   ierr = PetscViewerSetType(*viewer,PETSCVIEWERBINARY);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"
2052fbb9a2SHong Zhang PetscErrorCode TSTrajectorySet_Basic(TSTrajectory tj,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 
2922a86ba0SHong Zhang   PetscFunctionBegin;
30a8472d1cSHong Zhang   ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr);
31bc952696SBarry Smith   if (stepnum == 0) {
32bc952696SBarry Smith     PetscMPIInt rank;
33bc952696SBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);CHKERRQ(ierr);
34bc952696SBarry Smith     if (!rank) {
359085eb00SSatish Balay       ierr = PetscRMTree("SA-data");CHKERRQ(ierr);
369085eb00SSatish Balay       ierr = PetscMkdir("SA-data");CHKERRQ(ierr);
37bc952696SBarry Smith     }
3826656371SHong Zhang     ierr = PetscSNPrintf(filename,sizeof(filename),"SA-data/SA-%06d.bin",stepnum);CHKERRQ(ierr);
3926656371SHong Zhang     ierr = OutputBIN(filename,&viewer);CHKERRQ(ierr);
4026656371SHong Zhang     ierr = VecView(X,viewer);CHKERRQ(ierr);
4126656371SHong Zhang     ierr = PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
42c186a807SHong Zhang     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
43bc952696SBarry Smith     PetscFunctionReturn(0);
44bc952696SBarry Smith   }
45bc952696SBarry Smith   ierr = PetscSNPrintf(filename,sizeof(filename),"SA-data/SA-%06d.bin",stepnum);CHKERRQ(ierr);
46bc952696SBarry Smith   ierr = OutputBIN(filename,&viewer);CHKERRQ(ierr);
47bc952696SBarry Smith   ierr = VecView(X,viewer);CHKERRQ(ierr);
48c679fc15SHong Zhang   ierr = PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
49bc952696SBarry Smith 
50c679fc15SHong Zhang   ierr = TSGetStages(ts,&ns,&Y);CHKERRQ(ierr);
51bc952696SBarry Smith   for (i=0;i<ns;i++) {
52bc952696SBarry Smith     ierr = VecView(Y[i],viewer);CHKERRQ(ierr);
53bc952696SBarry Smith   }
54bc952696SBarry Smith 
55c679fc15SHong Zhang   ierr = TSGetPrevTime(ts,&tprev);CHKERRQ(ierr);
56c679fc15SHong Zhang   ierr = PetscViewerBinaryWrite(viewer,&tprev,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
57c679fc15SHong Zhang 
58bc952696SBarry Smith   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
59bc952696SBarry Smith   PetscFunctionReturn(0);
60bc952696SBarry Smith }
61bc952696SBarry Smith 
62bc952696SBarry Smith #undef __FUNCT__
63bc952696SBarry Smith #define __FUNCT__ "TSTrajectoryGet_Basic"
6452fbb9a2SHong Zhang PetscErrorCode TSTrajectoryGet_Basic(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal *t)
65bc952696SBarry Smith {
66bc952696SBarry Smith   Vec            Sol,*Y;
67060da220SMatthew G. Knepley   PetscInt       Nr,i;
68bc952696SBarry Smith   PetscViewer    viewer;
69bc952696SBarry Smith   PetscReal      timepre;
70bc952696SBarry Smith   char           filename[PETSC_MAX_PATH_LEN];
71bc952696SBarry Smith   PetscErrorCode ierr;
726cd0a353SHong Zhang 
7322a86ba0SHong Zhang   PetscFunctionBegin;
74a8472d1cSHong Zhang   ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr);
7526656371SHong Zhang   ierr = PetscSNPrintf(filename,sizeof filename,"SA-data/SA-%06d.bin",stepnum);CHKERRQ(ierr);
76bc952696SBarry Smith   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
77bc952696SBarry Smith 
78bc952696SBarry Smith   ierr = TSGetSolution(ts,&Sol);CHKERRQ(ierr);
79bc952696SBarry Smith   ierr = VecLoad(Sol,viewer);CHKERRQ(ierr);
80bc952696SBarry Smith 
81c679fc15SHong Zhang   ierr = PetscViewerBinaryRead(viewer,t,1,NULL,PETSC_REAL);CHKERRQ(ierr);
82bc952696SBarry Smith 
8326656371SHong Zhang   if (stepnum != 0) {
84bc952696SBarry Smith     ierr = TSGetStages(ts,&Nr,&Y);CHKERRQ(ierr);
85bc952696SBarry Smith     for (i=0;i<Nr ;i++) {
86bc952696SBarry Smith       ierr = VecLoad(Y[i],viewer);CHKERRQ(ierr);
87bc952696SBarry Smith     }
88c679fc15SHong Zhang     ierr = PetscViewerBinaryRead(viewer,&timepre,1,NULL,PETSC_REAL);CHKERRQ(ierr);
89acaebbb6SHong Zhang     ierr = TSSetTimeStep(ts,-(*t)+timepre);CHKERRQ(ierr);
9026656371SHong Zhang   }
91bc952696SBarry Smith 
9226656371SHong Zhang   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
93bc952696SBarry Smith   PetscFunctionReturn(0);
94bc952696SBarry Smith }
95bc952696SBarry Smith 
96bc952696SBarry Smith /*MC
971c8c567eSBarry Smith       TSTRAJECTORYBASIC - Stores each solution of the ODE/ADE in a file
98bc952696SBarry Smith 
99bc952696SBarry Smith   Level: intermediate
100bc952696SBarry Smith 
101bc952696SBarry Smith .seealso:  TSTrajectoryCreate(), TS, TSTrajectorySetType()
102bc952696SBarry Smith 
103bc952696SBarry Smith M*/
104bc952696SBarry Smith #undef __FUNCT__
105bc952696SBarry Smith #define __FUNCT__ "TSTrajectoryCreate_Basic"
106*c3df6c96SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Basic(TSTrajectory tj,TS ts)
107bc952696SBarry Smith {
108bc952696SBarry Smith   PetscFunctionBegin;
10952fbb9a2SHong Zhang   tj->ops->set  = TSTrajectorySet_Basic;
11052fbb9a2SHong Zhang   tj->ops->get  = TSTrajectoryGet_Basic;
111bc952696SBarry Smith   PetscFunctionReturn(0);
112bc952696SBarry Smith }
113