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