xref: /petsc/src/ts/trajectory/impls/visualization/trajvisualization.c (revision 78fbdcc8da47d2da69f91a746197bf2bceb22ad2)
12b043167SHong Zhang 
22b043167SHong Zhang #include <petsc/private/tsimpl.h>        /*I "petscts.h"  I*/
32b043167SHong Zhang 
4*78fbdcc8SBarry Smith static PetscErrorCode OutputBIN(MPI_Comm comm,const char *filename,PetscViewer *viewer)
52b043167SHong Zhang {
62b043167SHong Zhang   PetscErrorCode ierr;
72b043167SHong Zhang 
82b043167SHong Zhang   PetscFunctionBegin;
92b043167SHong Zhang   ierr = PetscViewerCreate(PETSC_COMM_WORLD,viewer);CHKERRQ(ierr);
102b043167SHong Zhang   ierr = PetscViewerSetType(*viewer,PETSCVIEWERBINARY);CHKERRQ(ierr);
112b043167SHong Zhang   ierr = PetscViewerFileSetMode(*viewer,FILE_MODE_WRITE);CHKERRQ(ierr);
122b043167SHong Zhang   ierr = PetscViewerFileSetName(*viewer,filename);CHKERRQ(ierr);
132b043167SHong Zhang   PetscFunctionReturn(0);
142b043167SHong Zhang }
152b043167SHong Zhang 
16560360afSLisandro Dalcin static PetscErrorCode TSTrajectorySet_Visualization(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X)
172b043167SHong Zhang {
182b043167SHong Zhang   PetscViewer    viewer;
192b043167SHong Zhang   char           filename[PETSC_MAX_PATH_LEN];
202b043167SHong Zhang   PetscReal      tprev;
212b043167SHong Zhang   PetscErrorCode ierr;
22*78fbdcc8SBarry Smith   MPI_Comm       comm;
232b043167SHong Zhang 
242b043167SHong Zhang   PetscFunctionBegin;
25*78fbdcc8SBarry Smith   ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
262b043167SHong Zhang   ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr);
272b043167SHong Zhang   if (stepnum == 0) {
282b043167SHong Zhang     PetscMPIInt rank;
29*78fbdcc8SBarry Smith     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
302b043167SHong Zhang     if (!rank) {
312b043167SHong Zhang       ierr = PetscRMTree("Visualization-data");CHKERRQ(ierr);
322b043167SHong Zhang       ierr = PetscMkdir("Visualization-data");CHKERRQ(ierr);
332b043167SHong Zhang     }
34*78fbdcc8SBarry Smith     if (tj->names) {
35*78fbdcc8SBarry Smith       PetscViewer bnames;
36*78fbdcc8SBarry Smith       ierr = PetscViewerBinaryOpen(comm,"Visualization-data/variablenames",FILE_MODE_WRITE,&bnames);CHKERRQ(ierr);
37*78fbdcc8SBarry Smith       ierr = PetscViewerBinaryWriteStringArray(bnames,(const char *const *)tj->names);CHKERRQ(ierr);
38*78fbdcc8SBarry Smith       ierr = PetscViewerDestroy(&bnames);CHKERRQ(ierr);
39*78fbdcc8SBarry Smith     }
402b043167SHong Zhang     ierr = PetscSNPrintf(filename,sizeof(filename),"Visualization-data/SA-%06d.bin",stepnum);CHKERRQ(ierr);
41*78fbdcc8SBarry Smith     ierr = OutputBIN(comm,filename,&viewer);CHKERRQ(ierr);
422b043167SHong Zhang     ierr = VecView(X,viewer);CHKERRQ(ierr);
432b043167SHong Zhang     ierr = PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
442b043167SHong Zhang     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
452b043167SHong Zhang     PetscFunctionReturn(0);
462b043167SHong Zhang   }
472b043167SHong Zhang   ierr = PetscSNPrintf(filename,sizeof(filename),"Visualization-data/SA-%06d.bin",stepnum);CHKERRQ(ierr);
48*78fbdcc8SBarry Smith   ierr = OutputBIN(comm,filename,&viewer);CHKERRQ(ierr);
492b043167SHong Zhang   ierr = VecView(X,viewer);CHKERRQ(ierr);
502b043167SHong Zhang   ierr = PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
512b043167SHong Zhang 
522b043167SHong Zhang   ierr = TSGetPrevTime(ts,&tprev);CHKERRQ(ierr);
532b043167SHong Zhang   ierr = PetscViewerBinaryWrite(viewer,&tprev,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
542b043167SHong Zhang 
552b043167SHong Zhang   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
562b043167SHong Zhang   PetscFunctionReturn(0);
572b043167SHong Zhang }
582b043167SHong Zhang 
592b043167SHong Zhang /*MC
605c9db295SHong Zhang       TSTRAJECTORYVISUALIZATION - Stores each solution of the ODE/DAE in a file
612b043167SHong Zhang 
62*78fbdcc8SBarry Smith       Saves each timestep into a seperate file in Visualization-data/SA-%06d.bin
63*78fbdcc8SBarry Smith 
64*78fbdcc8SBarry Smith       This version saves only the solutions at each timestep, it does not save the solution at each stage,
65*78fbdcc8SBarry Smith       see TSTRAJECTORYBASIC that saves all stages
66*78fbdcc8SBarry Smith 
67*78fbdcc8SBarry Smith       $PETSC_DIR/share/petsc/matlab/PetscReadBinaryTrajectory.m and $PETSC_DIR/bin/PetscBinaryIOTrajectory.py
68*78fbdcc8SBarry Smith       can read in files created with this format into MATLAB and Python.
69*78fbdcc8SBarry Smith 
702b043167SHong Zhang   Level: intermediate
712b043167SHong Zhang 
72*78fbdcc8SBarry Smith .seealso:  TSTrajectoryCreate(), TS, TSTrajectorySetType(), TSTrajectoryType, TSTrajectorySetVariableNames()
732b043167SHong Zhang 
742b043167SHong Zhang M*/
752b043167SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Visualization(TSTrajectory tj,TS ts)
762b043167SHong Zhang {
772b043167SHong Zhang   PetscFunctionBegin;
782b043167SHong Zhang   tj->ops->set  = TSTrajectorySet_Visualization;
792b043167SHong Zhang   PetscFunctionReturn(0);
802b043167SHong Zhang }
81