xref: /petsc/src/ts/trajectory/impls/visualization/trajvisualization.c (revision 2b04316782a66991320b63a468a551c6b8c983cd)
1*2b043167SHong Zhang 
2*2b043167SHong Zhang #include <petsc/private/tsimpl.h>        /*I "petscts.h"  I*/
3*2b043167SHong Zhang 
4*2b043167SHong Zhang #undef __FUNCT__
5*2b043167SHong Zhang #define __FUNCT__ "OutputBIN"
6*2b043167SHong Zhang static PetscErrorCode OutputBIN(const char *filename,PetscViewer *viewer)
7*2b043167SHong Zhang {
8*2b043167SHong Zhang   PetscErrorCode ierr;
9*2b043167SHong Zhang 
10*2b043167SHong Zhang   PetscFunctionBegin;
11*2b043167SHong Zhang   ierr = PetscViewerCreate(PETSC_COMM_WORLD,viewer);CHKERRQ(ierr);
12*2b043167SHong Zhang   ierr = PetscViewerSetType(*viewer,PETSCVIEWERBINARY);CHKERRQ(ierr);
13*2b043167SHong Zhang   ierr = PetscViewerFileSetMode(*viewer,FILE_MODE_WRITE);CHKERRQ(ierr);
14*2b043167SHong Zhang   ierr = PetscViewerFileSetName(*viewer,filename);CHKERRQ(ierr);
15*2b043167SHong Zhang   PetscFunctionReturn(0);
16*2b043167SHong Zhang }
17*2b043167SHong Zhang 
18*2b043167SHong Zhang #undef __FUNCT__
19*2b043167SHong Zhang #define __FUNCT__ "TSTrajectorySet_Visualization"
20*2b043167SHong Zhang PetscErrorCode TSTrajectorySet_Visualization(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X)
21*2b043167SHong Zhang {
22*2b043167SHong Zhang   PetscViewer    viewer;
23*2b043167SHong Zhang   char           filename[PETSC_MAX_PATH_LEN];
24*2b043167SHong Zhang   PetscReal      tprev;
25*2b043167SHong Zhang   PetscErrorCode ierr;
26*2b043167SHong Zhang 
27*2b043167SHong Zhang   PetscFunctionBegin;
28*2b043167SHong Zhang   ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr);
29*2b043167SHong Zhang   if (stepnum == 0) {
30*2b043167SHong Zhang     PetscMPIInt rank;
31*2b043167SHong Zhang     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);CHKERRQ(ierr);
32*2b043167SHong Zhang     if (!rank) {
33*2b043167SHong Zhang       ierr = PetscRMTree("Visualization-data");CHKERRQ(ierr);
34*2b043167SHong Zhang       ierr = PetscMkdir("Visualization-data");CHKERRQ(ierr);
35*2b043167SHong Zhang     }
36*2b043167SHong Zhang     ierr = PetscSNPrintf(filename,sizeof(filename),"Visualization-data/SA-%06d.bin",stepnum);CHKERRQ(ierr);
37*2b043167SHong Zhang     ierr = OutputBIN(filename,&viewer);CHKERRQ(ierr);
38*2b043167SHong Zhang     ierr = VecView(X,viewer);CHKERRQ(ierr);
39*2b043167SHong Zhang     ierr = PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
40*2b043167SHong Zhang     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
41*2b043167SHong Zhang     PetscFunctionReturn(0);
42*2b043167SHong Zhang   }
43*2b043167SHong Zhang   ierr = PetscSNPrintf(filename,sizeof(filename),"Visualization-data/SA-%06d.bin",stepnum);CHKERRQ(ierr);
44*2b043167SHong Zhang   ierr = OutputBIN(filename,&viewer);CHKERRQ(ierr);
45*2b043167SHong Zhang   ierr = VecView(X,viewer);CHKERRQ(ierr);
46*2b043167SHong Zhang   ierr = PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
47*2b043167SHong Zhang 
48*2b043167SHong Zhang   ierr = TSGetPrevTime(ts,&tprev);CHKERRQ(ierr);
49*2b043167SHong Zhang   ierr = PetscViewerBinaryWrite(viewer,&tprev,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
50*2b043167SHong Zhang 
51*2b043167SHong Zhang   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
52*2b043167SHong Zhang   PetscFunctionReturn(0);
53*2b043167SHong Zhang }
54*2b043167SHong Zhang 
55*2b043167SHong Zhang /*MC
56*2b043167SHong Zhang       TSTRAJECTORYBASIC - Stores each solution of the ODE/ADE in a file
57*2b043167SHong Zhang 
58*2b043167SHong Zhang   Level: intermediate
59*2b043167SHong Zhang 
60*2b043167SHong Zhang .seealso:  TSTrajectoryCreate(), TS, TSTrajectorySetType()
61*2b043167SHong Zhang 
62*2b043167SHong Zhang M*/
63*2b043167SHong Zhang #undef __FUNCT__
64*2b043167SHong Zhang #define __FUNCT__ "TSTrajectoryCreate_Visualization"
65*2b043167SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Visualization(TSTrajectory tj,TS ts)
66*2b043167SHong Zhang {
67*2b043167SHong Zhang   PetscFunctionBegin;
68*2b043167SHong Zhang   tj->ops->set  = TSTrajectorySet_Visualization;
69*2b043167SHong Zhang   PetscFunctionReturn(0);
70*2b043167SHong Zhang }
71