xref: /petsc/src/ts/trajectory/impls/visualization/trajvisualization.c (revision e8e8640d1cb9a3a2f50c0c0d7b26e5c4d521e2e4)
1 #include <petsc/private/tsimpl.h> /*I "petscts.h"  I*/
2 
OutputBIN(MPI_Comm comm,const char * filename,PetscViewer * viewer)3 static PetscErrorCode OutputBIN(MPI_Comm comm, const char *filename, PetscViewer *viewer)
4 {
5   PetscFunctionBegin;
6   PetscCall(PetscViewerCreate(comm, viewer));
7   PetscCall(PetscViewerSetType(*viewer, PETSCVIEWERBINARY));
8   PetscCall(PetscViewerFileSetMode(*viewer, FILE_MODE_WRITE));
9   PetscCall(PetscViewerFileSetName(*viewer, filename));
10   PetscFunctionReturn(PETSC_SUCCESS);
11 }
12 
TSTrajectorySet_Visualization(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X)13 static PetscErrorCode TSTrajectorySet_Visualization(TSTrajectory tj, TS ts, PetscInt stepnum, PetscReal time, Vec X)
14 {
15   PetscViewer viewer;
16   char        filename[PETSC_MAX_PATH_LEN];
17   PetscReal   tprev;
18   MPI_Comm    comm;
19 
20   PetscFunctionBegin;
21   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
22   if (stepnum == 0) {
23     PetscMPIInt rank;
24     PetscCallMPI(MPI_Comm_rank(comm, &rank));
25     if (rank == 0) {
26       PetscCall(PetscRMTree("Visualization-data"));
27       PetscCall(PetscMkdir("Visualization-data"));
28     }
29     if (tj->names) {
30       PetscViewer bnames;
31       PetscCall(PetscViewerBinaryOpen(comm, "Visualization-data/variablenames", FILE_MODE_WRITE, &bnames));
32       PetscCall(PetscViewerBinaryWriteStringArray(bnames, (const char *const *)tj->names));
33       PetscCall(PetscViewerDestroy(&bnames));
34     }
35     PetscCall(PetscSNPrintf(filename, sizeof(filename), "Visualization-data/SA-%06" PetscInt_FMT ".bin", stepnum));
36     PetscCall(OutputBIN(comm, filename, &viewer));
37     if (!tj->transform) {
38       PetscCall(VecView(X, viewer));
39     } else {
40       Vec XX;
41       PetscCall((*tj->transform)(tj->transformctx, X, &XX));
42       PetscCall(VecView(XX, viewer));
43       PetscCall(VecDestroy(&XX));
44     }
45     PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL));
46     PetscCall(PetscViewerDestroy(&viewer));
47     PetscFunctionReturn(PETSC_SUCCESS);
48   }
49   PetscCall(PetscSNPrintf(filename, sizeof(filename), "Visualization-data/SA-%06" PetscInt_FMT ".bin", stepnum));
50   PetscCall(OutputBIN(comm, filename, &viewer));
51   if (!tj->transform) {
52     PetscCall(VecView(X, viewer));
53   } else {
54     Vec XX;
55     PetscCall((*tj->transform)(tj->transformctx, X, &XX));
56     PetscCall(VecView(XX, viewer));
57     PetscCall(VecDestroy(&XX));
58   }
59   PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL));
60 
61   PetscCall(TSGetPrevTime(ts, &tprev));
62   PetscCall(PetscViewerBinaryWrite(viewer, &tprev, 1, PETSC_REAL));
63 
64   PetscCall(PetscViewerDestroy(&viewer));
65   PetscFunctionReturn(PETSC_SUCCESS);
66 }
67 
68 /*MC
69       TSTRAJECTORYVISUALIZATION - Stores each solution of the ODE/DAE in a file
70 
71       Saves each timestep into a separate file in Visualization-data/SA-%06d.bin
72 
73       This version saves only the solutions at each timestep, it does not save the solution at each stage,
74       see `TSTRAJECTORYBASIC` that saves all stages
75 
76       $PETSC_DIR/share/petsc/matlab/PetscReadBinaryTrajectory.m and $PETSC_DIR/lib/petsc/bin/PetscBinaryIOTrajectory.py
77       can read in files created with this format into MATLAB and Python.
78 
79   Level: intermediate
80 
81 .seealso: [](ch_ts), `TSTrajectoryCreate()`, `TS`, `TSTrajectorySetType()`, `TSTrajectoryType`, `TSTrajectorySetVariableNames()`,
82           `TSTrajectoryType`, `TSTrajectory`
83 M*/
TSTrajectoryCreate_Visualization(TSTrajectory tj,TS ts)84 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Visualization(TSTrajectory tj, TS ts)
85 {
86   PetscFunctionBegin;
87   tj->ops->set    = TSTrajectorySet_Visualization;
88   tj->setupcalled = PETSC_TRUE;
89   PetscFunctionReturn(PETSC_SUCCESS);
90 }
91