xref: /petsc/src/sys/classes/viewer/impls/binary/binv.c (revision bf31d2d34e67170cccc52032cb9ed4707f3cb2ab)
1af0996ceSBarry Smith #include <petsc/private/viewerimpl.h> /*I   "petscviewer.h"   I*/
25c6c1daeSBarry Smith 
376667918SBarry Smith /*
476667918SBarry Smith    This needs to start the same as PetscViewer_Socket.
576667918SBarry Smith */
65c6c1daeSBarry Smith typedef struct {
75c6c1daeSBarry Smith   int       fdes;        /* file descriptor, ignored if using MPI IO */
876667918SBarry Smith   PetscInt  flowcontrol; /* allow only <flowcontrol> messages outstanding at a time while doing IO */
976667918SBarry Smith   PetscBool skipheader;  /* don't write header, only raw data */
105c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
11bc196f7cSDave May   PetscBool  usempiio;
125c6c1daeSBarry Smith   MPI_File   mfdes; /* ignored unless using MPI IO */
13e8a65b7dSLisandro Dalcin   MPI_File   mfsub; /* subviewer support */
145c6c1daeSBarry Smith   MPI_Offset moff;
155c6c1daeSBarry Smith #endif
16cc843e7aSLisandro Dalcin   char         *filename;            /* file name */
17cc843e7aSLisandro Dalcin   PetscFileMode filemode;            /* read/write/append mode */
185c6c1daeSBarry Smith   FILE         *fdes_info;           /* optional file containing info on binary file*/
195c6c1daeSBarry Smith   PetscBool     storecompressed;     /* gzip the write binary file when closing it*/
20f90597f1SStefano Zampini   char         *ogzfilename;         /* gzip can be run after the filename has been updated */
215c6c1daeSBarry Smith   PetscBool     skipinfo;            /* Don't create info file for writing; don't use for reading */
225c6c1daeSBarry Smith   PetscBool     skipoptions;         /* don't use PETSc options database when loading */
23a261c58fSBarry Smith   PetscBool     matlabheaderwritten; /* if format is PETSC_VIEWER_BINARY_MATLAB has the MATLAB .info header been written yet */
24c98fd787SBarry Smith   PetscBool     setfromoptionscalled;
255c6c1daeSBarry Smith } PetscViewer_Binary;
265c6c1daeSBarry Smith 
27d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryClearFunctionList(PetscViewer v)
28d71ae5a4SJacob Faibussowitsch {
292e956fe4SStefano Zampini   PetscFunctionBegin;
302e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetFlowControl_C", NULL));
312e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetFlowControl_C", NULL));
322e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipHeader_C", NULL));
332e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipHeader_C", NULL));
342e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipOptions_C", NULL));
352e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipOptions_C", NULL));
362e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipInfo_C", NULL));
372e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipInfo_C", NULL));
382e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetInfoPointer_C", NULL));
392e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", NULL));
402e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", NULL));
412e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", NULL));
422e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", NULL));
432e956fe4SStefano Zampini #if defined(PETSC_HAVE_MPIIO)
442e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetUseMPIIO_C", NULL));
452e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetUseMPIIO_C", NULL));
462e956fe4SStefano Zampini #endif
473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
482e956fe4SStefano Zampini }
492e956fe4SStefano Zampini 
50cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
51d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySyncMPIIO(PetscViewer viewer)
52d71ae5a4SJacob Faibussowitsch {
53cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
54cc843e7aSLisandro Dalcin 
55cc843e7aSLisandro Dalcin   PetscFunctionBegin;
563ba16761SJacob Faibussowitsch   if (vbinary->filemode == FILE_MODE_READ) PetscFunctionReturn(PETSC_SUCCESS);
5748a46eb9SPierre Jolivet   if (vbinary->mfsub != MPI_FILE_NULL) PetscCallMPI(MPI_File_sync(vbinary->mfsub));
58cc843e7aSLisandro Dalcin   if (vbinary->mfdes != MPI_FILE_NULL) {
599566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer)));
609566063dSJacob Faibussowitsch     PetscCallMPI(MPI_File_sync(vbinary->mfdes));
61cc843e7aSLisandro Dalcin   }
623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
63cc843e7aSLisandro Dalcin }
64cc843e7aSLisandro Dalcin #endif
65cc843e7aSLisandro Dalcin 
66d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerGetSubViewer_Binary(PetscViewer viewer, MPI_Comm comm, PetscViewer *outviewer)
67d71ae5a4SJacob Faibussowitsch {
68e8a65b7dSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
69cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
705c6c1daeSBarry Smith 
715c6c1daeSBarry Smith   PetscFunctionBegin;
729566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
73e8a65b7dSLisandro Dalcin 
74e8a65b7dSLisandro Dalcin   /* Return subviewer in process zero */
759566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
76dd400576SPatrick Sanan   if (rank == 0) {
77e5afcf28SBarry Smith     PetscMPIInt flg;
78e5afcf28SBarry Smith 
799566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_compare(PETSC_COMM_SELF, comm, &flg));
80cc73adaaSBarry Smith     PetscCheck(flg == MPI_IDENT || flg == MPI_CONGRUENT, PETSC_COMM_SELF, PETSC_ERR_SUP, "PetscViewerGetSubViewer() for PETSCVIEWERBINARY requires a singleton MPI_Comm");
819566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(comm, outviewer));
829566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(*outviewer, PETSCVIEWERBINARY));
839566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy((*outviewer)->data, vbinary, sizeof(PetscViewer_Binary)));
8412f4c3a9SLisandro Dalcin     (*outviewer)->setupcalled = PETSC_TRUE;
8596509d17SLisandro Dalcin   } else {
8696509d17SLisandro Dalcin     *outviewer = NULL;
8796509d17SLisandro Dalcin   }
88e8a65b7dSLisandro Dalcin 
89e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
90e8a65b7dSLisandro Dalcin   if (vbinary->usempiio && *outviewer) {
91e8a65b7dSLisandro Dalcin     PetscViewer_Binary *obinary = (PetscViewer_Binary *)(*outviewer)->data;
92e8a65b7dSLisandro Dalcin     /* Parent viewer opens a new MPI file handle on PETSC_COMM_SELF and keeps track of it for future reuse */
93e8a65b7dSLisandro Dalcin     if (vbinary->mfsub == MPI_FILE_NULL) {
94e8a65b7dSLisandro Dalcin       int amode;
95cc843e7aSLisandro Dalcin       switch (vbinary->filemode) {
96d71ae5a4SJacob Faibussowitsch       case FILE_MODE_READ:
97d71ae5a4SJacob Faibussowitsch         amode = MPI_MODE_RDONLY;
98d71ae5a4SJacob Faibussowitsch         break;
99d71ae5a4SJacob Faibussowitsch       case FILE_MODE_WRITE:
100d71ae5a4SJacob Faibussowitsch         amode = MPI_MODE_WRONLY;
101d71ae5a4SJacob Faibussowitsch         break;
102d71ae5a4SJacob Faibussowitsch       case FILE_MODE_APPEND:
103d71ae5a4SJacob Faibussowitsch         amode = MPI_MODE_WRONLY;
104d71ae5a4SJacob Faibussowitsch         break;
105d71ae5a4SJacob Faibussowitsch       default:
106d71ae5a4SJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[vbinary->filemode]);
107e8a65b7dSLisandro Dalcin       }
1089566063dSJacob Faibussowitsch       PetscCallMPI(MPI_File_open(PETSC_COMM_SELF, vbinary->filename, amode, MPI_INFO_NULL, &vbinary->mfsub));
109e8a65b7dSLisandro Dalcin     }
110e8a65b7dSLisandro Dalcin     /* Subviewer gets the MPI file handle on PETSC_COMM_SELF */
111e8a65b7dSLisandro Dalcin     obinary->mfdes = vbinary->mfsub;
112e8a65b7dSLisandro Dalcin     obinary->mfsub = MPI_FILE_NULL;
113e8a65b7dSLisandro Dalcin     obinary->moff  = vbinary->moff;
114e8a65b7dSLisandro Dalcin   }
115e8a65b7dSLisandro Dalcin #endif
116cc843e7aSLisandro Dalcin 
117cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
1189566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinarySyncMPIIO(viewer));
119cc843e7aSLisandro Dalcin #endif
1203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1215c6c1daeSBarry Smith }
1225c6c1daeSBarry Smith 
123d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerRestoreSubViewer_Binary(PetscViewer viewer, MPI_Comm comm, PetscViewer *outviewer)
124d71ae5a4SJacob Faibussowitsch {
125e8a65b7dSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
126cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
127e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
128e8a65b7dSLisandro Dalcin   MPI_Offset moff = 0;
129e8a65b7dSLisandro Dalcin #endif
1305c6c1daeSBarry Smith 
1315c6c1daeSBarry Smith   PetscFunctionBegin;
1329566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
133c5853193SPierre Jolivet   PetscCheck(rank == 0 || !*outviewer, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Subviewer not obtained from viewer");
134e8a65b7dSLisandro Dalcin 
135e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
136e8a65b7dSLisandro Dalcin   if (vbinary->usempiio && *outviewer) {
137e8a65b7dSLisandro Dalcin     PetscViewer_Binary *obinary = (PetscViewer_Binary *)(*outviewer)->data;
13808401ef6SPierre Jolivet     PetscCheck(obinary->mfdes == vbinary->mfsub, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Subviewer not obtained from viewer");
1399566063dSJacob Faibussowitsch     if (obinary->mfsub != MPI_FILE_NULL) PetscCallMPI(MPI_File_close(&obinary->mfsub));
140e8a65b7dSLisandro Dalcin     moff = obinary->moff;
141e8a65b7dSLisandro Dalcin   }
142e8a65b7dSLisandro Dalcin #endif
143e8a65b7dSLisandro Dalcin 
144e8a65b7dSLisandro Dalcin   if (*outviewer) {
145e8a65b7dSLisandro Dalcin     PetscViewer_Binary *obinary = (PetscViewer_Binary *)(*outviewer)->data;
14608401ef6SPierre Jolivet     PetscCheck(obinary->fdes == vbinary->fdes, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Subviewer not obtained from viewer");
1479566063dSJacob Faibussowitsch     PetscCall(PetscFree((*outviewer)->data));
1482e956fe4SStefano Zampini     PetscCall(PetscViewerBinaryClearFunctionList(*outviewer));
1499566063dSJacob Faibussowitsch     PetscCall(PetscHeaderDestroy(outviewer));
1505c6c1daeSBarry Smith   }
151e8a65b7dSLisandro Dalcin 
152e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
153e8a65b7dSLisandro Dalcin   if (vbinary->usempiio) {
154e8a65b7dSLisandro Dalcin     PetscInt64 ioff = (PetscInt64)moff; /* We could use MPI_OFFSET datatype (requires MPI 2.2) */
1559566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&ioff, 1, MPIU_INT64, 0, PetscObjectComm((PetscObject)viewer)));
156e8a65b7dSLisandro Dalcin     vbinary->moff = (MPI_Offset)ioff;
157e8a65b7dSLisandro Dalcin   }
158e8a65b7dSLisandro Dalcin #endif
159cc843e7aSLisandro Dalcin 
160cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
1619566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinarySyncMPIIO(viewer));
162cc843e7aSLisandro Dalcin #endif
1633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1645c6c1daeSBarry Smith }
1655c6c1daeSBarry Smith 
1665c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
1675c6c1daeSBarry Smith /*@C
16885b8072bSPatrick Sanan     PetscViewerBinaryGetMPIIOOffset - Gets the current global offset that should be passed to `MPI_File_set_view()` or `MPI_File_{write|read}_at[_all]()`
1695c6c1daeSBarry Smith 
170cf53795eSBarry Smith     Not Collective; No Fortran Support
1715c6c1daeSBarry Smith 
1725c6c1daeSBarry Smith     Input Parameter:
17385b8072bSPatrick Sanan .   viewer - PetscViewer context, obtained from `PetscViewerBinaryOpen()`
1745c6c1daeSBarry Smith 
1755c6c1daeSBarry Smith     Output Parameter:
17622a8f86cSLisandro Dalcin .   off - the current global offset
1775c6c1daeSBarry Smith 
1785c6c1daeSBarry Smith     Level: advanced
1795c6c1daeSBarry Smith 
180811af0c4SBarry Smith     Note:
181811af0c4SBarry Smith     Use `PetscViewerBinaryAddMPIIOOffset()` to increase this value after you have written a view.
182811af0c4SBarry Smith 
183d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryAddMPIIOOffset()`
1845c6c1daeSBarry Smith @*/
185d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetMPIIOOffset(PetscViewer viewer, MPI_Offset *off)
186d71ae5a4SJacob Faibussowitsch {
18722a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
1885c6c1daeSBarry Smith 
1895c6c1daeSBarry Smith   PetscFunctionBegin;
19022a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
19122a8f86cSLisandro Dalcin   PetscValidPointer(off, 2);
19222a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
1935c6c1daeSBarry Smith   *off    = vbinary->moff;
1943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1955c6c1daeSBarry Smith }
1965c6c1daeSBarry Smith 
1975c6c1daeSBarry Smith /*@C
19822a8f86cSLisandro Dalcin     PetscViewerBinaryAddMPIIOOffset - Adds to the current global offset
1995c6c1daeSBarry Smith 
200cf53795eSBarry Smith     Logically Collective; No Fortran Support
2015c6c1daeSBarry Smith 
2025c6c1daeSBarry Smith     Input Parameters:
203811af0c4SBarry Smith +   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
20422a8f86cSLisandro Dalcin -   off - the addition to the global offset
2055c6c1daeSBarry Smith 
2065c6c1daeSBarry Smith     Level: advanced
2075c6c1daeSBarry Smith 
208811af0c4SBarry Smith     Note:
209811af0c4SBarry Smith     Use `PetscViewerBinaryGetMPIIOOffset()` to get the value that you should pass to `MPI_File_set_view()` or `MPI_File_{write|read}_at[_all]()`
210811af0c4SBarry Smith 
211d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
2125c6c1daeSBarry Smith @*/
213d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryAddMPIIOOffset(PetscViewer viewer, MPI_Offset off)
214d71ae5a4SJacob Faibussowitsch {
21522a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
2165c6c1daeSBarry Smith 
2175c6c1daeSBarry Smith   PetscFunctionBegin;
21822a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
21922a8f86cSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, (PetscInt)off, 2);
22022a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
2215c6c1daeSBarry Smith   vbinary->moff += off;
2223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2235c6c1daeSBarry Smith }
2245c6c1daeSBarry Smith 
2255c6c1daeSBarry Smith /*@C
226811af0c4SBarry Smith     PetscViewerBinaryGetMPIIODescriptor - Extracts the MPI IO file descriptor from a `PetscViewer`.
2275c6c1daeSBarry Smith 
228cf53795eSBarry Smith     Not Collective; No Fortran Support
2295c6c1daeSBarry Smith 
2305c6c1daeSBarry Smith     Input Parameter:
231811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
2325c6c1daeSBarry Smith 
2335c6c1daeSBarry Smith     Output Parameter:
2345c6c1daeSBarry Smith .   fdes - file descriptor
2355c6c1daeSBarry Smith 
2365c6c1daeSBarry Smith     Level: advanced
2375c6c1daeSBarry Smith 
238d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
2395c6c1daeSBarry Smith @*/
240d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetMPIIODescriptor(PetscViewer viewer, MPI_File *fdes)
241d71ae5a4SJacob Faibussowitsch {
24222a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
2435c6c1daeSBarry Smith 
2445c6c1daeSBarry Smith   PetscFunctionBegin;
24522a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
24622a8f86cSLisandro Dalcin   PetscValidPointer(fdes, 2);
2479566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
24822a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
2495c6c1daeSBarry Smith   *fdes   = vbinary->mfdes;
2503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2515c6c1daeSBarry Smith }
252cc843e7aSLisandro Dalcin #endif
2535c6c1daeSBarry Smith 
254cc843e7aSLisandro Dalcin /*@
255cc843e7aSLisandro Dalcin     PetscViewerBinarySetUseMPIIO - Sets a binary viewer to use MPI-IO for reading/writing. Must be called
256811af0c4SBarry Smith         before `PetscViewerFileSetName()`
257cc843e7aSLisandro Dalcin 
258c3339decSBarry Smith     Logically Collective
259cc843e7aSLisandro Dalcin 
260cc843e7aSLisandro Dalcin     Input Parameters:
261811af0c4SBarry Smith +   viewer - the `PetscViewer`; must be a `PETSCVIEWERBINARY`
262811af0c4SBarry Smith -   use - `PETSC_TRUE` means MPI-IO will be used
263cc843e7aSLisandro Dalcin 
264811af0c4SBarry Smith     Options Database Key:
265cc843e7aSLisandro Dalcin     -viewer_binary_mpiio : Flag for using MPI-IO
266cc843e7aSLisandro Dalcin 
267cc843e7aSLisandro Dalcin     Level: advanced
268cc843e7aSLisandro Dalcin 
269d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`,
270db781477SPatrick Sanan           `PetscViewerBinaryGetUseMPIIO()`
271cc843e7aSLisandro Dalcin @*/
272d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetUseMPIIO(PetscViewer viewer, PetscBool use)
273d71ae5a4SJacob Faibussowitsch {
274a8aae444SDave May   PetscFunctionBegin;
275cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
276cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, use, 2);
277cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetUseMPIIO_C", (PetscViewer, PetscBool), (viewer, use));
2783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
279cc843e7aSLisandro Dalcin }
280cc843e7aSLisandro Dalcin 
281cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
282d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetUseMPIIO_Binary(PetscViewer viewer, PetscBool use)
283d71ae5a4SJacob Faibussowitsch {
284cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
285cc843e7aSLisandro Dalcin   PetscFunctionBegin;
286cc73adaaSBarry Smith   PetscCheck(!viewer->setupcalled || vbinary->usempiio == use, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Cannot change MPIIO to %s after setup", PetscBools[use]);
287cc843e7aSLisandro Dalcin   vbinary->usempiio = use;
2883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
289a8aae444SDave May }
290a8aae444SDave May #endif
291a8aae444SDave May 
292cc843e7aSLisandro Dalcin /*@
2933f423023SBarry Smith     PetscViewerBinaryGetUseMPIIO - Returns `PETSC_TRUE` if the binary viewer uses MPI-IO.
2945c6c1daeSBarry Smith 
2955c6c1daeSBarry Smith     Not Collective
2965c6c1daeSBarry Smith 
2975c6c1daeSBarry Smith     Input Parameter:
298811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`; must be a `PETSCVIEWERBINARY`
2995c6c1daeSBarry Smith 
3005c6c1daeSBarry Smith     Output Parameter:
301811af0c4SBarry Smith .   use - `PETSC_TRUE` if MPI-IO is being used
3025c6c1daeSBarry Smith 
3035c6c1daeSBarry Smith     Level: advanced
3045c6c1daeSBarry Smith 
305bc196f7cSDave May     Note:
3063f423023SBarry Smith     If MPI-IO is not available, this function will always return `PETSC_FALSE`
307bc196f7cSDave May 
308d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
3095c6c1daeSBarry Smith @*/
310d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetUseMPIIO(PetscViewer viewer, PetscBool *use)
311d71ae5a4SJacob Faibussowitsch {
3125c6c1daeSBarry Smith   PetscFunctionBegin;
313cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
314cc843e7aSLisandro Dalcin   PetscValidBoolPointer(use, 2);
315cc843e7aSLisandro Dalcin   *use = PETSC_FALSE;
316cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinaryGetUseMPIIO_C", (PetscViewer, PetscBool *), (viewer, use));
3173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3185c6c1daeSBarry Smith }
3195c6c1daeSBarry Smith 
320cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
321d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetUseMPIIO_Binary(PetscViewer viewer, PetscBool *use)
322d71ae5a4SJacob Faibussowitsch {
3235c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
3245c6c1daeSBarry Smith 
3255c6c1daeSBarry Smith   PetscFunctionBegin;
326cc843e7aSLisandro Dalcin   *use = vbinary->usempiio;
3273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
328cc843e7aSLisandro Dalcin }
329cc843e7aSLisandro Dalcin #endif
330cc843e7aSLisandro Dalcin 
331cc843e7aSLisandro Dalcin /*@
332cc843e7aSLisandro Dalcin     PetscViewerBinarySetFlowControl - Sets how many messages are allowed to outstanding at the same time during parallel IO reads/writes
333cc843e7aSLisandro Dalcin 
334cc843e7aSLisandro Dalcin     Not Collective
335cc843e7aSLisandro Dalcin 
336d8d19677SJose E. Roman     Input Parameters:
337811af0c4SBarry Smith +   viewer - PetscViewer context, obtained from `PetscViewerBinaryOpen()`
338cc843e7aSLisandro Dalcin -   fc - the number of messages, defaults to 256 if this function was not called
339cc843e7aSLisandro Dalcin 
340cc843e7aSLisandro Dalcin     Level: advanced
341cc843e7aSLisandro Dalcin 
342d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetFlowControl()`
343cc843e7aSLisandro Dalcin @*/
344d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetFlowControl(PetscViewer viewer, PetscInt fc)
345d71ae5a4SJacob Faibussowitsch {
346cc843e7aSLisandro Dalcin   PetscFunctionBegin;
347cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
348cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, fc, 2);
349cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetFlowControl_C", (PetscViewer, PetscInt), (viewer, fc));
3503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3515c6c1daeSBarry Smith }
3525c6c1daeSBarry Smith 
353d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetFlowControl_Binary(PetscViewer viewer, PetscInt fc)
354d71ae5a4SJacob Faibussowitsch {
355cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
356cc843e7aSLisandro Dalcin 
357cc843e7aSLisandro Dalcin   PetscFunctionBegin;
35808401ef6SPierre Jolivet   PetscCheck(fc > 1, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_OUTOFRANGE, "Flow control count must be greater than 1, %" PetscInt_FMT " was set", fc);
359cc843e7aSLisandro Dalcin   vbinary->flowcontrol = fc;
3603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
361cc843e7aSLisandro Dalcin }
362cc843e7aSLisandro Dalcin 
363cc843e7aSLisandro Dalcin /*@
3645c6c1daeSBarry Smith     PetscViewerBinaryGetFlowControl - Returns how many messages are allowed to outstanding at the same time during parallel IO reads/writes
3655c6c1daeSBarry Smith 
3665c6c1daeSBarry Smith     Not Collective
3675c6c1daeSBarry Smith 
3685c6c1daeSBarry Smith     Input Parameter:
369811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
3705c6c1daeSBarry Smith 
3715c6c1daeSBarry Smith     Output Parameter:
3725c6c1daeSBarry Smith .   fc - the number of messages
3735c6c1daeSBarry Smith 
3745c6c1daeSBarry Smith     Level: advanced
3755c6c1daeSBarry Smith 
376d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinarySetFlowControl()`
3775c6c1daeSBarry Smith @*/
378d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetFlowControl(PetscViewer viewer, PetscInt *fc)
379d71ae5a4SJacob Faibussowitsch {
3805c6c1daeSBarry Smith   PetscFunctionBegin;
381cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
382cc843e7aSLisandro Dalcin   PetscValidIntPointer(fc, 2);
383cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetFlowControl_C", (PetscViewer, PetscInt *), (viewer, fc));
3843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3855c6c1daeSBarry Smith }
3865c6c1daeSBarry Smith 
38776667918SBarry Smith PETSC_INTERN PetscErrorCode PetscViewerBinaryGetFlowControl_Binary(PetscViewer viewer, PetscInt *fc)
388d71ae5a4SJacob Faibussowitsch {
3895c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
3905c6c1daeSBarry Smith 
3915c6c1daeSBarry Smith   PetscFunctionBegin;
392cc843e7aSLisandro Dalcin   *fc = vbinary->flowcontrol;
3933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3945c6c1daeSBarry Smith }
3955c6c1daeSBarry Smith 
3965c6c1daeSBarry Smith /*@C
397811af0c4SBarry Smith     PetscViewerBinaryGetDescriptor - Extracts the file descriptor from a `PetscViewer` of `PetscViewerType` `PETSCVIEWERBINARY`.
3985c6c1daeSBarry Smith 
399cd05f99aSJacob Faibussowitsch     Collective because it may trigger a `PetscViewerSetUp()` call; No Fortran Support
4005c6c1daeSBarry Smith 
4015c6c1daeSBarry Smith     Input Parameter:
402811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
4035c6c1daeSBarry Smith 
4045c6c1daeSBarry Smith     Output Parameter:
4055c6c1daeSBarry Smith .   fdes - file descriptor
4065c6c1daeSBarry Smith 
4075c6c1daeSBarry Smith     Level: advanced
4085c6c1daeSBarry Smith 
409811af0c4SBarry Smith     Note:
410811af0c4SBarry Smith       For writable binary `PetscViewer`s, the descriptor will only be valid for the
411811af0c4SBarry Smith     first processor in the communicator that shares the `PetscViewer`. For readable
4125c6c1daeSBarry Smith     files it will only be valid on nodes that have the file. If node 0 does not
4135c6c1daeSBarry Smith     have the file it generates an error even if another node does have the file.
4145c6c1daeSBarry Smith 
415d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`
4165c6c1daeSBarry Smith @*/
417d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetDescriptor(PetscViewer viewer, int *fdes)
418d71ae5a4SJacob Faibussowitsch {
41922a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
4205c6c1daeSBarry Smith 
4215c6c1daeSBarry Smith   PetscFunctionBegin;
42222a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
42322a8f86cSLisandro Dalcin   PetscValidPointer(fdes, 2);
4249566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
42522a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
4265c6c1daeSBarry Smith   *fdes   = vbinary->fdes;
4273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4285c6c1daeSBarry Smith }
4295c6c1daeSBarry Smith 
4305c6c1daeSBarry Smith /*@
4315c6c1daeSBarry Smith     PetscViewerBinarySkipInfo - Binary file will not have .info file created with it
4325c6c1daeSBarry Smith 
4335c6c1daeSBarry Smith     Not Collective
4345c6c1daeSBarry Smith 
435fd292e60Sprj-     Input Parameter:
436811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerCreate()`
4375c6c1daeSBarry Smith 
4385c6c1daeSBarry Smith     Options Database Key:
43910699b91SBarry Smith .   -viewer_binary_skip_info - true indicates do not generate .info file
4405c6c1daeSBarry Smith 
4415c6c1daeSBarry Smith     Level: advanced
4425c6c1daeSBarry Smith 
44395452b02SPatrick Sanan     Notes:
444811af0c4SBarry Smith     This must be called after `PetscViewerSetType()`. If you use `PetscViewerBinaryOpen()` then
4453f423023SBarry Smith     you can only skip the info file with the `-viewer_binary_skip_info` flag. To use the function you must open the
446811af0c4SBarry Smith     viewer with `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinarySkipInfo()`.
4475c6c1daeSBarry Smith 
4485c6c1daeSBarry Smith     The .info contains meta information about the data in the binary file, for example the block size if it was
4495c6c1daeSBarry Smith     set for a vector or matrix.
4505c6c1daeSBarry Smith 
451811af0c4SBarry Smith     This routine is deprecated, use `PetscViewerBinarySetSkipInfo()`
452811af0c4SBarry Smith 
453d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySetSkipOptions()`,
454db781477SPatrick Sanan           `PetscViewerBinaryGetSkipOptions()`, `PetscViewerBinaryGetSkipInfo()`
4555c6c1daeSBarry Smith @*/
456d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySkipInfo(PetscViewer viewer)
457d71ae5a4SJacob Faibussowitsch {
4585c6c1daeSBarry Smith   PetscFunctionBegin;
4599566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinarySetSkipInfo(viewer, PETSC_TRUE));
4603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
461807ea322SDave May }
462807ea322SDave May 
463807ea322SDave May /*@
464807ea322SDave May     PetscViewerBinarySetSkipInfo - Binary file will not have .info file created with it
465807ea322SDave May 
466807ea322SDave May     Not Collective
467807ea322SDave May 
468d8d19677SJose E. Roman     Input Parameters:
4693f423023SBarry Smith +   viewer - PetscViewer context, obtained from `PetscViewerCreate()`
4703f423023SBarry Smith -   skip - `PETSC_TRUE` implies the .info file will not be generated
471807ea322SDave May 
472807ea322SDave May     Options Database Key:
47310699b91SBarry Smith .   -viewer_binary_skip_info - true indicates do not generate .info file
474807ea322SDave May 
475807ea322SDave May     Level: advanced
476807ea322SDave May 
477d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySetSkipOptions()`,
478db781477SPatrick Sanan           `PetscViewerBinaryGetSkipOptions()`, `PetscViewerBinaryGetSkipInfo()`
479807ea322SDave May @*/
480d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetSkipInfo(PetscViewer viewer, PetscBool skip)
481d71ae5a4SJacob Faibussowitsch {
482807ea322SDave May   PetscFunctionBegin;
483cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
484cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, skip, 2);
485cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetSkipInfo_C", (PetscViewer, PetscBool), (viewer, skip));
4863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
487807ea322SDave May }
488807ea322SDave May 
489d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetSkipInfo_Binary(PetscViewer viewer, PetscBool skip)
490d71ae5a4SJacob Faibussowitsch {
491807ea322SDave May   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
492807ea322SDave May 
493807ea322SDave May   PetscFunctionBegin;
494cc843e7aSLisandro Dalcin   vbinary->skipinfo = skip;
4953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
496807ea322SDave May }
497807ea322SDave May 
498807ea322SDave May /*@
499807ea322SDave May     PetscViewerBinaryGetSkipInfo - check if viewer wrote a .info file
500807ea322SDave May 
501807ea322SDave May     Not Collective
502807ea322SDave May 
503807ea322SDave May     Input Parameter:
504811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
505807ea322SDave May 
506807ea322SDave May     Output Parameter:
507811af0c4SBarry Smith .   skip - `PETSC_TRUE` implies the .info file was not generated
508807ea322SDave May 
509807ea322SDave May     Level: advanced
510807ea322SDave May 
511811af0c4SBarry Smith     Note:
512811af0c4SBarry Smith     This must be called after `PetscViewerSetType()`
513807ea322SDave May 
514d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`,
515db781477SPatrick Sanan           `PetscViewerBinarySetSkipOptions()`, `PetscViewerBinarySetSkipInfo()`
516807ea322SDave May @*/
517d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetSkipInfo(PetscViewer viewer, PetscBool *skip)
518d71ae5a4SJacob Faibussowitsch {
519807ea322SDave May   PetscFunctionBegin;
520cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
521cc843e7aSLisandro Dalcin   PetscValidBoolPointer(skip, 2);
522cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetSkipInfo_C", (PetscViewer, PetscBool *), (viewer, skip));
5233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
524807ea322SDave May }
525807ea322SDave May 
526d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetSkipInfo_Binary(PetscViewer viewer, PetscBool *skip)
527d71ae5a4SJacob Faibussowitsch {
528807ea322SDave May   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
529807ea322SDave May 
530807ea322SDave May   PetscFunctionBegin;
531cc843e7aSLisandro Dalcin   *skip = vbinary->skipinfo;
5323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
533807ea322SDave May }
534807ea322SDave May 
5355c6c1daeSBarry Smith /*@
5365c6c1daeSBarry Smith     PetscViewerBinarySetSkipOptions - do not use the PETSc options database when loading objects
5375c6c1daeSBarry Smith 
5385c6c1daeSBarry Smith     Not Collective
5395c6c1daeSBarry Smith 
5405c6c1daeSBarry Smith     Input Parameters:
541811af0c4SBarry Smith +   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
542811af0c4SBarry Smith -   skip - `PETSC_TRUE` means do not use the options from the options database
5435c6c1daeSBarry Smith 
5445c6c1daeSBarry Smith     Options Database Key:
545811af0c4SBarry Smith .   -viewer_binary_skip_options <true or false> - true means do not use the options from the options database
5465c6c1daeSBarry Smith 
5475c6c1daeSBarry Smith     Level: advanced
5485c6c1daeSBarry Smith 
549811af0c4SBarry Smith     Note:
550811af0c4SBarry Smith     This must be called after `PetscViewerSetType()`
5515c6c1daeSBarry Smith 
552d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
553db781477SPatrick Sanan           `PetscViewerBinaryGetSkipOptions()`
5545c6c1daeSBarry Smith @*/
555d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetSkipOptions(PetscViewer viewer, PetscBool skip)
556d71ae5a4SJacob Faibussowitsch {
5575c6c1daeSBarry Smith   PetscFunctionBegin;
558cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
559cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, skip, 2);
560cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetSkipOptions_C", (PetscViewer, PetscBool), (viewer, skip));
5613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
562807ea322SDave May }
563807ea322SDave May 
564d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetSkipOptions_Binary(PetscViewer viewer, PetscBool skip)
565d71ae5a4SJacob Faibussowitsch {
566cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
567807ea322SDave May 
568807ea322SDave May   PetscFunctionBegin;
569cc843e7aSLisandro Dalcin   vbinary->skipoptions = skip;
5703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5715c6c1daeSBarry Smith }
5725c6c1daeSBarry Smith 
5735c6c1daeSBarry Smith /*@
5745c6c1daeSBarry Smith     PetscViewerBinaryGetSkipOptions - checks if viewer uses the PETSc options database when loading objects
5755c6c1daeSBarry Smith 
5765c6c1daeSBarry Smith     Not Collective
5775c6c1daeSBarry Smith 
5785c6c1daeSBarry Smith     Input Parameter:
579811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
5805c6c1daeSBarry Smith 
5815c6c1daeSBarry Smith     Output Parameter:
582811af0c4SBarry Smith .   skip - `PETSC_TRUE` means do not use
5835c6c1daeSBarry Smith 
5845c6c1daeSBarry Smith     Level: advanced
5855c6c1daeSBarry Smith 
586811af0c4SBarry Smith     Note:
587811af0c4SBarry Smith     This must be called after `PetscViewerSetType()`
5885c6c1daeSBarry Smith 
589d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
590db781477SPatrick Sanan           `PetscViewerBinarySetSkipOptions()`
5915c6c1daeSBarry Smith @*/
592d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetSkipOptions(PetscViewer viewer, PetscBool *skip)
593d71ae5a4SJacob Faibussowitsch {
5945c6c1daeSBarry Smith   PetscFunctionBegin;
595cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
596cc843e7aSLisandro Dalcin   PetscValidBoolPointer(skip, 2);
597cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetSkipOptions_C", (PetscViewer, PetscBool *), (viewer, skip));
5983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5995c6c1daeSBarry Smith }
6005c6c1daeSBarry Smith 
601d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetSkipOptions_Binary(PetscViewer viewer, PetscBool *skip)
602d71ae5a4SJacob Faibussowitsch {
603d21b9a37SPierre Jolivet   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
6045c6c1daeSBarry Smith 
6055c6c1daeSBarry Smith   PetscFunctionBegin;
606cc843e7aSLisandro Dalcin   *skip = vbinary->skipoptions;
6073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6085c6c1daeSBarry Smith }
6095c6c1daeSBarry Smith 
6105c6c1daeSBarry Smith /*@
6115c6c1daeSBarry Smith     PetscViewerBinarySetSkipHeader - do not write a header with size information on output, just raw data
6125c6c1daeSBarry Smith 
6135c6c1daeSBarry Smith     Not Collective
6145c6c1daeSBarry Smith 
6155c6c1daeSBarry Smith     Input Parameters:
616811af0c4SBarry Smith +   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
617811af0c4SBarry Smith -   skip - `PETSC_TRUE` means do not write header
6185c6c1daeSBarry Smith 
6195c6c1daeSBarry Smith     Options Database Key:
620811af0c4SBarry Smith .   -viewer_binary_skip_header <true or false> - true means do not write header
6215c6c1daeSBarry Smith 
6225c6c1daeSBarry Smith     Level: advanced
6235c6c1daeSBarry Smith 
6243f423023SBarry Smith     Note:
625811af0c4SBarry Smith       This must be called after `PetscViewerSetType()`
6265c6c1daeSBarry Smith 
627d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
628db781477SPatrick Sanan           `PetscViewerBinaryGetSkipHeader()`
6295c6c1daeSBarry Smith @*/
630d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetSkipHeader(PetscViewer viewer, PetscBool skip)
631d71ae5a4SJacob Faibussowitsch {
6325c6c1daeSBarry Smith   PetscFunctionBegin;
633cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
634cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, skip, 2);
635cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetSkipHeader_C", (PetscViewer, PetscBool), (viewer, skip));
6363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6375c6c1daeSBarry Smith }
6385c6c1daeSBarry Smith 
639d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetSkipHeader_Binary(PetscViewer viewer, PetscBool skip)
640d71ae5a4SJacob Faibussowitsch {
6415c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
6425c6c1daeSBarry Smith 
6435c6c1daeSBarry Smith   PetscFunctionBegin;
644cc843e7aSLisandro Dalcin   vbinary->skipheader = skip;
6453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6465c6c1daeSBarry Smith }
6475c6c1daeSBarry Smith 
6485c6c1daeSBarry Smith /*@
6495c6c1daeSBarry Smith     PetscViewerBinaryGetSkipHeader - checks whether to write a header with size information on output, or just raw data
6505c6c1daeSBarry Smith 
6515c6c1daeSBarry Smith     Not Collective
6525c6c1daeSBarry Smith 
6535c6c1daeSBarry Smith     Input Parameter:
654811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
6555c6c1daeSBarry Smith 
6565c6c1daeSBarry Smith     Output Parameter:
657811af0c4SBarry Smith .   skip - `PETSC_TRUE` means do not write header
6585c6c1daeSBarry Smith 
6595c6c1daeSBarry Smith     Level: advanced
6605c6c1daeSBarry Smith 
66195452b02SPatrick Sanan     Notes:
66295452b02SPatrick Sanan     This must be called after PetscViewerSetType()
6635c6c1daeSBarry Smith 
664811af0c4SBarry Smith     Returns `PETSC_FALSE` for `PETSCSOCKETVIEWER`, you cannot skip the header for it.
6655c6c1daeSBarry Smith 
666d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
667db781477SPatrick Sanan           `PetscViewerBinarySetSkipHeader()`
6685c6c1daeSBarry Smith @*/
669d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetSkipHeader(PetscViewer viewer, PetscBool *skip)
670d71ae5a4SJacob Faibussowitsch {
6715c6c1daeSBarry Smith   PetscFunctionBegin;
672cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
673cc843e7aSLisandro Dalcin   PetscValidBoolPointer(skip, 2);
674cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetSkipHeader_C", (PetscViewer, PetscBool *), (viewer, skip));
6753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6765c6c1daeSBarry Smith }
6775c6c1daeSBarry Smith 
678d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetSkipHeader_Binary(PetscViewer viewer, PetscBool *skip)
679d71ae5a4SJacob Faibussowitsch {
680f3b211e4SSatish Balay   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
6815c6c1daeSBarry Smith 
6825c6c1daeSBarry Smith   PetscFunctionBegin;
683cc843e7aSLisandro Dalcin   *skip = vbinary->skipheader;
6843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6855c6c1daeSBarry Smith }
6865c6c1daeSBarry Smith 
6875c6c1daeSBarry Smith /*@C
6885c6c1daeSBarry Smith     PetscViewerBinaryGetInfoPointer - Extracts the file pointer for the ASCII
6895c6c1daeSBarry Smith           info file associated with a binary file.
6905c6c1daeSBarry Smith 
691cf53795eSBarry Smith     Not Collective; No Fortran Support
6925c6c1daeSBarry Smith 
6935c6c1daeSBarry Smith     Input Parameter:
694811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
6955c6c1daeSBarry Smith 
6965c6c1daeSBarry Smith     Output Parameter:
6970298fd71SBarry Smith .   file - file pointer  Always returns NULL if not a binary viewer
6985c6c1daeSBarry Smith 
6995c6c1daeSBarry Smith     Level: advanced
7005c6c1daeSBarry Smith 
701811af0c4SBarry Smith     Note:
702811af0c4SBarry Smith       For writable binary `PetscViewer`s, the file pointer will only be valid for the
703811af0c4SBarry Smith     first processor in the communicator that shares the `PetscViewer`.
7045c6c1daeSBarry Smith 
705d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`
7065c6c1daeSBarry Smith @*/
707d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetInfoPointer(PetscViewer viewer, FILE **file)
708d71ae5a4SJacob Faibussowitsch {
7095c6c1daeSBarry Smith   PetscFunctionBegin;
710cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
711cc843e7aSLisandro Dalcin   PetscValidPointer(file, 2);
7120298fd71SBarry Smith   *file = NULL;
713cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinaryGetInfoPointer_C", (PetscViewer, FILE **), (viewer, file));
7143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7155c6c1daeSBarry Smith }
7165c6c1daeSBarry Smith 
717d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetInfoPointer_Binary(PetscViewer viewer, FILE **file)
718d71ae5a4SJacob Faibussowitsch {
719cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
7205c6c1daeSBarry Smith 
7215c6c1daeSBarry Smith   PetscFunctionBegin;
7229566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
723cc843e7aSLisandro Dalcin   *file = vbinary->fdes_info;
724cc843e7aSLisandro Dalcin   if (viewer->format == PETSC_VIEWER_BINARY_MATLAB && !vbinary->matlabheaderwritten) {
7255c6c1daeSBarry Smith     if (vbinary->fdes_info) {
726cc843e7aSLisandro Dalcin       FILE *info = vbinary->fdes_info;
7279566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- begin code written by PetscViewerBinary for MATLAB format ---#\n"));
7289566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#$$ Set.filename = '%s';\n", vbinary->filename));
7299566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#$$ fd = PetscOpenFile(Set.filename);\n"));
7309566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n"));
731cc843e7aSLisandro Dalcin     }
732cc843e7aSLisandro Dalcin     vbinary->matlabheaderwritten = PETSC_TRUE;
7335c6c1daeSBarry Smith   }
7343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7355c6c1daeSBarry Smith }
7365c6c1daeSBarry Smith 
737e0385b85SDave May #if defined(PETSC_HAVE_MPIIO)
738d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_BinaryMPIIO(PetscViewer v)
739d71ae5a4SJacob Faibussowitsch {
740e0385b85SDave May   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
741e0385b85SDave May 
742e0385b85SDave May   PetscFunctionBegin;
74348a46eb9SPierre Jolivet   if (vbinary->mfdes != MPI_FILE_NULL) PetscCallMPI(MPI_File_close(&vbinary->mfdes));
74448a46eb9SPierre Jolivet   if (vbinary->mfsub != MPI_FILE_NULL) PetscCallMPI(MPI_File_close(&vbinary->mfsub));
745cc843e7aSLisandro Dalcin   vbinary->moff = 0;
7463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
747e0385b85SDave May }
748e0385b85SDave May #endif
749e0385b85SDave May 
750d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_BinarySTDIO(PetscViewer v)
751d71ae5a4SJacob Faibussowitsch {
752cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
753cc843e7aSLisandro Dalcin 
754cc843e7aSLisandro Dalcin   PetscFunctionBegin;
755cc843e7aSLisandro Dalcin   if (vbinary->fdes != -1) {
7569566063dSJacob Faibussowitsch     PetscCall(PetscBinaryClose(vbinary->fdes));
757cc843e7aSLisandro Dalcin     vbinary->fdes = -1;
758cc843e7aSLisandro Dalcin     if (vbinary->storecompressed) {
759cc843e7aSLisandro Dalcin       char        cmd[8 + PETSC_MAX_PATH_LEN], out[64 + PETSC_MAX_PATH_LEN] = "";
760cc843e7aSLisandro Dalcin       const char *gzfilename = vbinary->ogzfilename ? vbinary->ogzfilename : vbinary->filename;
761cc843e7aSLisandro Dalcin       /* compress the file */
7629566063dSJacob Faibussowitsch       PetscCall(PetscStrncpy(cmd, "gzip -f ", sizeof(cmd)));
7639566063dSJacob Faibussowitsch       PetscCall(PetscStrlcat(cmd, gzfilename, sizeof(cmd)));
764cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_POPEN)
765cc843e7aSLisandro Dalcin       {
766cc843e7aSLisandro Dalcin         FILE *fp;
7679566063dSJacob Faibussowitsch         PetscCall(PetscPOpen(PETSC_COMM_SELF, NULL, cmd, "r", &fp));
768cc73adaaSBarry Smith         PetscCheck(!fgets(out, (int)(sizeof(out) - 1), fp), PETSC_COMM_SELF, PETSC_ERR_LIB, "Error from command %s\n%s", cmd, out);
7699566063dSJacob Faibussowitsch         PetscCall(PetscPClose(PETSC_COMM_SELF, fp));
770cc843e7aSLisandro Dalcin       }
771cc843e7aSLisandro Dalcin #endif
772cc843e7aSLisandro Dalcin     }
773cc843e7aSLisandro Dalcin   }
7749566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary->ogzfilename));
7753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
776cc843e7aSLisandro Dalcin }
777cc843e7aSLisandro Dalcin 
778d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_BinaryInfo(PetscViewer v)
779d71ae5a4SJacob Faibussowitsch {
780cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
781cc843e7aSLisandro Dalcin 
782cc843e7aSLisandro Dalcin   PetscFunctionBegin;
783cc843e7aSLisandro Dalcin   if (v->format == PETSC_VIEWER_BINARY_MATLAB && vbinary->matlabheaderwritten) {
784cc843e7aSLisandro Dalcin     if (vbinary->fdes_info) {
785cc843e7aSLisandro Dalcin       FILE *info = vbinary->fdes_info;
7869566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- begin code written by PetscViewerBinary for MATLAB format ---#\n"));
7879566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#$$ close(fd);\n"));
7889566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n"));
789cc843e7aSLisandro Dalcin     }
790cc843e7aSLisandro Dalcin   }
791cc843e7aSLisandro Dalcin   if (vbinary->fdes_info) {
792cc843e7aSLisandro Dalcin     FILE *info         = vbinary->fdes_info;
793cc843e7aSLisandro Dalcin     vbinary->fdes_info = NULL;
794cc73adaaSBarry Smith     PetscCheck(!fclose(info), PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file");
795cc843e7aSLisandro Dalcin   }
7963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
797cc843e7aSLisandro Dalcin }
798cc843e7aSLisandro Dalcin 
799d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_Binary(PetscViewer v)
800d71ae5a4SJacob Faibussowitsch {
801cc843e7aSLisandro Dalcin   PetscFunctionBegin;
802cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
8039566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_BinaryMPIIO(v));
804cc843e7aSLisandro Dalcin #endif
8059566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_BinarySTDIO(v));
8069566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_BinaryInfo(v));
8073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
808cc843e7aSLisandro Dalcin }
809cc843e7aSLisandro Dalcin 
810d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerDestroy_Binary(PetscViewer v)
811d71ae5a4SJacob Faibussowitsch {
8125c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
8135c6c1daeSBarry Smith 
8145c6c1daeSBarry Smith   PetscFunctionBegin;
8159566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_Binary(v));
8169566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary->filename));
8179566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary));
8182e956fe4SStefano Zampini   PetscCall(PetscViewerBinaryClearFunctionList(v));
8193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
820e0385b85SDave May }
8215c6c1daeSBarry Smith 
8225c6c1daeSBarry Smith /*@C
8235c6c1daeSBarry Smith    PetscViewerBinaryOpen - Opens a file for binary input/output.
8245c6c1daeSBarry Smith 
825d083f849SBarry Smith    Collective
8265c6c1daeSBarry Smith 
8275c6c1daeSBarry Smith    Input Parameters:
8285c6c1daeSBarry Smith +  comm - MPI communicator
8295c6c1daeSBarry Smith .  name - name of file
830cc843e7aSLisandro Dalcin -  mode - open mode of file
8313f423023SBarry Smith .vb
8323f423023SBarry Smith     FILE_MODE_WRITE - create new file for binary output
8333f423023SBarry Smith     FILE_MODE_READ - open existing file for binary input
8343f423023SBarry Smith     FILE_MODE_APPEND - open existing file for binary output
8353f423023SBarry Smith .ve
8365c6c1daeSBarry Smith 
8375c6c1daeSBarry Smith    Output Parameter:
838cc843e7aSLisandro Dalcin .  viewer - PetscViewer for binary input/output to use with the specified file
8395c6c1daeSBarry Smith 
8405c6c1daeSBarry Smith     Options Database Keys:
841811af0c4SBarry Smith +    -viewer_binary_filename <name> - name of file to use
842811af0c4SBarry Smith .    -viewer_binary_skip_info - true to skip opening an info file
843811af0c4SBarry Smith .    -viewer_binary_skip_options - true to not use options database while creating viewer
844811af0c4SBarry Smith .    -viewer_binary_skip_header - true to skip output object headers to the file
845811af0c4SBarry Smith -    -viewer_binary_mpiio - true to use MPI-IO for input and output to the file (more scalable for large problems)
8465c6c1daeSBarry Smith 
8475c6c1daeSBarry Smith    Level: beginner
8485c6c1daeSBarry Smith 
8495c6c1daeSBarry Smith    Note:
850811af0c4SBarry Smith    This `PetscViewer` should be destroyed with `PetscViewerDestroy()`.
8515c6c1daeSBarry Smith 
8525c6c1daeSBarry Smith     For reading files, the filename may begin with ftp:// or http:// and/or
8535c6c1daeSBarry Smith     end with .gz; in this case file is brought over and uncompressed.
8545c6c1daeSBarry Smith 
8555c6c1daeSBarry Smith     For creating files, if the file name ends with .gz it is automatically
8565c6c1daeSBarry Smith     compressed when closed.
8575c6c1daeSBarry Smith 
858d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
859db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
860db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`, `PetscViewerBinarySetUseMPIIO()`,
861db781477SPatrick Sanan           `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
8625c6c1daeSBarry Smith @*/
863d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryOpen(MPI_Comm comm, const char name[], PetscFileMode mode, PetscViewer *viewer)
864d71ae5a4SJacob Faibussowitsch {
8655c6c1daeSBarry Smith   PetscFunctionBegin;
8669566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(comm, viewer));
8679566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(*viewer, PETSCVIEWERBINARY));
8689566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(*viewer, mode));
8699566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetName(*viewer, name));
8709566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions(*viewer));
8713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8725c6c1daeSBarry Smith }
8735c6c1daeSBarry Smith 
8745c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
875d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryWriteReadMPIIO(PetscViewer viewer, void *data, PetscInt num, PetscInt *count, PetscDataType dtype, PetscBool write)
876d71ae5a4SJacob Faibussowitsch {
87722a8f86cSLisandro Dalcin   MPI_Comm            comm    = PetscObjectComm((PetscObject)viewer);
8785c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
87922a8f86cSLisandro Dalcin   MPI_File            mfdes   = vbinary->mfdes;
8805c6c1daeSBarry Smith   MPI_Datatype        mdtype;
88169e10bbaSLisandro Dalcin   PetscMPIInt         rank, cnt;
8825c6c1daeSBarry Smith   MPI_Status          status;
8835c6c1daeSBarry Smith   MPI_Aint            ul, dsize;
8845c6c1daeSBarry Smith 
8855c6c1daeSBarry Smith   PetscFunctionBegin;
8869566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
8879566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(num, &cnt));
8889566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToMPIDataType(dtype, &mdtype));
8895c6c1daeSBarry Smith   if (write) {
89048a46eb9SPierre Jolivet     if (rank == 0) PetscCall(MPIU_File_write_at(mfdes, vbinary->moff, data, cnt, mdtype, &status));
8915c6c1daeSBarry Smith   } else {
892dd400576SPatrick Sanan     if (rank == 0) {
8939566063dSJacob Faibussowitsch       PetscCall(MPIU_File_read_at(mfdes, vbinary->moff, data, cnt, mdtype, &status));
8949566063dSJacob Faibussowitsch       if (cnt > 0) PetscCallMPI(MPI_Get_count(&status, mdtype, &cnt));
8955c6c1daeSBarry Smith     }
8969566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&cnt, 1, MPI_INT, 0, comm));
8979566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(data, cnt, mdtype, 0, comm));
89869e10bbaSLisandro Dalcin   }
8999566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_get_extent(mdtype, &ul, &dsize));
9005c6c1daeSBarry Smith   vbinary->moff += dsize * cnt;
9019860990eSLisandro Dalcin   if (count) *count = cnt;
9023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9035c6c1daeSBarry Smith }
9045c6c1daeSBarry Smith #endif
9055c6c1daeSBarry Smith 
9065c6c1daeSBarry Smith /*@C
9075c6c1daeSBarry Smith    PetscViewerBinaryRead - Reads from a binary file, all processors get the same result
9085c6c1daeSBarry Smith 
909d083f849SBarry Smith    Collective
9105c6c1daeSBarry Smith 
9115c6c1daeSBarry Smith    Input Parameters:
9125c6c1daeSBarry Smith +  viewer - the binary viewer
913907376e6SBarry Smith .  data - location of the data to be written
914060da220SMatthew G. Knepley .  num - number of items of data to read
915907376e6SBarry Smith -  dtype - type of data to read
9165c6c1daeSBarry Smith 
9172fe279fdSBarry Smith    Output Parameter:
9183f423023SBarry Smith .  count - number of items of data actually read, or `NULL`.
919f8e4bde8SMatthew G. Knepley 
9205c6c1daeSBarry Smith    Level: beginner
9215c6c1daeSBarry Smith 
922d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
923db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
924db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
9255c6c1daeSBarry Smith @*/
926d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryRead(PetscViewer viewer, void *data, PetscInt num, PetscInt *count, PetscDataType dtype)
927d71ae5a4SJacob Faibussowitsch {
92822a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
9295c6c1daeSBarry Smith 
93022a8f86cSLisandro Dalcin   PetscFunctionBegin;
93122a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
93222a8f86cSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, num, 3);
9339566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
93422a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
9355c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
936bc196f7cSDave May   if (vbinary->usempiio) {
9379566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWriteReadMPIIO(viewer, data, num, count, dtype, PETSC_FALSE));
9385c6c1daeSBarry Smith   } else {
9395c6c1daeSBarry Smith #endif
9409566063dSJacob Faibussowitsch     PetscCall(PetscBinarySynchronizedRead(PetscObjectComm((PetscObject)viewer), vbinary->fdes, data, num, count, dtype));
9415c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
9425c6c1daeSBarry Smith   }
9435c6c1daeSBarry Smith #endif
9443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9455c6c1daeSBarry Smith }
9465c6c1daeSBarry Smith 
9475c6c1daeSBarry Smith /*@C
948811af0c4SBarry Smith    PetscViewerBinaryWrite - writes to a binary file, only from the first MPI rank
9495c6c1daeSBarry Smith 
950d083f849SBarry Smith    Collective
9515c6c1daeSBarry Smith 
9525c6c1daeSBarry Smith    Input Parameters:
9535c6c1daeSBarry Smith +  viewer - the binary viewer
9545c6c1daeSBarry Smith .  data - location of data
9555824c9d0SPatrick Sanan .  count - number of items of data to write
956f253e43cSLisandro Dalcin -  dtype - type of data to write
9575c6c1daeSBarry Smith 
9585c6c1daeSBarry Smith    Level: beginner
9595c6c1daeSBarry Smith 
960d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
961db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`, `PetscDataType`
962db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
9635c6c1daeSBarry Smith @*/
964d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryWrite(PetscViewer viewer, const void *data, PetscInt count, PetscDataType dtype)
965d71ae5a4SJacob Faibussowitsch {
96622a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
9675c6c1daeSBarry Smith 
9685c6c1daeSBarry Smith   PetscFunctionBegin;
96922a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
97022a8f86cSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, count, 3);
9719566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
97222a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
9735c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
974bc196f7cSDave May   if (vbinary->usempiio) {
9759566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWriteReadMPIIO(viewer, (void *)data, count, NULL, dtype, PETSC_TRUE));
9765c6c1daeSBarry Smith   } else {
9775c6c1daeSBarry Smith #endif
9789566063dSJacob Faibussowitsch     PetscCall(PetscBinarySynchronizedWrite(PetscObjectComm((PetscObject)viewer), vbinary->fdes, data, count, dtype));
9795c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
9805c6c1daeSBarry Smith   }
9815c6c1daeSBarry Smith #endif
9823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9835c6c1daeSBarry Smith }
9845c6c1daeSBarry Smith 
985d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryWriteReadAll(PetscViewer viewer, PetscBool write, void *data, PetscInt count, PetscInt start, PetscInt total, PetscDataType dtype)
986d71ae5a4SJacob Faibussowitsch {
9875972f5f3SLisandro Dalcin   MPI_Comm              comm = PetscObjectComm((PetscObject)viewer);
9885972f5f3SLisandro Dalcin   PetscMPIInt           size, rank;
9895972f5f3SLisandro Dalcin   MPI_Datatype          mdtype;
990ec4bef21SJose E. Roman   PETSC_UNUSED MPI_Aint lb;
991ec4bef21SJose E. Roman   MPI_Aint              dsize;
9925972f5f3SLisandro Dalcin   PetscBool             useMPIIO;
9935972f5f3SLisandro Dalcin 
9945972f5f3SLisandro Dalcin   PetscFunctionBegin;
9955972f5f3SLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
996064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveBool(viewer, ((start >= 0) || (start == PETSC_DETERMINE)), 5);
997064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveBool(viewer, ((total >= 0) || (total == PETSC_DETERMINE)), 6);
998064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveInt(viewer, total, 6);
9999566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
10005972f5f3SLisandro Dalcin 
10019566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToMPIDataType(dtype, &mdtype));
10029566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_get_extent(mdtype, &lb, &dsize));
10039566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
10049566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
10055972f5f3SLisandro Dalcin 
10069566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetUseMPIIO(viewer, &useMPIIO));
10075972f5f3SLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
10085972f5f3SLisandro Dalcin   if (useMPIIO) {
10095972f5f3SLisandro Dalcin     MPI_File    mfdes;
10105972f5f3SLisandro Dalcin     MPI_Offset  off;
10115972f5f3SLisandro Dalcin     PetscMPIInt cnt;
10125972f5f3SLisandro Dalcin 
10135972f5f3SLisandro Dalcin     if (start == PETSC_DETERMINE) {
10149566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Scan(&count, &start, 1, MPIU_INT, MPI_SUM, comm));
10155972f5f3SLisandro Dalcin       start -= count;
10165972f5f3SLisandro Dalcin     }
10175972f5f3SLisandro Dalcin     if (total == PETSC_DETERMINE) {
10185972f5f3SLisandro Dalcin       total = start + count;
10199566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Bcast(&total, 1, MPIU_INT, size - 1, comm));
10205972f5f3SLisandro Dalcin     }
10219566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(count, &cnt));
10229566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryGetMPIIODescriptor(viewer, &mfdes));
10239566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryGetMPIIOOffset(viewer, &off));
10245972f5f3SLisandro Dalcin     off += (MPI_Offset)(start * dsize);
10255972f5f3SLisandro Dalcin     if (write) {
10269566063dSJacob Faibussowitsch       PetscCall(MPIU_File_write_at_all(mfdes, off, data, cnt, mdtype, MPI_STATUS_IGNORE));
10275972f5f3SLisandro Dalcin     } else {
10289566063dSJacob Faibussowitsch       PetscCall(MPIU_File_read_at_all(mfdes, off, data, cnt, mdtype, MPI_STATUS_IGNORE));
10295972f5f3SLisandro Dalcin     }
10305972f5f3SLisandro Dalcin     off = (MPI_Offset)(total * dsize);
10319566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryAddMPIIOOffset(viewer, off));
10323ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
10335972f5f3SLisandro Dalcin   }
10345972f5f3SLisandro Dalcin #endif
10355972f5f3SLisandro Dalcin   {
10365972f5f3SLisandro Dalcin     int         fdes;
10375972f5f3SLisandro Dalcin     char       *workbuf = NULL;
1038dd400576SPatrick Sanan     PetscInt    tcount = rank == 0 ? 0 : count, maxcount = 0, message_count, flowcontrolcount;
10395972f5f3SLisandro Dalcin     PetscMPIInt tag, cnt, maxcnt, scnt = 0, rcnt = 0, j;
10405972f5f3SLisandro Dalcin     MPI_Status  status;
10415972f5f3SLisandro Dalcin 
10429566063dSJacob Faibussowitsch     PetscCall(PetscCommGetNewTag(comm, &tag));
10439566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Reduce(&tcount, &maxcount, 1, MPIU_INT, MPI_MAX, 0, comm));
10449566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(maxcount, &maxcnt));
10459566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(count, &cnt));
10465972f5f3SLisandro Dalcin 
10479566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryGetDescriptor(viewer, &fdes));
10489566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlowControlStart(viewer, &message_count, &flowcontrolcount));
1049dd400576SPatrick Sanan     if (rank == 0) {
10509566063dSJacob Faibussowitsch       PetscCall(PetscMalloc(maxcnt * dsize, &workbuf));
10515972f5f3SLisandro Dalcin       if (write) {
10529566063dSJacob Faibussowitsch         PetscCall(PetscBinaryWrite(fdes, data, cnt, dtype));
10535972f5f3SLisandro Dalcin       } else {
10549566063dSJacob Faibussowitsch         PetscCall(PetscBinaryRead(fdes, data, cnt, NULL, dtype));
10555972f5f3SLisandro Dalcin       }
10565972f5f3SLisandro Dalcin       for (j = 1; j < size; j++) {
10579566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlowControlStepMain(viewer, j, &message_count, flowcontrolcount));
10585972f5f3SLisandro Dalcin         if (write) {
10599566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Recv(workbuf, maxcnt, mdtype, j, tag, comm, &status));
10609566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Get_count(&status, mdtype, &rcnt));
10619566063dSJacob Faibussowitsch           PetscCall(PetscBinaryWrite(fdes, workbuf, rcnt, dtype));
10625972f5f3SLisandro Dalcin         } else {
10639566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Recv(&scnt, 1, MPI_INT, j, tag, comm, MPI_STATUS_IGNORE));
10649566063dSJacob Faibussowitsch           PetscCall(PetscBinaryRead(fdes, workbuf, scnt, NULL, dtype));
10659566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Send(workbuf, scnt, mdtype, j, tag, comm));
10665972f5f3SLisandro Dalcin         }
10675972f5f3SLisandro Dalcin       }
10689566063dSJacob Faibussowitsch       PetscCall(PetscFree(workbuf));
10699566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlowControlEndMain(viewer, &message_count));
10705972f5f3SLisandro Dalcin     } else {
10719566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlowControlStepWorker(viewer, rank, &message_count));
10725972f5f3SLisandro Dalcin       if (write) {
10739566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(data, cnt, mdtype, 0, tag, comm));
10745972f5f3SLisandro Dalcin       } else {
10759566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(&cnt, 1, MPI_INT, 0, tag, comm));
10769566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Recv(data, cnt, mdtype, 0, tag, comm, MPI_STATUS_IGNORE));
10775972f5f3SLisandro Dalcin       }
10789566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlowControlEndWorker(viewer, &message_count));
10795972f5f3SLisandro Dalcin     }
10805972f5f3SLisandro Dalcin   }
10813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10825972f5f3SLisandro Dalcin }
10835972f5f3SLisandro Dalcin 
10845972f5f3SLisandro Dalcin /*@C
1085811af0c4SBarry Smith    PetscViewerBinaryReadAll - reads from a binary file from all MPI ranks, each rank receives its own portion of the data
10865972f5f3SLisandro Dalcin 
10875972f5f3SLisandro Dalcin    Collective
10885972f5f3SLisandro Dalcin 
10895972f5f3SLisandro Dalcin    Input Parameters:
10905972f5f3SLisandro Dalcin +  viewer - the binary viewer
10915972f5f3SLisandro Dalcin .  data - location of data
10925972f5f3SLisandro Dalcin .  count - local number of items of data to read
1093811af0c4SBarry Smith .  start - local start, can be `PETSC_DETERMINE`
1094811af0c4SBarry Smith .  total - global number of items of data to read, can be `PETSC_DETERMINE`
10955972f5f3SLisandro Dalcin -  dtype - type of data to read
10965972f5f3SLisandro Dalcin 
10975972f5f3SLisandro Dalcin    Level: advanced
10985972f5f3SLisandro Dalcin 
1099d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryRead()`, `PetscViewerBinaryWriteAll()`
11005972f5f3SLisandro Dalcin @*/
1101d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryReadAll(PetscViewer viewer, void *data, PetscInt count, PetscInt start, PetscInt total, PetscDataType dtype)
1102d71ae5a4SJacob Faibussowitsch {
11035972f5f3SLisandro Dalcin   PetscFunctionBegin;
11049566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWriteReadAll(viewer, PETSC_FALSE, data, count, start, total, dtype));
11053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11065972f5f3SLisandro Dalcin }
11075972f5f3SLisandro Dalcin 
11085972f5f3SLisandro Dalcin /*@C
1109811af0c4SBarry Smith    PetscViewerBinaryWriteAll - writes to a binary file from all MPI ranks, each rank writes its own portion of the data
11105972f5f3SLisandro Dalcin 
11115972f5f3SLisandro Dalcin    Collective
11125972f5f3SLisandro Dalcin 
11135972f5f3SLisandro Dalcin    Input Parameters:
11145972f5f3SLisandro Dalcin +  viewer - the binary viewer
11155972f5f3SLisandro Dalcin .  data - location of data
11165972f5f3SLisandro Dalcin .  count - local number of items of data to write
1117811af0c4SBarry Smith .  start - local start, can be `PETSC_DETERMINE`
1118811af0c4SBarry Smith .  total - global number of items of data to write, can be `PETSC_DETERMINE`
11195972f5f3SLisandro Dalcin -  dtype - type of data to write
11205972f5f3SLisandro Dalcin 
11215972f5f3SLisandro Dalcin    Level: advanced
11225972f5f3SLisandro Dalcin 
1123d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryWriteAll()`, `PetscViewerBinaryReadAll()`
11245972f5f3SLisandro Dalcin @*/
1125d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryWriteAll(PetscViewer viewer, const void *data, PetscInt count, PetscInt start, PetscInt total, PetscDataType dtype)
1126d71ae5a4SJacob Faibussowitsch {
11275972f5f3SLisandro Dalcin   PetscFunctionBegin;
11289566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWriteReadAll(viewer, PETSC_TRUE, (void *)data, count, start, total, dtype));
11293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11305972f5f3SLisandro Dalcin }
11315972f5f3SLisandro Dalcin 
11325c6c1daeSBarry Smith /*@C
1133811af0c4SBarry Smith    PetscViewerBinaryWriteStringArray - writes to a binary file, only from the first MPI rank, an array of strings
11345c6c1daeSBarry Smith 
1135d083f849SBarry Smith    Collective
11365c6c1daeSBarry Smith 
11375c6c1daeSBarry Smith    Input Parameters:
11385c6c1daeSBarry Smith +  viewer - the binary viewer
11395c6c1daeSBarry Smith -  data - location of the array of strings
11405c6c1daeSBarry Smith 
11415c6c1daeSBarry Smith    Level: intermediate
11425c6c1daeSBarry Smith 
1143811af0c4SBarry Smith     Note:
11443f423023SBarry Smith     The array of strings must be `NULL` terminated
11455c6c1daeSBarry Smith 
1146d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
1147db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
1148db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
11495c6c1daeSBarry Smith @*/
1150d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryWriteStringArray(PetscViewer viewer, const char *const *data)
1151d71ae5a4SJacob Faibussowitsch {
11525c6c1daeSBarry Smith   PetscInt i, n = 0, *sizes;
1153cc843e7aSLisandro Dalcin   size_t   len;
11545c6c1daeSBarry Smith 
115522a8f86cSLisandro Dalcin   PetscFunctionBegin;
11569566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
11575c6c1daeSBarry Smith   /* count number of strings */
11589371c9d4SSatish Balay   while (data[n++])
11599371c9d4SSatish Balay     ;
11605c6c1daeSBarry Smith   n--;
11619566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n + 1, &sizes));
11625c6c1daeSBarry Smith   sizes[0] = n;
11635c6c1daeSBarry Smith   for (i = 0; i < n; i++) {
11649566063dSJacob Faibussowitsch     PetscCall(PetscStrlen(data[i], &len));
1165cc843e7aSLisandro Dalcin     sizes[i + 1] = (PetscInt)len + 1; /* size includes space for the null terminator */
11665c6c1daeSBarry Smith   }
11679566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWrite(viewer, sizes, n + 1, PETSC_INT));
116848a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall(PetscViewerBinaryWrite(viewer, (void *)data[i], sizes[i + 1], PETSC_CHAR));
11699566063dSJacob Faibussowitsch   PetscCall(PetscFree(sizes));
11703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11715c6c1daeSBarry Smith }
11725c6c1daeSBarry Smith 
11735c6c1daeSBarry Smith /*@C
1174811af0c4SBarry Smith    PetscViewerBinaryReadStringArray - reads a binary file an array of strings to all MPI ranks
11755c6c1daeSBarry Smith 
1176d083f849SBarry Smith    Collective
11775c6c1daeSBarry Smith 
11785c6c1daeSBarry Smith    Input Parameter:
11795c6c1daeSBarry Smith .  viewer - the binary viewer
11805c6c1daeSBarry Smith 
11815c6c1daeSBarry Smith    Output Parameter:
11825c6c1daeSBarry Smith .  data - location of the array of strings
11835c6c1daeSBarry Smith 
11845c6c1daeSBarry Smith    Level: intermediate
11855c6c1daeSBarry Smith 
1186811af0c4SBarry Smith     Note:
11873f423023SBarry Smith     The array of strings must `NULL` terminated
11885c6c1daeSBarry Smith 
1189d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
1190db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
1191db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
11925c6c1daeSBarry Smith @*/
1193d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryReadStringArray(PetscViewer viewer, char ***data)
1194d71ae5a4SJacob Faibussowitsch {
1195060da220SMatthew G. Knepley   PetscInt i, n, *sizes, N = 0;
11965c6c1daeSBarry Smith 
119722a8f86cSLisandro Dalcin   PetscFunctionBegin;
11989566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
11995c6c1daeSBarry Smith   /* count number of strings */
12009566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &n, 1, NULL, PETSC_INT));
12019566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &sizes));
12029566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, sizes, n, NULL, PETSC_INT));
1203a297a907SKarl Rupp   for (i = 0; i < n; i++) N += sizes[i];
12049566063dSJacob Faibussowitsch   PetscCall(PetscMalloc((n + 1) * sizeof(char *) + N * sizeof(char), data));
12055c6c1daeSBarry Smith   (*data)[0] = (char *)((*data) + n + 1);
1206a297a907SKarl Rupp   for (i = 1; i < n; i++) (*data)[i] = (*data)[i - 1] + sizes[i - 1];
12079566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, (*data)[0], N, NULL, PETSC_CHAR));
120802c9f0b5SLisandro Dalcin   (*data)[n] = NULL;
12099566063dSJacob Faibussowitsch   PetscCall(PetscFree(sizes));
12103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12115c6c1daeSBarry Smith }
12125c6c1daeSBarry Smith 
1213cc843e7aSLisandro Dalcin /*@C
1214cc843e7aSLisandro Dalcin      PetscViewerFileSetMode - Sets the open mode of file
1215cc843e7aSLisandro Dalcin 
1216c3339decSBarry Smith     Logically Collective
1217cc843e7aSLisandro Dalcin 
1218cc843e7aSLisandro Dalcin   Input Parameters:
1219811af0c4SBarry Smith +  viewer - the `PetscViewer`; must be a a `PETSCVIEWERBINARY`, `PETSCVIEWERMATLAB`, `PETSCVIEWERHDF5`, or `PETSCVIEWERASCII`  `PetscViewer`
1220cc843e7aSLisandro Dalcin -  mode - open mode of file
1221cf53795eSBarry Smith .vb
1222cf53795eSBarry Smith     FILE_MODE_WRITE - create new file for output
1223cf53795eSBarry Smith     FILE_MODE_READ - open existing file for input
1224cf53795eSBarry Smith     FILE_MODE_APPEND - open existing file for output
1225cf53795eSBarry Smith .ve
1226cc843e7aSLisandro Dalcin 
1227cc843e7aSLisandro Dalcin   Level: advanced
1228cc843e7aSLisandro Dalcin 
1229d1f92df0SBarry Smith .seealso: [](sec_viewers), `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
1230cc843e7aSLisandro Dalcin @*/
1231d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerFileSetMode(PetscViewer viewer, PetscFileMode mode)
1232d71ae5a4SJacob Faibussowitsch {
1233cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1234cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1235cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveEnum(viewer, mode, 2);
123608401ef6SPierre Jolivet   PetscCheck(mode != FILE_MODE_UNDEFINED, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Cannot set FILE_MODE_UNDEFINED");
1237f7d195e4SLawrence Mitchell   PetscCheck(mode >= FILE_MODE_UNDEFINED && mode <= FILE_MODE_APPEND_UPDATE, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_OUTOFRANGE, "Invalid file mode %d", (int)mode);
1238cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerFileSetMode_C", (PetscViewer, PetscFileMode), (viewer, mode));
12393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1240cc843e7aSLisandro Dalcin }
1241cc843e7aSLisandro Dalcin 
1242d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetMode_Binary(PetscViewer viewer, PetscFileMode mode)
1243d71ae5a4SJacob Faibussowitsch {
1244cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1245cc843e7aSLisandro Dalcin 
1246cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1247cc73adaaSBarry Smith   PetscCheck(!viewer->setupcalled || vbinary->filemode == mode, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Cannot change mode to %s after setup", PetscFileModes[mode]);
1248cc843e7aSLisandro Dalcin   vbinary->filemode = mode;
12493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1250cc843e7aSLisandro Dalcin }
1251cc843e7aSLisandro Dalcin 
1252cc843e7aSLisandro Dalcin /*@C
1253cc843e7aSLisandro Dalcin      PetscViewerFileGetMode - Gets the open mode of file
1254cc843e7aSLisandro Dalcin 
1255cc843e7aSLisandro Dalcin     Not Collective
1256cc843e7aSLisandro Dalcin 
1257cc843e7aSLisandro Dalcin   Input Parameter:
1258811af0c4SBarry Smith .  viewer - the `PetscViewer`; must be a `PETSCVIEWERBINARY`, `PETSCVIEWERMATLAB`, `PETSCVIEWERHDF5`, or `PETSCVIEWERASCII`  `PetscViewer`
1259cc843e7aSLisandro Dalcin 
1260cc843e7aSLisandro Dalcin   Output Parameter:
1261cc843e7aSLisandro Dalcin .  mode - open mode of file
1262cf53795eSBarry Smith .vb
1263cf53795eSBarry Smith     FILE_MODE_WRITE - create new file for binary output
1264cf53795eSBarry Smith     FILE_MODE_READ - open existing file for binary input
1265cf53795eSBarry Smith     FILE_MODE_APPEND - open existing file for binary output
1266cf53795eSBarry Smith .ve
1267cc843e7aSLisandro Dalcin 
1268cc843e7aSLisandro Dalcin   Level: advanced
1269cc843e7aSLisandro Dalcin 
1270d1f92df0SBarry Smith .seealso: [](sec_viewers), `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
1271cc843e7aSLisandro Dalcin @*/
1272d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerFileGetMode(PetscViewer viewer, PetscFileMode *mode)
1273d71ae5a4SJacob Faibussowitsch {
1274cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1275cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1276cc843e7aSLisandro Dalcin   PetscValidPointer(mode, 2);
1277cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerFileGetMode_C", (PetscViewer, PetscFileMode *), (viewer, mode));
12783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1279cc843e7aSLisandro Dalcin }
1280cc843e7aSLisandro Dalcin 
1281d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetMode_Binary(PetscViewer viewer, PetscFileMode *mode)
1282d71ae5a4SJacob Faibussowitsch {
1283cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1284cc843e7aSLisandro Dalcin 
1285cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1286cc843e7aSLisandro Dalcin   *mode = vbinary->filemode;
12873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1288cc843e7aSLisandro Dalcin }
1289cc843e7aSLisandro Dalcin 
1290d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetName_Binary(PetscViewer viewer, const char name[])
1291d71ae5a4SJacob Faibussowitsch {
1292cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1293cc843e7aSLisandro Dalcin 
1294cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1295cc843e7aSLisandro Dalcin   if (viewer->setupcalled && vbinary->filename) {
1296cc843e7aSLisandro Dalcin     /* gzip can be run after the file with the previous filename has been closed */
12979566063dSJacob Faibussowitsch     PetscCall(PetscFree(vbinary->ogzfilename));
12989566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(vbinary->filename, &vbinary->ogzfilename));
1299cc843e7aSLisandro Dalcin   }
13009566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary->filename));
13019566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &vbinary->filename));
1302cc843e7aSLisandro Dalcin   viewer->setupcalled = PETSC_FALSE;
13033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1304cc843e7aSLisandro Dalcin }
1305cc843e7aSLisandro Dalcin 
1306d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetName_Binary(PetscViewer viewer, const char **name)
1307d71ae5a4SJacob Faibussowitsch {
13085c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
13095c6c1daeSBarry Smith 
13105c6c1daeSBarry Smith   PetscFunctionBegin;
13115c6c1daeSBarry Smith   *name = vbinary->filename;
13123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13135c6c1daeSBarry Smith }
13145c6c1daeSBarry Smith 
13155c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
1316d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetUp_BinaryMPIIO(PetscViewer viewer)
1317d71ae5a4SJacob Faibussowitsch {
13185c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1319e8a65b7dSLisandro Dalcin   int                 amode;
13205c6c1daeSBarry Smith 
13215c6c1daeSBarry Smith   PetscFunctionBegin;
1322a297a907SKarl Rupp   vbinary->storecompressed = PETSC_FALSE;
13235c6c1daeSBarry Smith 
1324cc843e7aSLisandro Dalcin   vbinary->moff = 0;
1325cc843e7aSLisandro Dalcin   switch (vbinary->filemode) {
1326d71ae5a4SJacob Faibussowitsch   case FILE_MODE_READ:
1327d71ae5a4SJacob Faibussowitsch     amode = MPI_MODE_RDONLY;
1328d71ae5a4SJacob Faibussowitsch     break;
1329d71ae5a4SJacob Faibussowitsch   case FILE_MODE_WRITE:
1330d71ae5a4SJacob Faibussowitsch     amode = MPI_MODE_WRONLY | MPI_MODE_CREATE;
1331d71ae5a4SJacob Faibussowitsch     break;
1332d71ae5a4SJacob Faibussowitsch   case FILE_MODE_APPEND:
1333d71ae5a4SJacob Faibussowitsch     amode = MPI_MODE_WRONLY | MPI_MODE_CREATE | MPI_MODE_APPEND;
1334d71ae5a4SJacob Faibussowitsch     break;
1335d71ae5a4SJacob Faibussowitsch   case FILE_MODE_UNDEFINED:
1336d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerSetUp()");
1337d71ae5a4SJacob Faibussowitsch   default:
1338d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[vbinary->filemode]);
13395c6c1daeSBarry Smith   }
13409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_File_open(PetscObjectComm((PetscObject)viewer), vbinary->filename, amode, MPI_INFO_NULL, &vbinary->mfdes));
134122a8f86cSLisandro Dalcin   /*
134222a8f86cSLisandro Dalcin       The MPI standard does not have MPI_MODE_TRUNCATE. We emulate this behavior by setting the file size to zero.
134322a8f86cSLisandro Dalcin   */
13449566063dSJacob Faibussowitsch   if (vbinary->filemode == FILE_MODE_WRITE) PetscCallMPI(MPI_File_set_size(vbinary->mfdes, 0));
134522a8f86cSLisandro Dalcin   /*
134622a8f86cSLisandro Dalcin       Initially, all processes view the file as a linear byte stream. Therefore, for files opened with MPI_MODE_APPEND,
134722a8f86cSLisandro Dalcin       MPI_File_get_position[_shared](fh, &offset) returns the absolute byte position at the end of file.
134822a8f86cSLisandro Dalcin       Otherwise, we would need to call MPI_File_get_byte_offset(fh, offset, &byte_offset) to convert
134922a8f86cSLisandro Dalcin       the offset in etype units to an absolute byte position.
135022a8f86cSLisandro Dalcin    */
13519566063dSJacob Faibussowitsch   if (vbinary->filemode == FILE_MODE_APPEND) PetscCallMPI(MPI_File_get_position(vbinary->mfdes, &vbinary->moff));
13523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1353cc843e7aSLisandro Dalcin }
1354cc843e7aSLisandro Dalcin #endif
13555c6c1daeSBarry Smith 
1356d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetUp_BinarySTDIO(PetscViewer viewer)
1357d71ae5a4SJacob Faibussowitsch {
1358cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1359cc843e7aSLisandro Dalcin   const char         *fname;
1360bbcf679cSJacob Faibussowitsch   char                bname[PETSC_MAX_PATH_LEN], *gz = NULL;
1361cc843e7aSLisandro Dalcin   PetscBool           found;
1362cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
13635c6c1daeSBarry Smith 
1364cc843e7aSLisandro Dalcin   PetscFunctionBegin;
13659566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
13665c6c1daeSBarry Smith 
1367cc843e7aSLisandro Dalcin   /* if file name ends in .gz strip that off and note user wants file compressed */
1368cc843e7aSLisandro Dalcin   vbinary->storecompressed = PETSC_FALSE;
1369cc843e7aSLisandro Dalcin   if (vbinary->filemode == FILE_MODE_WRITE) {
13709566063dSJacob Faibussowitsch     PetscCall(PetscStrstr(vbinary->filename, ".gz", &gz));
13719371c9d4SSatish Balay     if (gz && gz[3] == 0) {
13729371c9d4SSatish Balay       *gz                      = 0;
13739371c9d4SSatish Balay       vbinary->storecompressed = PETSC_TRUE;
13749371c9d4SSatish Balay     }
1375cc843e7aSLisandro Dalcin   }
1376cc843e7aSLisandro Dalcin #if !defined(PETSC_HAVE_POPEN)
137728b400f6SJacob Faibussowitsch   PetscCheck(!vbinary->storecompressed, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP_SYS, "Cannot run gzip on this machine");
1378cc843e7aSLisandro Dalcin #endif
1379cc843e7aSLisandro Dalcin 
1380cc843e7aSLisandro Dalcin   fname = vbinary->filename;
1381cc843e7aSLisandro Dalcin   if (vbinary->filemode == FILE_MODE_READ) { /* possibly get the file from remote site or compressed file */
13829566063dSJacob Faibussowitsch     PetscCall(PetscFileRetrieve(PetscObjectComm((PetscObject)viewer), fname, bname, PETSC_MAX_PATH_LEN, &found));
138328b400f6SJacob Faibussowitsch     PetscCheck(found, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_OPEN, "Cannot locate file: %s", fname);
1384cc843e7aSLisandro Dalcin     fname = bname;
13855c6c1daeSBarry Smith   }
13865c6c1daeSBarry Smith 
1387cc843e7aSLisandro Dalcin   vbinary->fdes = -1;
1388dd400576SPatrick Sanan   if (rank == 0) { /* only first processor opens file*/
1389cc843e7aSLisandro Dalcin     PetscFileMode mode = vbinary->filemode;
1390cc843e7aSLisandro Dalcin     if (mode == FILE_MODE_APPEND) {
1391cc843e7aSLisandro Dalcin       /* check if asked to append to a non-existing file */
13929566063dSJacob Faibussowitsch       PetscCall(PetscTestFile(fname, '\0', &found));
1393cc843e7aSLisandro Dalcin       if (!found) mode = FILE_MODE_WRITE;
1394cc843e7aSLisandro Dalcin     }
13959566063dSJacob Faibussowitsch     PetscCall(PetscBinaryOpen(fname, mode, &vbinary->fdes));
1396cc843e7aSLisandro Dalcin   }
13973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1398cc843e7aSLisandro Dalcin }
1399cc843e7aSLisandro Dalcin 
1400*bf31d2d3SBarry Smith #include <errno.h>
1401d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetUp_BinaryInfo(PetscViewer viewer)
1402d71ae5a4SJacob Faibussowitsch {
1403cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1404cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
1405cc843e7aSLisandro Dalcin   PetscBool           found;
1406cc843e7aSLisandro Dalcin 
1407cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1408cc843e7aSLisandro Dalcin   vbinary->fdes_info = NULL;
14099566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
1410dd400576SPatrick Sanan   if (!vbinary->skipinfo && (vbinary->filemode == FILE_MODE_READ || rank == 0)) {
1411cc843e7aSLisandro Dalcin     char infoname[PETSC_MAX_PATH_LEN], iname[PETSC_MAX_PATH_LEN], *gz;
1412cc843e7aSLisandro Dalcin 
14139566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(infoname, vbinary->filename, sizeof(infoname)));
1414cc843e7aSLisandro Dalcin     /* remove .gz if it ends file name */
14159566063dSJacob Faibussowitsch     PetscCall(PetscStrstr(infoname, ".gz", &gz));
1416cc843e7aSLisandro Dalcin     if (gz && gz[3] == 0) *gz = 0;
1417cc843e7aSLisandro Dalcin 
14189566063dSJacob Faibussowitsch     PetscCall(PetscStrlcat(infoname, ".info", sizeof(infoname)));
1419cc843e7aSLisandro Dalcin     if (vbinary->filemode == FILE_MODE_READ) {
14209566063dSJacob Faibussowitsch       PetscCall(PetscFixFilename(infoname, iname));
14219566063dSJacob Faibussowitsch       PetscCall(PetscFileRetrieve(PetscObjectComm((PetscObject)viewer), iname, infoname, PETSC_MAX_PATH_LEN, &found));
14229566063dSJacob Faibussowitsch       if (found) PetscCall(PetscOptionsInsertFile(PetscObjectComm((PetscObject)viewer), ((PetscObject)viewer)->options, infoname, PETSC_FALSE));
1423dd400576SPatrick Sanan     } else if (rank == 0) { /* write or append */
1424cc843e7aSLisandro Dalcin       const char *omode  = (vbinary->filemode == FILE_MODE_APPEND) ? "a" : "w";
1425cc843e7aSLisandro Dalcin       vbinary->fdes_info = fopen(infoname, omode);
1426*bf31d2d3SBarry Smith       PetscCheck(vbinary->fdes_info, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open .info file %s for writing due to \"%s\"", infoname, strerror(errno));
14275c6c1daeSBarry Smith     }
14285c6c1daeSBarry Smith   }
14293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14305c6c1daeSBarry Smith }
14315c6c1daeSBarry Smith 
1432d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerSetUp_Binary(PetscViewer viewer)
1433d71ae5a4SJacob Faibussowitsch {
14345c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1435cc843e7aSLisandro Dalcin   PetscBool           usempiio;
1436cc843e7aSLisandro Dalcin 
14375c6c1daeSBarry Smith   PetscFunctionBegin;
14389566063dSJacob Faibussowitsch   if (!vbinary->setfromoptionscalled) PetscCall(PetscViewerSetFromOptions(viewer));
143928b400f6SJacob Faibussowitsch   PetscCheck(vbinary->filename, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetName()");
144008401ef6SPierre Jolivet   PetscCheck(vbinary->filemode != (PetscFileMode)-1, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode()");
14419566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_Binary(viewer));
1442cc843e7aSLisandro Dalcin 
14439566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetUseMPIIO(viewer, &usempiio));
1444cc843e7aSLisandro Dalcin   if (usempiio) {
1445cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
14469566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetUp_BinaryMPIIO(viewer));
1447cc843e7aSLisandro Dalcin #endif
1448cc843e7aSLisandro Dalcin   } else {
14499566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetUp_BinarySTDIO(viewer));
1450cc843e7aSLisandro Dalcin   }
14519566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetUp_BinaryInfo(viewer));
1452cc843e7aSLisandro Dalcin 
14539566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectState((PetscObject)viewer, "File: %s", vbinary->filename));
14543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14555c6c1daeSBarry Smith }
14565c6c1daeSBarry Smith 
1457d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerView_Binary(PetscViewer v, PetscViewer viewer)
1458d71ae5a4SJacob Faibussowitsch {
1459cb6ad94fSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
1460cb6ad94fSLisandro Dalcin   const char         *fname   = vbinary->filename ? vbinary->filename : "not yet set";
1461cc843e7aSLisandro Dalcin   const char         *fmode   = vbinary->filemode != (PetscFileMode)-1 ? PetscFileModes[vbinary->filemode] : "not yet set";
1462cb6ad94fSLisandro Dalcin   PetscBool           usempiio;
14632bf49c77SBarry Smith 
14642bf49c77SBarry Smith   PetscFunctionBegin;
14659566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetUseMPIIO(v, &usempiio));
14669566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Filename: %s\n", fname));
14679566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Mode: %s (%s)\n", fmode, usempiio ? "mpiio" : "stdio"));
14683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14692bf49c77SBarry Smith }
14702bf49c77SBarry Smith 
1471d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerSetFromOptions_Binary(PetscViewer viewer, PetscOptionItems *PetscOptionsObject)
1472d71ae5a4SJacob Faibussowitsch {
147322a8f86cSLisandro Dalcin   PetscViewer_Binary *binary = (PetscViewer_Binary *)viewer->data;
1474d22fd6bcSDave May   char                defaultname[PETSC_MAX_PATH_LEN];
147503a1d59fSDave May   PetscBool           flg;
1476e0385b85SDave May 
147703a1d59fSDave May   PetscFunctionBegin;
14783ba16761SJacob Faibussowitsch   if (viewer->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
1479d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Binary PetscViewer Options");
14809566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(defaultname, PETSC_MAX_PATH_LEN - 1, "binaryoutput"));
14819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-viewer_binary_filename", "Specify filename", "PetscViewerFileSetName", defaultname, defaultname, sizeof(defaultname), &flg));
14829566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscViewerFileSetName_Binary(viewer, defaultname));
14839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_skip_info", "Skip writing/reading .info file", "PetscViewerBinarySetSkipInfo", binary->skipinfo, &binary->skipinfo, NULL));
14849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_skip_options", "Skip parsing Vec/Mat load options", "PetscViewerBinarySetSkipOptions", binary->skipoptions, &binary->skipoptions, NULL));
14859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_skip_header", "Skip writing/reading header information", "PetscViewerBinarySetSkipHeader", binary->skipheader, &binary->skipheader, NULL));
148603a1d59fSDave May #if defined(PETSC_HAVE_MPIIO)
14879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_mpiio", "Use MPI-IO functionality to write/read binary file", "PetscViewerBinarySetUseMPIIO", binary->usempiio, &binary->usempiio, NULL));
14885972f5f3SLisandro Dalcin #else
14899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_mpiio", "Use MPI-IO functionality to write/read binary file (NOT AVAILABLE)", "PetscViewerBinarySetUseMPIIO", PETSC_FALSE, NULL, NULL));
149003a1d59fSDave May #endif
1491d0609cedSBarry Smith   PetscOptionsHeadEnd();
1492bc196f7cSDave May   binary->setfromoptionscalled = PETSC_TRUE;
14933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
149403a1d59fSDave May }
149503a1d59fSDave May 
14968556b5ebSBarry Smith /*MC
14978556b5ebSBarry Smith    PETSCVIEWERBINARY - A viewer that saves to binary files
14988556b5ebSBarry Smith 
14991b266c99SBarry Smith   Level: beginner
15001b266c99SBarry Smith 
1501d1f92df0SBarry Smith .seealso: [](sec_viewers), `PetscViewerBinaryOpen()`, `PETSC_VIEWER_STDOUT_()`, `PETSC_VIEWER_STDOUT_SELF`, `PETSC_VIEWER_STDOUT_WORLD`, `PetscViewerCreate()`, `PetscViewerASCIIOpen()`,
1502811af0c4SBarry Smith           `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`, `PETSCVIEWERDRAW`, `PETSCVIEWERSOCKET`
1503811af0c4SBarry Smith           `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()`,
1504811af0c4SBarry Smith           `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`
15058556b5ebSBarry Smith M*/
15068556b5ebSBarry Smith 
1507d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscViewerCreate_Binary(PetscViewer v)
1508d71ae5a4SJacob Faibussowitsch {
15095c6c1daeSBarry Smith   PetscViewer_Binary *vbinary;
15105c6c1daeSBarry Smith 
15115c6c1daeSBarry Smith   PetscFunctionBegin;
15124dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&vbinary));
15135c6c1daeSBarry Smith   v->data = (void *)vbinary;
1514cc843e7aSLisandro Dalcin 
151503a1d59fSDave May   v->ops->setfromoptions   = PetscViewerSetFromOptions_Binary;
15165c6c1daeSBarry Smith   v->ops->destroy          = PetscViewerDestroy_Binary;
15172bf49c77SBarry Smith   v->ops->view             = PetscViewerView_Binary;
1518c98fd787SBarry Smith   v->ops->setup            = PetscViewerSetUp_Binary;
1519cc843e7aSLisandro Dalcin   v->ops->flush            = NULL; /* Should we support Flush() ? */
1520cc843e7aSLisandro Dalcin   v->ops->getsubviewer     = PetscViewerGetSubViewer_Binary;
1521cc843e7aSLisandro Dalcin   v->ops->restoresubviewer = PetscViewerRestoreSubViewer_Binary;
1522cc843e7aSLisandro Dalcin   v->ops->read             = PetscViewerBinaryRead;
1523cc843e7aSLisandro Dalcin 
1524cc843e7aSLisandro Dalcin   vbinary->fdes = -1;
1525e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
1526cc843e7aSLisandro Dalcin   vbinary->usempiio = PETSC_FALSE;
1527e8a65b7dSLisandro Dalcin   vbinary->mfdes    = MPI_FILE_NULL;
1528e8a65b7dSLisandro Dalcin   vbinary->mfsub    = MPI_FILE_NULL;
1529e8a65b7dSLisandro Dalcin #endif
1530cc843e7aSLisandro Dalcin   vbinary->filename        = NULL;
15317e4fd573SVaclav Hapla   vbinary->filemode        = FILE_MODE_UNDEFINED;
153202c9f0b5SLisandro Dalcin   vbinary->fdes_info       = NULL;
15335c6c1daeSBarry Smith   vbinary->skipinfo        = PETSC_FALSE;
15345c6c1daeSBarry Smith   vbinary->skipoptions     = PETSC_TRUE;
15355c6c1daeSBarry Smith   vbinary->skipheader      = PETSC_FALSE;
15365c6c1daeSBarry Smith   vbinary->storecompressed = PETSC_FALSE;
1537f90597f1SStefano Zampini   vbinary->ogzfilename     = NULL;
15385c6c1daeSBarry Smith   vbinary->flowcontrol     = 256; /* seems a good number for Cray XT-5 */
15395c6c1daeSBarry Smith 
1540cc843e7aSLisandro Dalcin   vbinary->setfromoptionscalled = PETSC_FALSE;
1541cc843e7aSLisandro Dalcin 
15429566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetFlowControl_C", PetscViewerBinaryGetFlowControl_Binary));
15439566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetFlowControl_C", PetscViewerBinarySetFlowControl_Binary));
15449566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipHeader_C", PetscViewerBinaryGetSkipHeader_Binary));
15459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipHeader_C", PetscViewerBinarySetSkipHeader_Binary));
15469566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipOptions_C", PetscViewerBinaryGetSkipOptions_Binary));
15479566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipOptions_C", PetscViewerBinarySetSkipOptions_Binary));
15489566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipInfo_C", PetscViewerBinaryGetSkipInfo_Binary));
15499566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipInfo_C", PetscViewerBinarySetSkipInfo_Binary));
15509566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetInfoPointer_C", PetscViewerBinaryGetInfoPointer_Binary));
15519566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_Binary));
15529566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_Binary));
15539566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_Binary));
15549566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_Binary));
15555c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
15569566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetUseMPIIO_C", PetscViewerBinaryGetUseMPIIO_Binary));
15579566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetUseMPIIO_C", PetscViewerBinarySetUseMPIIO_Binary));
15585c6c1daeSBarry Smith #endif
15593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15605c6c1daeSBarry Smith }
15615c6c1daeSBarry Smith 
15625c6c1daeSBarry Smith /*
15635c6c1daeSBarry Smith     The variable Petsc_Viewer_Binary_keyval is used to indicate an MPI attribute that
15645c6c1daeSBarry Smith   is attached to a communicator, in this case the attribute is a PetscViewer.
15655c6c1daeSBarry Smith */
1566d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_Binary_keyval = MPI_KEYVAL_INVALID;
15675c6c1daeSBarry Smith 
15685c6c1daeSBarry Smith /*@C
1569811af0c4SBarry Smith      PETSC_VIEWER_BINARY_ - Creates a `PETSCVIEWERBINARY` `PetscViewer` shared by all processors
15705c6c1daeSBarry Smith                      in a communicator.
15715c6c1daeSBarry Smith 
1572d083f849SBarry Smith      Collective
15735c6c1daeSBarry Smith 
15745c6c1daeSBarry Smith      Input Parameter:
1575811af0c4SBarry Smith .    comm - the MPI communicator to share the `PETSCVIEWERBINARY`
15765c6c1daeSBarry Smith 
15775c6c1daeSBarry Smith      Level: intermediate
15785c6c1daeSBarry Smith 
15795c6c1daeSBarry Smith    Options Database Keys:
158010699b91SBarry Smith +    -viewer_binary_filename <name> - filename in which to store the binary data, defaults to binaryoutput
158110699b91SBarry Smith .    -viewer_binary_skip_info - true means do not create .info file for this viewer
158210699b91SBarry Smith .    -viewer_binary_skip_options - true means do not use the options database for this viewer
158310699b91SBarry Smith .    -viewer_binary_skip_header - true means do not store the usual header information in the binary file
158410699b91SBarry Smith -    -viewer_binary_mpiio - true means use the file via MPI-IO, maybe faster for large files and many MPI ranks
15855c6c1daeSBarry Smith 
1586811af0c4SBarry Smith    Environmental variable:
158710699b91SBarry Smith -   PETSC_VIEWER_BINARY_FILENAME - filename in which to store the binary data, defaults to binaryoutput
15885c6c1daeSBarry Smith 
1589811af0c4SBarry Smith      Note:
1590811af0c4SBarry Smith      Unlike almost all other PETSc routines, `PETSC_VIEWER_BINARY_` does not return
15915c6c1daeSBarry Smith      an error code.  The binary PetscViewer is usually used in the form
15925c6c1daeSBarry Smith $       XXXView(XXX object,PETSC_VIEWER_BINARY_(comm));
15935c6c1daeSBarry Smith 
1594d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PETSC_VIEWER_BINARY_WORLD`, `PETSC_VIEWER_BINARY_SELF`, `PetscViewerBinaryOpen()`, `PetscViewerCreate()`,
1595db781477SPatrick Sanan           `PetscViewerDestroy()`
15965c6c1daeSBarry Smith @*/
1597d71ae5a4SJacob Faibussowitsch PetscViewer PETSC_VIEWER_BINARY_(MPI_Comm comm)
1598d71ae5a4SJacob Faibussowitsch {
15995c6c1daeSBarry Smith   PetscErrorCode ierr;
16003ba16761SJacob Faibussowitsch   PetscMPIInt    mpi_ierr;
16015c6c1daeSBarry Smith   PetscBool      flg;
16025c6c1daeSBarry Smith   PetscViewer    viewer;
16035c6c1daeSBarry Smith   char           fname[PETSC_MAX_PATH_LEN];
16045c6c1daeSBarry Smith   MPI_Comm       ncomm;
16055c6c1daeSBarry Smith 
16065c6c1daeSBarry Smith   PetscFunctionBegin;
16079371c9d4SSatish Balay   ierr = PetscCommDuplicate(comm, &ncomm, NULL);
16089371c9d4SSatish Balay   if (ierr) {
16093ba16761SJacob Faibussowitsch     ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16109371c9d4SSatish Balay     PetscFunctionReturn(NULL);
16119371c9d4SSatish Balay   }
16125c6c1daeSBarry Smith   if (Petsc_Viewer_Binary_keyval == MPI_KEYVAL_INVALID) {
16133ba16761SJacob Faibussowitsch     mpi_ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Viewer_Binary_keyval, NULL);
16143ba16761SJacob Faibussowitsch     if (mpi_ierr) {
16153ba16761SJacob Faibussowitsch       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16169371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16179371c9d4SSatish Balay     }
16185c6c1daeSBarry Smith   }
16193ba16761SJacob Faibussowitsch   mpi_ierr = MPI_Comm_get_attr(ncomm, Petsc_Viewer_Binary_keyval, (void **)&viewer, (int *)&flg);
16203ba16761SJacob Faibussowitsch   if (mpi_ierr) {
16213ba16761SJacob Faibussowitsch     ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16229371c9d4SSatish Balay     PetscFunctionReturn(NULL);
16239371c9d4SSatish Balay   }
16245c6c1daeSBarry Smith   if (!flg) { /* PetscViewer not yet created */
16255c6c1daeSBarry Smith     ierr = PetscOptionsGetenv(ncomm, "PETSC_VIEWER_BINARY_FILENAME", fname, PETSC_MAX_PATH_LEN, &flg);
16269371c9d4SSatish Balay     if (ierr) {
16273ba16761SJacob Faibussowitsch       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16289371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16299371c9d4SSatish Balay     }
16305c6c1daeSBarry Smith     if (!flg) {
1631c6a7a370SJeremy L Thompson       ierr = PetscStrncpy(fname, "binaryoutput", sizeof(fname));
16329371c9d4SSatish Balay       if (ierr) {
16333ba16761SJacob Faibussowitsch         ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16349371c9d4SSatish Balay         PetscFunctionReturn(NULL);
16359371c9d4SSatish Balay       }
16365c6c1daeSBarry Smith     }
16375c6c1daeSBarry Smith     ierr = PetscViewerBinaryOpen(ncomm, fname, FILE_MODE_WRITE, &viewer);
16389371c9d4SSatish Balay     if (ierr) {
16393ba16761SJacob Faibussowitsch       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16409371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16419371c9d4SSatish Balay     }
16425c6c1daeSBarry Smith     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
16439371c9d4SSatish Balay     if (ierr) {
16443ba16761SJacob Faibussowitsch       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16459371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16469371c9d4SSatish Balay     }
16473ba16761SJacob Faibussowitsch     mpi_ierr = MPI_Comm_set_attr(ncomm, Petsc_Viewer_Binary_keyval, (void *)viewer);
16483ba16761SJacob Faibussowitsch     if (mpi_ierr) {
16493ba16761SJacob Faibussowitsch       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16509371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16519371c9d4SSatish Balay     }
16525c6c1daeSBarry Smith   }
16535c6c1daeSBarry Smith   ierr = PetscCommDestroy(&ncomm);
16549371c9d4SSatish Balay   if (ierr) {
16553ba16761SJacob Faibussowitsch     ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16569371c9d4SSatish Balay     PetscFunctionReturn(NULL);
16579371c9d4SSatish Balay   }
16585c6c1daeSBarry Smith   PetscFunctionReturn(viewer);
16595c6c1daeSBarry Smith }
1660