xref: /petsc/src/ts/trajectory/impls/visualization/trajvisualization.c (revision 9566063d113dddea24716c546802770db7481bc0)
12b043167SHong Zhang 
22b043167SHong Zhang #include <petsc/private/tsimpl.h>        /*I "petscts.h"  I*/
32b043167SHong Zhang 
478fbdcc8SBarry Smith static PetscErrorCode OutputBIN(MPI_Comm comm,const char *filename,PetscViewer *viewer)
52b043167SHong Zhang {
62b043167SHong Zhang   PetscFunctionBegin;
7*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(comm,viewer));
8*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(*viewer,PETSCVIEWERBINARY));
9*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(*viewer,FILE_MODE_WRITE));
10*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetName(*viewer,filename));
112b043167SHong Zhang   PetscFunctionReturn(0);
122b043167SHong Zhang }
132b043167SHong Zhang 
14560360afSLisandro Dalcin static PetscErrorCode TSTrajectorySet_Visualization(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X)
152b043167SHong Zhang {
162b043167SHong Zhang   PetscViewer    viewer;
172b043167SHong Zhang   char           filename[PETSC_MAX_PATH_LEN];
182b043167SHong Zhang   PetscReal      tprev;
1978fbdcc8SBarry Smith   MPI_Comm       comm;
202b043167SHong Zhang 
212b043167SHong Zhang   PetscFunctionBegin;
22*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ts,&comm));
232b043167SHong Zhang   if (stepnum == 0) {
242b043167SHong Zhang     PetscMPIInt rank;
25*9566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm,&rank));
26dd400576SPatrick Sanan     if (rank == 0) {
27*9566063dSJacob Faibussowitsch       PetscCall(PetscRMTree("Visualization-data"));
28*9566063dSJacob Faibussowitsch       PetscCall(PetscMkdir("Visualization-data"));
292b043167SHong Zhang     }
3078fbdcc8SBarry Smith     if (tj->names) {
3178fbdcc8SBarry Smith       PetscViewer bnames;
32*9566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryOpen(comm,"Visualization-data/variablenames",FILE_MODE_WRITE,&bnames));
33*9566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryWriteStringArray(bnames,(const char *const *)tj->names));
34*9566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&bnames));
3578fbdcc8SBarry Smith     }
36*9566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(filename,sizeof(filename),"Visualization-data/SA-%06d.bin",stepnum));
37*9566063dSJacob Faibussowitsch     PetscCall(OutputBIN(comm,filename,&viewer));
3808347785SBarry Smith     if (!tj->transform) {
39*9566063dSJacob Faibussowitsch       PetscCall(VecView(X,viewer));
4008347785SBarry Smith     } else {
4108347785SBarry Smith       Vec XX;
42*9566063dSJacob Faibussowitsch       PetscCall((*tj->transform)(tj->transformctx,X,&XX));
43*9566063dSJacob Faibussowitsch       PetscCall(VecView(XX,viewer));
44*9566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&XX));
4508347785SBarry Smith     }
46*9566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL));
47*9566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
482b043167SHong Zhang     PetscFunctionReturn(0);
492b043167SHong Zhang   }
50*9566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(filename,sizeof(filename),"Visualization-data/SA-%06d.bin",stepnum));
51*9566063dSJacob Faibussowitsch   PetscCall(OutputBIN(comm,filename,&viewer));
5208347785SBarry Smith   if (!tj->transform) {
53*9566063dSJacob Faibussowitsch     PetscCall(VecView(X,viewer));
5408347785SBarry Smith   } else {
5508347785SBarry Smith     Vec XX;
56*9566063dSJacob Faibussowitsch     PetscCall((*tj->transform)(tj->transformctx,X,&XX));
57*9566063dSJacob Faibussowitsch     PetscCall(VecView(XX,viewer));
58*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&XX));
5908347785SBarry Smith   }
60*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL));
612b043167SHong Zhang 
62*9566063dSJacob Faibussowitsch   PetscCall(TSGetPrevTime(ts,&tprev));
63*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWrite(viewer,&tprev,1,PETSC_REAL));
642b043167SHong Zhang 
65*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
662b043167SHong Zhang   PetscFunctionReturn(0);
672b043167SHong Zhang }
682b043167SHong Zhang 
692b043167SHong Zhang /*MC
705c9db295SHong Zhang       TSTRAJECTORYVISUALIZATION - Stores each solution of the ODE/DAE in a file
712b043167SHong Zhang 
728e5aa403SBarry Smith       Saves each timestep into a separate file in Visualization-data/SA-%06d.bin
7378fbdcc8SBarry Smith 
7478fbdcc8SBarry Smith       This version saves only the solutions at each timestep, it does not save the solution at each stage,
7578fbdcc8SBarry Smith       see TSTRAJECTORYBASIC that saves all stages
7678fbdcc8SBarry Smith 
77c3a89c15SBarry Smith       $PETSC_DIR/share/petsc/matlab/PetscReadBinaryTrajectory.m and $PETSC_DIR/lib/petsc/bin/PetscBinaryIOTrajectory.py
7878fbdcc8SBarry Smith       can read in files created with this format into MATLAB and Python.
7978fbdcc8SBarry Smith 
802b043167SHong Zhang   Level: intermediate
812b043167SHong Zhang 
8278fbdcc8SBarry Smith .seealso:  TSTrajectoryCreate(), TS, TSTrajectorySetType(), TSTrajectoryType, TSTrajectorySetVariableNames()
832b043167SHong Zhang 
842b043167SHong Zhang M*/
852b043167SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Visualization(TSTrajectory tj,TS ts)
862b043167SHong Zhang {
872b043167SHong Zhang   PetscFunctionBegin;
882b043167SHong Zhang   tj->ops->set    = TSTrajectorySet_Visualization;
891aebe4aeSStefano Zampini   tj->setupcalled = PETSC_TRUE;
902b043167SHong Zhang   PetscFunctionReturn(0);
912b043167SHong Zhang }
92