xref: /petsc/src/sys/classes/viewer/impls/binary/binv.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
1af0996ceSBarry Smith #include <petsc/private/viewerimpl.h> /*I   "petscviewer.h"   I*/
25c6c1daeSBarry Smith 
35c6c1daeSBarry Smith typedef struct {
45c6c1daeSBarry Smith   int fdes; /* file descriptor, ignored if using MPI IO */
55c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
6bc196f7cSDave May   PetscBool  usempiio;
75c6c1daeSBarry Smith   MPI_File   mfdes; /* ignored unless using MPI IO */
8e8a65b7dSLisandro Dalcin   MPI_File   mfsub; /* subviewer support */
95c6c1daeSBarry Smith   MPI_Offset moff;
105c6c1daeSBarry Smith #endif
11cc843e7aSLisandro Dalcin   char         *filename;            /* file name */
12cc843e7aSLisandro Dalcin   PetscFileMode filemode;            /* read/write/append mode */
135c6c1daeSBarry Smith   FILE         *fdes_info;           /* optional file containing info on binary file*/
145c6c1daeSBarry Smith   PetscBool     storecompressed;     /* gzip the write binary file when closing it*/
15f90597f1SStefano Zampini   char         *ogzfilename;         /* gzip can be run after the filename has been updated */
165c6c1daeSBarry Smith   PetscBool     skipinfo;            /* Don't create info file for writing; don't use for reading */
175c6c1daeSBarry Smith   PetscBool     skipoptions;         /* don't use PETSc options database when loading */
185c6c1daeSBarry Smith   PetscInt      flowcontrol;         /* allow only <flowcontrol> messages outstanding at a time while doing IO */
195c6c1daeSBarry Smith   PetscBool     skipheader;          /* don't write header, only raw data */
20a261c58fSBarry Smith   PetscBool     matlabheaderwritten; /* if format is PETSC_VIEWER_BINARY_MATLAB has the MATLAB .info header been written yet */
21c98fd787SBarry Smith   PetscBool     setfromoptionscalled;
225c6c1daeSBarry Smith } PetscViewer_Binary;
235c6c1daeSBarry Smith 
249371c9d4SSatish Balay static PetscErrorCode PetscViewerBinaryClearFunctionList(PetscViewer v) {
252e956fe4SStefano Zampini   PetscFunctionBegin;
262e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetFlowControl_C", NULL));
272e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetFlowControl_C", NULL));
282e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipHeader_C", NULL));
292e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipHeader_C", NULL));
302e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipOptions_C", NULL));
312e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipOptions_C", NULL));
322e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipInfo_C", NULL));
332e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipInfo_C", NULL));
342e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetInfoPointer_C", NULL));
352e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", NULL));
362e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", NULL));
372e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", NULL));
382e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", NULL));
392e956fe4SStefano Zampini #if defined(PETSC_HAVE_MPIIO)
402e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetUseMPIIO_C", NULL));
412e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetUseMPIIO_C", NULL));
422e956fe4SStefano Zampini #endif
432e956fe4SStefano Zampini   PetscFunctionReturn(0);
442e956fe4SStefano Zampini }
452e956fe4SStefano Zampini 
46cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
479371c9d4SSatish Balay static PetscErrorCode PetscViewerBinarySyncMPIIO(PetscViewer viewer) {
48cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
49cc843e7aSLisandro Dalcin 
50cc843e7aSLisandro Dalcin   PetscFunctionBegin;
51cc843e7aSLisandro Dalcin   if (vbinary->filemode == FILE_MODE_READ) PetscFunctionReturn(0);
5248a46eb9SPierre Jolivet   if (vbinary->mfsub != MPI_FILE_NULL) PetscCallMPI(MPI_File_sync(vbinary->mfsub));
53cc843e7aSLisandro Dalcin   if (vbinary->mfdes != MPI_FILE_NULL) {
549566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer)));
559566063dSJacob Faibussowitsch     PetscCallMPI(MPI_File_sync(vbinary->mfdes));
56cc843e7aSLisandro Dalcin   }
57cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
58cc843e7aSLisandro Dalcin }
59cc843e7aSLisandro Dalcin #endif
60cc843e7aSLisandro Dalcin 
619371c9d4SSatish Balay static PetscErrorCode PetscViewerGetSubViewer_Binary(PetscViewer viewer, MPI_Comm comm, PetscViewer *outviewer) {
62e8a65b7dSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
63cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
645c6c1daeSBarry Smith 
655c6c1daeSBarry Smith   PetscFunctionBegin;
669566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
67e8a65b7dSLisandro Dalcin 
68e8a65b7dSLisandro Dalcin   /* Return subviewer in process zero */
699566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
70dd400576SPatrick Sanan   if (rank == 0) {
71e5afcf28SBarry Smith     PetscMPIInt flg;
72e5afcf28SBarry Smith 
739566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_compare(PETSC_COMM_SELF, comm, &flg));
74cc73adaaSBarry Smith     PetscCheck(flg == MPI_IDENT || flg == MPI_CONGRUENT, PETSC_COMM_SELF, PETSC_ERR_SUP, "PetscViewerGetSubViewer() for PETSCVIEWERBINARY requires a singleton MPI_Comm");
759566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(comm, outviewer));
769566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(*outviewer, PETSCVIEWERBINARY));
779566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy((*outviewer)->data, vbinary, sizeof(PetscViewer_Binary)));
7812f4c3a9SLisandro Dalcin     (*outviewer)->setupcalled = PETSC_TRUE;
7996509d17SLisandro Dalcin   } else {
8096509d17SLisandro Dalcin     *outviewer = NULL;
8196509d17SLisandro Dalcin   }
82e8a65b7dSLisandro Dalcin 
83e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
84e8a65b7dSLisandro Dalcin   if (vbinary->usempiio && *outviewer) {
85e8a65b7dSLisandro Dalcin     PetscViewer_Binary *obinary = (PetscViewer_Binary *)(*outviewer)->data;
86e8a65b7dSLisandro Dalcin     /* Parent viewer opens a new MPI file handle on PETSC_COMM_SELF and keeps track of it for future reuse */
87e8a65b7dSLisandro Dalcin     if (vbinary->mfsub == MPI_FILE_NULL) {
88e8a65b7dSLisandro Dalcin       int amode;
89cc843e7aSLisandro Dalcin       switch (vbinary->filemode) {
90e8a65b7dSLisandro Dalcin       case FILE_MODE_READ: amode = MPI_MODE_RDONLY; break;
91e8a65b7dSLisandro Dalcin       case FILE_MODE_WRITE: amode = MPI_MODE_WRONLY; break;
9222a8f86cSLisandro Dalcin       case FILE_MODE_APPEND: amode = MPI_MODE_WRONLY; break;
9398921bdaSJacob Faibussowitsch       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[vbinary->filemode]);
94e8a65b7dSLisandro Dalcin       }
959566063dSJacob Faibussowitsch       PetscCallMPI(MPI_File_open(PETSC_COMM_SELF, vbinary->filename, amode, MPI_INFO_NULL, &vbinary->mfsub));
96e8a65b7dSLisandro Dalcin     }
97e8a65b7dSLisandro Dalcin     /* Subviewer gets the MPI file handle on PETSC_COMM_SELF */
98e8a65b7dSLisandro Dalcin     obinary->mfdes = vbinary->mfsub;
99e8a65b7dSLisandro Dalcin     obinary->mfsub = MPI_FILE_NULL;
100e8a65b7dSLisandro Dalcin     obinary->moff  = vbinary->moff;
101e8a65b7dSLisandro Dalcin   }
102e8a65b7dSLisandro Dalcin #endif
103cc843e7aSLisandro Dalcin 
104cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
1059566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinarySyncMPIIO(viewer));
106cc843e7aSLisandro Dalcin #endif
1075c6c1daeSBarry Smith   PetscFunctionReturn(0);
1085c6c1daeSBarry Smith }
1095c6c1daeSBarry Smith 
1109371c9d4SSatish Balay static PetscErrorCode PetscViewerRestoreSubViewer_Binary(PetscViewer viewer, MPI_Comm comm, PetscViewer *outviewer) {
111e8a65b7dSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
112cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
113e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
114e8a65b7dSLisandro Dalcin   MPI_Offset moff = 0;
115e8a65b7dSLisandro Dalcin #endif
1165c6c1daeSBarry Smith 
1175c6c1daeSBarry Smith   PetscFunctionBegin;
1189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
119c5853193SPierre Jolivet   PetscCheck(rank == 0 || !*outviewer, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Subviewer not obtained from viewer");
120e8a65b7dSLisandro Dalcin 
121e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
122e8a65b7dSLisandro Dalcin   if (vbinary->usempiio && *outviewer) {
123e8a65b7dSLisandro Dalcin     PetscViewer_Binary *obinary = (PetscViewer_Binary *)(*outviewer)->data;
12408401ef6SPierre Jolivet     PetscCheck(obinary->mfdes == vbinary->mfsub, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Subviewer not obtained from viewer");
1259566063dSJacob Faibussowitsch     if (obinary->mfsub != MPI_FILE_NULL) PetscCallMPI(MPI_File_close(&obinary->mfsub));
126e8a65b7dSLisandro Dalcin     moff = obinary->moff;
127e8a65b7dSLisandro Dalcin   }
128e8a65b7dSLisandro Dalcin #endif
129e8a65b7dSLisandro Dalcin 
130e8a65b7dSLisandro Dalcin   if (*outviewer) {
131e8a65b7dSLisandro Dalcin     PetscViewer_Binary *obinary = (PetscViewer_Binary *)(*outviewer)->data;
13208401ef6SPierre Jolivet     PetscCheck(obinary->fdes == vbinary->fdes, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Subviewer not obtained from viewer");
1339566063dSJacob Faibussowitsch     PetscCall(PetscFree((*outviewer)->data));
1342e956fe4SStefano Zampini     PetscCall(PetscViewerBinaryClearFunctionList(*outviewer));
1359566063dSJacob Faibussowitsch     PetscCall(PetscHeaderDestroy(outviewer));
1365c6c1daeSBarry Smith   }
137e8a65b7dSLisandro Dalcin 
138e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
139e8a65b7dSLisandro Dalcin   if (vbinary->usempiio) {
140e8a65b7dSLisandro Dalcin     PetscInt64 ioff = (PetscInt64)moff; /* We could use MPI_OFFSET datatype (requires MPI 2.2) */
1419566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&ioff, 1, MPIU_INT64, 0, PetscObjectComm((PetscObject)viewer)));
142e8a65b7dSLisandro Dalcin     vbinary->moff = (MPI_Offset)ioff;
143e8a65b7dSLisandro Dalcin   }
144e8a65b7dSLisandro Dalcin #endif
145cc843e7aSLisandro Dalcin 
146cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
1479566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinarySyncMPIIO(viewer));
148cc843e7aSLisandro Dalcin #endif
1495c6c1daeSBarry Smith   PetscFunctionReturn(0);
1505c6c1daeSBarry Smith }
1515c6c1daeSBarry Smith 
1525c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
1535c6c1daeSBarry Smith /*@C
15485b8072bSPatrick Sanan     PetscViewerBinaryGetMPIIOOffset - Gets the current global offset that should be passed to `MPI_File_set_view()` or `MPI_File_{write|read}_at[_all]()`
1555c6c1daeSBarry Smith 
1565c6c1daeSBarry Smith     Not Collective
1575c6c1daeSBarry Smith 
1585c6c1daeSBarry Smith     Input Parameter:
15985b8072bSPatrick Sanan .   viewer - PetscViewer context, obtained from `PetscViewerBinaryOpen()`
1605c6c1daeSBarry Smith 
1615c6c1daeSBarry Smith     Output Parameter:
16222a8f86cSLisandro Dalcin .   off - the current global offset
1635c6c1daeSBarry Smith 
1645c6c1daeSBarry Smith     Level: advanced
1655c6c1daeSBarry Smith 
166811af0c4SBarry Smith     Note:
167811af0c4SBarry Smith     Use `PetscViewerBinaryAddMPIIOOffset()` to increase this value after you have written a view.
168811af0c4SBarry Smith 
1695c6c1daeSBarry Smith     Fortran Note:
1705c6c1daeSBarry Smith     This routine is not supported in Fortran.
1715c6c1daeSBarry Smith 
172811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryAddMPIIOOffset()`
1735c6c1daeSBarry Smith @*/
1749371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryGetMPIIOOffset(PetscViewer viewer, MPI_Offset *off) {
17522a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
1765c6c1daeSBarry Smith 
1775c6c1daeSBarry Smith   PetscFunctionBegin;
17822a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
17922a8f86cSLisandro Dalcin   PetscValidPointer(off, 2);
18022a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
1815c6c1daeSBarry Smith   *off    = vbinary->moff;
1825c6c1daeSBarry Smith   PetscFunctionReturn(0);
1835c6c1daeSBarry Smith }
1845c6c1daeSBarry Smith 
1855c6c1daeSBarry Smith /*@C
18622a8f86cSLisandro Dalcin     PetscViewerBinaryAddMPIIOOffset - Adds to the current global offset
1875c6c1daeSBarry Smith 
18822a8f86cSLisandro Dalcin     Logically Collective
1895c6c1daeSBarry Smith 
1905c6c1daeSBarry Smith     Input Parameters:
191811af0c4SBarry Smith +   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
19222a8f86cSLisandro Dalcin -   off - the addition to the global offset
1935c6c1daeSBarry Smith 
1945c6c1daeSBarry Smith     Level: advanced
1955c6c1daeSBarry Smith 
196811af0c4SBarry Smith     Note:
197811af0c4SBarry Smith     Use `PetscViewerBinaryGetMPIIOOffset()` to get the value that you should pass to `MPI_File_set_view()` or `MPI_File_{write|read}_at[_all]()`
198811af0c4SBarry Smith 
1995c6c1daeSBarry Smith     Fortran Note:
2005c6c1daeSBarry Smith     This routine is not supported in Fortran.
2015c6c1daeSBarry Smith 
202811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
2035c6c1daeSBarry Smith @*/
2049371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryAddMPIIOOffset(PetscViewer viewer, MPI_Offset off) {
20522a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
2065c6c1daeSBarry Smith 
2075c6c1daeSBarry Smith   PetscFunctionBegin;
20822a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
20922a8f86cSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, (PetscInt)off, 2);
21022a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
2115c6c1daeSBarry Smith   vbinary->moff += off;
2125c6c1daeSBarry Smith   PetscFunctionReturn(0);
2135c6c1daeSBarry Smith }
2145c6c1daeSBarry Smith 
2155c6c1daeSBarry Smith /*@C
216811af0c4SBarry Smith     PetscViewerBinaryGetMPIIODescriptor - Extracts the MPI IO file descriptor from a `PetscViewer`.
2175c6c1daeSBarry Smith 
2185c6c1daeSBarry Smith     Not Collective
2195c6c1daeSBarry Smith 
2205c6c1daeSBarry Smith     Input Parameter:
221811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
2225c6c1daeSBarry Smith 
2235c6c1daeSBarry Smith     Output Parameter:
2245c6c1daeSBarry Smith .   fdes - file descriptor
2255c6c1daeSBarry Smith 
2265c6c1daeSBarry Smith     Level: advanced
2275c6c1daeSBarry Smith 
2285c6c1daeSBarry Smith     Fortran Note:
2295c6c1daeSBarry Smith     This routine is not supported in Fortran.
2305c6c1daeSBarry Smith 
231811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
2325c6c1daeSBarry Smith @*/
2339371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryGetMPIIODescriptor(PetscViewer viewer, MPI_File *fdes) {
23422a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
2355c6c1daeSBarry Smith 
2365c6c1daeSBarry Smith   PetscFunctionBegin;
23722a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
23822a8f86cSLisandro Dalcin   PetscValidPointer(fdes, 2);
2399566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
24022a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
2415c6c1daeSBarry Smith   *fdes   = vbinary->mfdes;
2425c6c1daeSBarry Smith   PetscFunctionReturn(0);
2435c6c1daeSBarry Smith }
244cc843e7aSLisandro Dalcin #endif
2455c6c1daeSBarry Smith 
246cc843e7aSLisandro Dalcin /*@
247cc843e7aSLisandro Dalcin     PetscViewerBinarySetUseMPIIO - Sets a binary viewer to use MPI-IO for reading/writing. Must be called
248811af0c4SBarry Smith         before `PetscViewerFileSetName()`
249cc843e7aSLisandro Dalcin 
250811af0c4SBarry Smith     Logically Collective on viewer
251cc843e7aSLisandro Dalcin 
252cc843e7aSLisandro Dalcin     Input Parameters:
253811af0c4SBarry Smith +   viewer - the `PetscViewer`; must be a `PETSCVIEWERBINARY`
254811af0c4SBarry Smith -   use - `PETSC_TRUE` means MPI-IO will be used
255cc843e7aSLisandro Dalcin 
256811af0c4SBarry Smith     Options Database Key:
257cc843e7aSLisandro Dalcin     -viewer_binary_mpiio : Flag for using MPI-IO
258cc843e7aSLisandro Dalcin 
259cc843e7aSLisandro Dalcin     Level: advanced
260cc843e7aSLisandro Dalcin 
261811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`,
262db781477SPatrick Sanan           `PetscViewerBinaryGetUseMPIIO()`
263cc843e7aSLisandro Dalcin @*/
2649371c9d4SSatish Balay PetscErrorCode PetscViewerBinarySetUseMPIIO(PetscViewer viewer, PetscBool use) {
265a8aae444SDave May   PetscFunctionBegin;
266cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
267cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, use, 2);
268cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetUseMPIIO_C", (PetscViewer, PetscBool), (viewer, use));
269cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
270cc843e7aSLisandro Dalcin }
271cc843e7aSLisandro Dalcin 
272cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
2739371c9d4SSatish Balay static PetscErrorCode PetscViewerBinarySetUseMPIIO_Binary(PetscViewer viewer, PetscBool use) {
274cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
275cc843e7aSLisandro Dalcin   PetscFunctionBegin;
276cc73adaaSBarry Smith   PetscCheck(!viewer->setupcalled || vbinary->usempiio == use, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Cannot change MPIIO to %s after setup", PetscBools[use]);
277cc843e7aSLisandro Dalcin   vbinary->usempiio = use;
278a8aae444SDave May   PetscFunctionReturn(0);
279a8aae444SDave May }
280a8aae444SDave May #endif
281a8aae444SDave May 
282cc843e7aSLisandro Dalcin /*@
283bc196f7cSDave May     PetscViewerBinaryGetUseMPIIO - Returns PETSC_TRUE if the binary viewer uses MPI-IO.
2845c6c1daeSBarry Smith 
2855c6c1daeSBarry Smith     Not Collective
2865c6c1daeSBarry Smith 
2875c6c1daeSBarry Smith     Input Parameter:
288811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`; must be a `PETSCVIEWERBINARY`
2895c6c1daeSBarry Smith 
2905c6c1daeSBarry Smith     Output Parameter:
291811af0c4SBarry Smith .   use - `PETSC_TRUE` if MPI-IO is being used
2925c6c1daeSBarry Smith 
2935c6c1daeSBarry Smith     Level: advanced
2945c6c1daeSBarry Smith 
295bc196f7cSDave May     Note:
296bc196f7cSDave May     If MPI-IO is not available, this function will always return PETSC_FALSE
297bc196f7cSDave May 
298811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
2995c6c1daeSBarry Smith @*/
3009371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryGetUseMPIIO(PetscViewer viewer, PetscBool *use) {
3015c6c1daeSBarry Smith   PetscFunctionBegin;
302cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
303cc843e7aSLisandro Dalcin   PetscValidBoolPointer(use, 2);
304cc843e7aSLisandro Dalcin   *use = PETSC_FALSE;
305cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinaryGetUseMPIIO_C", (PetscViewer, PetscBool *), (viewer, use));
3065c6c1daeSBarry Smith   PetscFunctionReturn(0);
3075c6c1daeSBarry Smith }
3085c6c1daeSBarry Smith 
309cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
3109371c9d4SSatish Balay static PetscErrorCode PetscViewerBinaryGetUseMPIIO_Binary(PetscViewer viewer, PetscBool *use) {
3115c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
3125c6c1daeSBarry Smith 
3135c6c1daeSBarry Smith   PetscFunctionBegin;
314cc843e7aSLisandro Dalcin   *use = vbinary->usempiio;
315cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
316cc843e7aSLisandro Dalcin }
317cc843e7aSLisandro Dalcin #endif
318cc843e7aSLisandro Dalcin 
319cc843e7aSLisandro Dalcin /*@
320cc843e7aSLisandro Dalcin     PetscViewerBinarySetFlowControl - Sets how many messages are allowed to outstanding at the same time during parallel IO reads/writes
321cc843e7aSLisandro Dalcin 
322cc843e7aSLisandro Dalcin     Not Collective
323cc843e7aSLisandro Dalcin 
324d8d19677SJose E. Roman     Input Parameters:
325811af0c4SBarry Smith +   viewer - PetscViewer context, obtained from `PetscViewerBinaryOpen()`
326cc843e7aSLisandro Dalcin -   fc - the number of messages, defaults to 256 if this function was not called
327cc843e7aSLisandro Dalcin 
328cc843e7aSLisandro Dalcin     Level: advanced
329cc843e7aSLisandro Dalcin 
330811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetFlowControl()`
331cc843e7aSLisandro Dalcin @*/
3329371c9d4SSatish Balay PetscErrorCode PetscViewerBinarySetFlowControl(PetscViewer viewer, PetscInt fc) {
333cc843e7aSLisandro Dalcin   PetscFunctionBegin;
334cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
335cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, fc, 2);
336cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetFlowControl_C", (PetscViewer, PetscInt), (viewer, fc));
3375c6c1daeSBarry Smith   PetscFunctionReturn(0);
3385c6c1daeSBarry Smith }
3395c6c1daeSBarry Smith 
3409371c9d4SSatish Balay static PetscErrorCode PetscViewerBinarySetFlowControl_Binary(PetscViewer viewer, PetscInt fc) {
341cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
342cc843e7aSLisandro Dalcin 
343cc843e7aSLisandro Dalcin   PetscFunctionBegin;
34408401ef6SPierre Jolivet   PetscCheck(fc > 1, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_OUTOFRANGE, "Flow control count must be greater than 1, %" PetscInt_FMT " was set", fc);
345cc843e7aSLisandro Dalcin   vbinary->flowcontrol = fc;
346cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
347cc843e7aSLisandro Dalcin }
348cc843e7aSLisandro Dalcin 
349cc843e7aSLisandro Dalcin /*@
3505c6c1daeSBarry Smith     PetscViewerBinaryGetFlowControl - Returns how many messages are allowed to outstanding at the same time during parallel IO reads/writes
3515c6c1daeSBarry Smith 
3525c6c1daeSBarry Smith     Not Collective
3535c6c1daeSBarry Smith 
3545c6c1daeSBarry Smith     Input Parameter:
355811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
3565c6c1daeSBarry Smith 
3575c6c1daeSBarry Smith     Output Parameter:
3585c6c1daeSBarry Smith .   fc - the number of messages
3595c6c1daeSBarry Smith 
3605c6c1daeSBarry Smith     Level: advanced
3615c6c1daeSBarry Smith 
362811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinarySetFlowControl()`
3635c6c1daeSBarry Smith @*/
3649371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryGetFlowControl(PetscViewer viewer, PetscInt *fc) {
3655c6c1daeSBarry Smith   PetscFunctionBegin;
366cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
367cc843e7aSLisandro Dalcin   PetscValidIntPointer(fc, 2);
368cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetFlowControl_C", (PetscViewer, PetscInt *), (viewer, fc));
3695c6c1daeSBarry Smith   PetscFunctionReturn(0);
3705c6c1daeSBarry Smith }
3715c6c1daeSBarry Smith 
3729371c9d4SSatish Balay static PetscErrorCode PetscViewerBinaryGetFlowControl_Binary(PetscViewer viewer, PetscInt *fc) {
3735c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
3745c6c1daeSBarry Smith 
3755c6c1daeSBarry Smith   PetscFunctionBegin;
376cc843e7aSLisandro Dalcin   *fc = vbinary->flowcontrol;
3775c6c1daeSBarry Smith   PetscFunctionReturn(0);
3785c6c1daeSBarry Smith }
3795c6c1daeSBarry Smith 
3805c6c1daeSBarry Smith /*@C
381811af0c4SBarry Smith     PetscViewerBinaryGetDescriptor - Extracts the file descriptor from a `PetscViewer` of `PetscViewerType` `PETSCVIEWERBINARY`.
3825c6c1daeSBarry Smith 
383811af0c4SBarry Smith     Collective on viewer because it may trigger a `PetscViewerSetUp()` call
3845c6c1daeSBarry Smith 
3855c6c1daeSBarry Smith     Input Parameter:
386811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
3875c6c1daeSBarry Smith 
3885c6c1daeSBarry Smith     Output Parameter:
3895c6c1daeSBarry Smith .   fdes - file descriptor
3905c6c1daeSBarry Smith 
3915c6c1daeSBarry Smith     Level: advanced
3925c6c1daeSBarry Smith 
393811af0c4SBarry Smith     Note:
394811af0c4SBarry Smith       For writable binary `PetscViewer`s, the descriptor will only be valid for the
395811af0c4SBarry Smith     first processor in the communicator that shares the `PetscViewer`. For readable
3965c6c1daeSBarry Smith     files it will only be valid on nodes that have the file. If node 0 does not
3975c6c1daeSBarry Smith     have the file it generates an error even if another node does have the file.
3985c6c1daeSBarry Smith 
3995c6c1daeSBarry Smith     Fortran Note:
4005c6c1daeSBarry Smith     This routine is not supported in Fortran.
4015c6c1daeSBarry Smith 
402811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`
4035c6c1daeSBarry Smith @*/
4049371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryGetDescriptor(PetscViewer viewer, int *fdes) {
40522a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
4065c6c1daeSBarry Smith 
4075c6c1daeSBarry Smith   PetscFunctionBegin;
40822a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
40922a8f86cSLisandro Dalcin   PetscValidPointer(fdes, 2);
4109566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
41122a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
4125c6c1daeSBarry Smith   *fdes   = vbinary->fdes;
4135c6c1daeSBarry Smith   PetscFunctionReturn(0);
4145c6c1daeSBarry Smith }
4155c6c1daeSBarry Smith 
4165c6c1daeSBarry Smith /*@
4175c6c1daeSBarry Smith     PetscViewerBinarySkipInfo - Binary file will not have .info file created with it
4185c6c1daeSBarry Smith 
4195c6c1daeSBarry Smith     Not Collective
4205c6c1daeSBarry Smith 
421fd292e60Sprj-     Input Parameter:
422811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerCreate()`
4235c6c1daeSBarry Smith 
4245c6c1daeSBarry Smith     Options Database Key:
42510699b91SBarry Smith .   -viewer_binary_skip_info - true indicates do not generate .info file
4265c6c1daeSBarry Smith 
4275c6c1daeSBarry Smith     Level: advanced
4285c6c1daeSBarry Smith 
42995452b02SPatrick Sanan     Notes:
430811af0c4SBarry Smith     This must be called after `PetscViewerSetType()`. If you use `PetscViewerBinaryOpen()` then
4315c6c1daeSBarry Smith     you can only skip the info file with the -viewer_binary_skip_info flag. To use the function you must open the
432811af0c4SBarry Smith     viewer with `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinarySkipInfo()`.
4335c6c1daeSBarry Smith 
4345c6c1daeSBarry Smith     The .info contains meta information about the data in the binary file, for example the block size if it was
4355c6c1daeSBarry Smith     set for a vector or matrix.
4365c6c1daeSBarry Smith 
437811af0c4SBarry Smith     This routine is deprecated, use `PetscViewerBinarySetSkipInfo()`
438811af0c4SBarry Smith 
439811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySetSkipOptions()`,
440db781477SPatrick Sanan           `PetscViewerBinaryGetSkipOptions()`, `PetscViewerBinaryGetSkipInfo()`
4415c6c1daeSBarry Smith @*/
4429371c9d4SSatish Balay PetscErrorCode PetscViewerBinarySkipInfo(PetscViewer viewer) {
4435c6c1daeSBarry Smith   PetscFunctionBegin;
4449566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinarySetSkipInfo(viewer, PETSC_TRUE));
445807ea322SDave May   PetscFunctionReturn(0);
446807ea322SDave May }
447807ea322SDave May 
448807ea322SDave May /*@
449807ea322SDave May     PetscViewerBinarySetSkipInfo - Binary file will not have .info file created with it
450807ea322SDave May 
451807ea322SDave May     Not Collective
452807ea322SDave May 
453d8d19677SJose E. Roman     Input Parameters:
454cc843e7aSLisandro Dalcin +   viewer - PetscViewer context, obtained from PetscViewerCreate()
455cc843e7aSLisandro Dalcin -   skip - PETSC_TRUE implies the .info file will not be generated
456807ea322SDave May 
457807ea322SDave May     Options Database Key:
45810699b91SBarry Smith .   -viewer_binary_skip_info - true indicates do not generate .info file
459807ea322SDave May 
460807ea322SDave May     Level: advanced
461807ea322SDave May 
462811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySetSkipOptions()`,
463db781477SPatrick Sanan           `PetscViewerBinaryGetSkipOptions()`, `PetscViewerBinaryGetSkipInfo()`
464807ea322SDave May @*/
4659371c9d4SSatish Balay PetscErrorCode PetscViewerBinarySetSkipInfo(PetscViewer viewer, PetscBool skip) {
466807ea322SDave May   PetscFunctionBegin;
467cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
468cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, skip, 2);
469cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetSkipInfo_C", (PetscViewer, PetscBool), (viewer, skip));
470807ea322SDave May   PetscFunctionReturn(0);
471807ea322SDave May }
472807ea322SDave May 
4739371c9d4SSatish Balay static PetscErrorCode PetscViewerBinarySetSkipInfo_Binary(PetscViewer viewer, PetscBool skip) {
474807ea322SDave May   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
475807ea322SDave May 
476807ea322SDave May   PetscFunctionBegin;
477cc843e7aSLisandro Dalcin   vbinary->skipinfo = skip;
478807ea322SDave May   PetscFunctionReturn(0);
479807ea322SDave May }
480807ea322SDave May 
481807ea322SDave May /*@
482807ea322SDave May     PetscViewerBinaryGetSkipInfo - check if viewer wrote a .info file
483807ea322SDave May 
484807ea322SDave May     Not Collective
485807ea322SDave May 
486807ea322SDave May     Input Parameter:
487811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
488807ea322SDave May 
489807ea322SDave May     Output Parameter:
490811af0c4SBarry Smith .   skip - `PETSC_TRUE` implies the .info file was not generated
491807ea322SDave May 
492807ea322SDave May     Level: advanced
493807ea322SDave May 
494811af0c4SBarry Smith     Note:
495811af0c4SBarry Smith     This must be called after `PetscViewerSetType()`
496807ea322SDave May 
497811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`,
498db781477SPatrick Sanan           `PetscViewerBinarySetSkipOptions()`, `PetscViewerBinarySetSkipInfo()`
499807ea322SDave May @*/
5009371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryGetSkipInfo(PetscViewer viewer, PetscBool *skip) {
501807ea322SDave May   PetscFunctionBegin;
502cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
503cc843e7aSLisandro Dalcin   PetscValidBoolPointer(skip, 2);
504cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetSkipInfo_C", (PetscViewer, PetscBool *), (viewer, skip));
505807ea322SDave May   PetscFunctionReturn(0);
506807ea322SDave May }
507807ea322SDave May 
5089371c9d4SSatish Balay static PetscErrorCode PetscViewerBinaryGetSkipInfo_Binary(PetscViewer viewer, PetscBool *skip) {
509807ea322SDave May   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
510807ea322SDave May 
511807ea322SDave May   PetscFunctionBegin;
512cc843e7aSLisandro Dalcin   *skip = vbinary->skipinfo;
513807ea322SDave May   PetscFunctionReturn(0);
514807ea322SDave May }
515807ea322SDave May 
5165c6c1daeSBarry Smith /*@
5175c6c1daeSBarry Smith     PetscViewerBinarySetSkipOptions - do not use the PETSc options database when loading objects
5185c6c1daeSBarry Smith 
5195c6c1daeSBarry Smith     Not Collective
5205c6c1daeSBarry Smith 
5215c6c1daeSBarry Smith     Input Parameters:
522811af0c4SBarry Smith +   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
523811af0c4SBarry Smith -   skip - `PETSC_TRUE` means do not use the options from the options database
5245c6c1daeSBarry Smith 
5255c6c1daeSBarry Smith     Options Database Key:
526811af0c4SBarry Smith .   -viewer_binary_skip_options <true or false> - true means do not use the options from the options database
5275c6c1daeSBarry Smith 
5285c6c1daeSBarry Smith     Level: advanced
5295c6c1daeSBarry Smith 
530811af0c4SBarry Smith     Note:
531811af0c4SBarry Smith     This must be called after `PetscViewerSetType()`
5325c6c1daeSBarry Smith 
533811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
534db781477SPatrick Sanan           `PetscViewerBinaryGetSkipOptions()`
5355c6c1daeSBarry Smith @*/
5369371c9d4SSatish Balay PetscErrorCode PetscViewerBinarySetSkipOptions(PetscViewer viewer, PetscBool skip) {
5375c6c1daeSBarry Smith   PetscFunctionBegin;
538cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
539cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, skip, 2);
540cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetSkipOptions_C", (PetscViewer, PetscBool), (viewer, skip));
541807ea322SDave May   PetscFunctionReturn(0);
542807ea322SDave May }
543807ea322SDave May 
5449371c9d4SSatish Balay static PetscErrorCode PetscViewerBinarySetSkipOptions_Binary(PetscViewer viewer, PetscBool skip) {
545cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
546807ea322SDave May 
547807ea322SDave May   PetscFunctionBegin;
548cc843e7aSLisandro Dalcin   vbinary->skipoptions = skip;
5495c6c1daeSBarry Smith   PetscFunctionReturn(0);
5505c6c1daeSBarry Smith }
5515c6c1daeSBarry Smith 
5525c6c1daeSBarry Smith /*@
5535c6c1daeSBarry Smith     PetscViewerBinaryGetSkipOptions - checks if viewer uses the PETSc options database when loading objects
5545c6c1daeSBarry Smith 
5555c6c1daeSBarry Smith     Not Collective
5565c6c1daeSBarry Smith 
5575c6c1daeSBarry Smith     Input Parameter:
558811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
5595c6c1daeSBarry Smith 
5605c6c1daeSBarry Smith     Output Parameter:
561811af0c4SBarry Smith .   skip - `PETSC_TRUE` means do not use
5625c6c1daeSBarry Smith 
5635c6c1daeSBarry Smith     Level: advanced
5645c6c1daeSBarry Smith 
565811af0c4SBarry Smith     Note:
566811af0c4SBarry Smith     This must be called after `PetscViewerSetType()`
5675c6c1daeSBarry Smith 
568811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
569db781477SPatrick Sanan           `PetscViewerBinarySetSkipOptions()`
5705c6c1daeSBarry Smith @*/
5719371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryGetSkipOptions(PetscViewer viewer, PetscBool *skip) {
5725c6c1daeSBarry Smith   PetscFunctionBegin;
573cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
574cc843e7aSLisandro Dalcin   PetscValidBoolPointer(skip, 2);
575cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetSkipOptions_C", (PetscViewer, PetscBool *), (viewer, skip));
5765c6c1daeSBarry Smith   PetscFunctionReturn(0);
5775c6c1daeSBarry Smith }
5785c6c1daeSBarry Smith 
5799371c9d4SSatish Balay static PetscErrorCode PetscViewerBinaryGetSkipOptions_Binary(PetscViewer viewer, PetscBool *skip) {
580d21b9a37SPierre Jolivet   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
5815c6c1daeSBarry Smith 
5825c6c1daeSBarry Smith   PetscFunctionBegin;
583cc843e7aSLisandro Dalcin   *skip = vbinary->skipoptions;
5845c6c1daeSBarry Smith   PetscFunctionReturn(0);
5855c6c1daeSBarry Smith }
5865c6c1daeSBarry Smith 
5875c6c1daeSBarry Smith /*@
5885c6c1daeSBarry Smith     PetscViewerBinarySetSkipHeader - do not write a header with size information on output, just raw data
5895c6c1daeSBarry Smith 
5905c6c1daeSBarry Smith     Not Collective
5915c6c1daeSBarry Smith 
5925c6c1daeSBarry Smith     Input Parameters:
593811af0c4SBarry Smith +   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
594811af0c4SBarry Smith -   skip - `PETSC_TRUE `means do not write header
5955c6c1daeSBarry Smith 
5965c6c1daeSBarry Smith     Options Database Key:
597811af0c4SBarry Smith .   -viewer_binary_skip_header <true or false> - true means do not write header
5985c6c1daeSBarry Smith 
5995c6c1daeSBarry Smith     Level: advanced
6005c6c1daeSBarry Smith 
60195452b02SPatrick Sanan     Notes:
602811af0c4SBarry Smith       This must be called after `PetscViewerSetType()`
6035c6c1daeSBarry Smith 
60410699b91SBarry Smith       Is ignored on anything but a binary viewer
6055c6c1daeSBarry Smith 
606811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
607db781477SPatrick Sanan           `PetscViewerBinaryGetSkipHeader()`
6085c6c1daeSBarry Smith @*/
6099371c9d4SSatish Balay PetscErrorCode PetscViewerBinarySetSkipHeader(PetscViewer viewer, PetscBool skip) {
6105c6c1daeSBarry Smith   PetscFunctionBegin;
611cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
612cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, skip, 2);
613cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetSkipHeader_C", (PetscViewer, PetscBool), (viewer, skip));
6145c6c1daeSBarry Smith   PetscFunctionReturn(0);
6155c6c1daeSBarry Smith }
6165c6c1daeSBarry Smith 
6179371c9d4SSatish Balay static PetscErrorCode PetscViewerBinarySetSkipHeader_Binary(PetscViewer viewer, PetscBool skip) {
6185c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
6195c6c1daeSBarry Smith 
6205c6c1daeSBarry Smith   PetscFunctionBegin;
621cc843e7aSLisandro Dalcin   vbinary->skipheader = skip;
6225c6c1daeSBarry Smith   PetscFunctionReturn(0);
6235c6c1daeSBarry Smith }
6245c6c1daeSBarry Smith 
6255c6c1daeSBarry Smith /*@
6265c6c1daeSBarry Smith     PetscViewerBinaryGetSkipHeader - checks whether to write a header with size information on output, or just raw data
6275c6c1daeSBarry Smith 
6285c6c1daeSBarry Smith     Not Collective
6295c6c1daeSBarry Smith 
6305c6c1daeSBarry Smith     Input Parameter:
631811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
6325c6c1daeSBarry Smith 
6335c6c1daeSBarry Smith     Output Parameter:
634811af0c4SBarry Smith .   skip - `PETSC_TRUE` means do not write header
6355c6c1daeSBarry Smith 
6365c6c1daeSBarry Smith     Level: advanced
6375c6c1daeSBarry Smith 
63895452b02SPatrick Sanan     Notes:
63995452b02SPatrick Sanan     This must be called after PetscViewerSetType()
6405c6c1daeSBarry Smith 
641811af0c4SBarry Smith     Returns `PETSC_FALSE` for `PETSCSOCKETVIEWER`, you cannot skip the header for it.
6425c6c1daeSBarry Smith 
643811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
644db781477SPatrick Sanan           `PetscViewerBinarySetSkipHeader()`
6455c6c1daeSBarry Smith @*/
6469371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryGetSkipHeader(PetscViewer viewer, PetscBool *skip) {
6475c6c1daeSBarry Smith   PetscFunctionBegin;
648cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
649cc843e7aSLisandro Dalcin   PetscValidBoolPointer(skip, 2);
650cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetSkipHeader_C", (PetscViewer, PetscBool *), (viewer, skip));
6515c6c1daeSBarry Smith   PetscFunctionReturn(0);
6525c6c1daeSBarry Smith }
6535c6c1daeSBarry Smith 
6549371c9d4SSatish Balay static PetscErrorCode PetscViewerBinaryGetSkipHeader_Binary(PetscViewer viewer, PetscBool *skip) {
655f3b211e4SSatish Balay   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
6565c6c1daeSBarry Smith 
6575c6c1daeSBarry Smith   PetscFunctionBegin;
658cc843e7aSLisandro Dalcin   *skip = vbinary->skipheader;
6595c6c1daeSBarry Smith   PetscFunctionReturn(0);
6605c6c1daeSBarry Smith }
6615c6c1daeSBarry Smith 
6625c6c1daeSBarry Smith /*@C
6635c6c1daeSBarry Smith     PetscViewerBinaryGetInfoPointer - Extracts the file pointer for the ASCII
6645c6c1daeSBarry Smith           info file associated with a binary file.
6655c6c1daeSBarry Smith 
6665c6c1daeSBarry Smith     Not Collective
6675c6c1daeSBarry Smith 
6685c6c1daeSBarry Smith     Input Parameter:
669811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
6705c6c1daeSBarry Smith 
6715c6c1daeSBarry Smith     Output Parameter:
6720298fd71SBarry Smith .   file - file pointer  Always returns NULL if not a binary viewer
6735c6c1daeSBarry Smith 
6745c6c1daeSBarry Smith     Level: advanced
6755c6c1daeSBarry Smith 
676811af0c4SBarry Smith     Note:
677811af0c4SBarry Smith       For writable binary `PetscViewer`s, the file pointer will only be valid for the
678811af0c4SBarry Smith     first processor in the communicator that shares the `PetscViewer`.
6795c6c1daeSBarry Smith 
6805c6c1daeSBarry Smith     Fortran Note:
6815c6c1daeSBarry Smith     This routine is not supported in Fortran.
6825c6c1daeSBarry Smith 
683811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`
6845c6c1daeSBarry Smith @*/
6859371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryGetInfoPointer(PetscViewer viewer, FILE **file) {
6865c6c1daeSBarry Smith   PetscFunctionBegin;
687cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
688cc843e7aSLisandro Dalcin   PetscValidPointer(file, 2);
6890298fd71SBarry Smith   *file = NULL;
690cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinaryGetInfoPointer_C", (PetscViewer, FILE **), (viewer, file));
6915c6c1daeSBarry Smith   PetscFunctionReturn(0);
6925c6c1daeSBarry Smith }
6935c6c1daeSBarry Smith 
6949371c9d4SSatish Balay static PetscErrorCode PetscViewerBinaryGetInfoPointer_Binary(PetscViewer viewer, FILE **file) {
695cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
6965c6c1daeSBarry Smith 
6975c6c1daeSBarry Smith   PetscFunctionBegin;
6989566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
699cc843e7aSLisandro Dalcin   *file = vbinary->fdes_info;
700cc843e7aSLisandro Dalcin   if (viewer->format == PETSC_VIEWER_BINARY_MATLAB && !vbinary->matlabheaderwritten) {
7015c6c1daeSBarry Smith     if (vbinary->fdes_info) {
702cc843e7aSLisandro Dalcin       FILE *info = vbinary->fdes_info;
7039566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- begin code written by PetscViewerBinary for MATLAB format ---#\n"));
7049566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#$$ Set.filename = '%s';\n", vbinary->filename));
7059566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#$$ fd = PetscOpenFile(Set.filename);\n"));
7069566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n"));
707cc843e7aSLisandro Dalcin     }
708cc843e7aSLisandro Dalcin     vbinary->matlabheaderwritten = PETSC_TRUE;
7095c6c1daeSBarry Smith   }
7105c6c1daeSBarry Smith   PetscFunctionReturn(0);
7115c6c1daeSBarry Smith }
7125c6c1daeSBarry Smith 
713e0385b85SDave May #if defined(PETSC_HAVE_MPIIO)
7149371c9d4SSatish Balay static PetscErrorCode PetscViewerFileClose_BinaryMPIIO(PetscViewer v) {
715e0385b85SDave May   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
716e0385b85SDave May 
717e0385b85SDave May   PetscFunctionBegin;
71848a46eb9SPierre Jolivet   if (vbinary->mfdes != MPI_FILE_NULL) PetscCallMPI(MPI_File_close(&vbinary->mfdes));
71948a46eb9SPierre Jolivet   if (vbinary->mfsub != MPI_FILE_NULL) PetscCallMPI(MPI_File_close(&vbinary->mfsub));
720cc843e7aSLisandro Dalcin   vbinary->moff = 0;
721e0385b85SDave May   PetscFunctionReturn(0);
722e0385b85SDave May }
723e0385b85SDave May #endif
724e0385b85SDave May 
7259371c9d4SSatish Balay static PetscErrorCode PetscViewerFileClose_BinarySTDIO(PetscViewer v) {
726cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
727cc843e7aSLisandro Dalcin 
728cc843e7aSLisandro Dalcin   PetscFunctionBegin;
729cc843e7aSLisandro Dalcin   if (vbinary->fdes != -1) {
7309566063dSJacob Faibussowitsch     PetscCall(PetscBinaryClose(vbinary->fdes));
731cc843e7aSLisandro Dalcin     vbinary->fdes = -1;
732cc843e7aSLisandro Dalcin     if (vbinary->storecompressed) {
733cc843e7aSLisandro Dalcin       char        cmd[8 + PETSC_MAX_PATH_LEN], out[64 + PETSC_MAX_PATH_LEN] = "";
734cc843e7aSLisandro Dalcin       const char *gzfilename = vbinary->ogzfilename ? vbinary->ogzfilename : vbinary->filename;
735cc843e7aSLisandro Dalcin       /* compress the file */
7369566063dSJacob Faibussowitsch       PetscCall(PetscStrncpy(cmd, "gzip -f ", sizeof(cmd)));
7379566063dSJacob Faibussowitsch       PetscCall(PetscStrlcat(cmd, gzfilename, sizeof(cmd)));
738cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_POPEN)
739cc843e7aSLisandro Dalcin       {
740cc843e7aSLisandro Dalcin         FILE *fp;
7419566063dSJacob Faibussowitsch         PetscCall(PetscPOpen(PETSC_COMM_SELF, NULL, cmd, "r", &fp));
742cc73adaaSBarry Smith         PetscCheck(!fgets(out, (int)(sizeof(out) - 1), fp), PETSC_COMM_SELF, PETSC_ERR_LIB, "Error from command %s\n%s", cmd, out);
7439566063dSJacob Faibussowitsch         PetscCall(PetscPClose(PETSC_COMM_SELF, fp));
744cc843e7aSLisandro Dalcin       }
745cc843e7aSLisandro Dalcin #endif
746cc843e7aSLisandro Dalcin     }
747cc843e7aSLisandro Dalcin   }
7489566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary->ogzfilename));
749cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
750cc843e7aSLisandro Dalcin }
751cc843e7aSLisandro Dalcin 
7529371c9d4SSatish Balay static PetscErrorCode PetscViewerFileClose_BinaryInfo(PetscViewer v) {
753cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
754cc843e7aSLisandro Dalcin 
755cc843e7aSLisandro Dalcin   PetscFunctionBegin;
756cc843e7aSLisandro Dalcin   if (v->format == PETSC_VIEWER_BINARY_MATLAB && vbinary->matlabheaderwritten) {
757cc843e7aSLisandro Dalcin     if (vbinary->fdes_info) {
758cc843e7aSLisandro Dalcin       FILE *info = vbinary->fdes_info;
7599566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- begin code written by PetscViewerBinary for MATLAB format ---#\n"));
7609566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#$$ close(fd);\n"));
7619566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n"));
762cc843e7aSLisandro Dalcin     }
763cc843e7aSLisandro Dalcin   }
764cc843e7aSLisandro Dalcin   if (vbinary->fdes_info) {
765cc843e7aSLisandro Dalcin     FILE *info         = vbinary->fdes_info;
766cc843e7aSLisandro Dalcin     vbinary->fdes_info = NULL;
767cc73adaaSBarry Smith     PetscCheck(!fclose(info), PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file");
768cc843e7aSLisandro Dalcin   }
769cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
770cc843e7aSLisandro Dalcin }
771cc843e7aSLisandro Dalcin 
7729371c9d4SSatish Balay static PetscErrorCode PetscViewerFileClose_Binary(PetscViewer v) {
773cc843e7aSLisandro Dalcin   PetscFunctionBegin;
774cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
7759566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_BinaryMPIIO(v));
776cc843e7aSLisandro Dalcin #endif
7779566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_BinarySTDIO(v));
7789566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_BinaryInfo(v));
779cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
780cc843e7aSLisandro Dalcin }
781cc843e7aSLisandro Dalcin 
7829371c9d4SSatish Balay static PetscErrorCode PetscViewerDestroy_Binary(PetscViewer v) {
7835c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
7845c6c1daeSBarry Smith 
7855c6c1daeSBarry Smith   PetscFunctionBegin;
7869566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_Binary(v));
7879566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary->filename));
7889566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary));
7892e956fe4SStefano Zampini   PetscCall(PetscViewerBinaryClearFunctionList(v));
790e0385b85SDave May   PetscFunctionReturn(0);
791e0385b85SDave May }
7925c6c1daeSBarry Smith 
7935c6c1daeSBarry Smith /*@C
7945c6c1daeSBarry Smith    PetscViewerBinaryOpen - Opens a file for binary input/output.
7955c6c1daeSBarry Smith 
796d083f849SBarry Smith    Collective
7975c6c1daeSBarry Smith 
7985c6c1daeSBarry Smith    Input Parameters:
7995c6c1daeSBarry Smith +  comm - MPI communicator
8005c6c1daeSBarry Smith .  name - name of file
801cc843e7aSLisandro Dalcin -  mode - open mode of file
8025c6c1daeSBarry Smith $    FILE_MODE_WRITE - create new file for binary output
8035c6c1daeSBarry Smith $    FILE_MODE_READ - open existing file for binary input
8045c6c1daeSBarry Smith $    FILE_MODE_APPEND - open existing file for binary output
8055c6c1daeSBarry Smith 
8065c6c1daeSBarry Smith    Output Parameter:
807cc843e7aSLisandro Dalcin .  viewer - PetscViewer for binary input/output to use with the specified file
8085c6c1daeSBarry Smith 
8095c6c1daeSBarry Smith     Options Database Keys:
810811af0c4SBarry Smith +    -viewer_binary_filename <name> - name of file to use
811811af0c4SBarry Smith .    -viewer_binary_skip_info - true to skip opening an info file
812811af0c4SBarry Smith .    -viewer_binary_skip_options - true to not use options database while creating viewer
813811af0c4SBarry Smith .    -viewer_binary_skip_header - true to skip output object headers to the file
814811af0c4SBarry Smith -    -viewer_binary_mpiio - true to use MPI-IO for input and output to the file (more scalable for large problems)
8155c6c1daeSBarry Smith 
8165c6c1daeSBarry Smith    Level: beginner
8175c6c1daeSBarry Smith 
8185c6c1daeSBarry Smith    Note:
819811af0c4SBarry Smith    This `PetscViewer` should be destroyed with `PetscViewerDestroy()`.
8205c6c1daeSBarry Smith 
8215c6c1daeSBarry Smith     For reading files, the filename may begin with ftp:// or http:// and/or
8225c6c1daeSBarry Smith     end with .gz; in this case file is brought over and uncompressed.
8235c6c1daeSBarry Smith 
8245c6c1daeSBarry Smith     For creating files, if the file name ends with .gz it is automatically
8255c6c1daeSBarry Smith     compressed when closed.
8265c6c1daeSBarry Smith 
827811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
828db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
829db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`, `PetscViewerBinarySetUseMPIIO()`,
830db781477SPatrick Sanan           `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
8315c6c1daeSBarry Smith @*/
8329371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryOpen(MPI_Comm comm, const char name[], PetscFileMode mode, PetscViewer *viewer) {
8335c6c1daeSBarry Smith   PetscFunctionBegin;
8349566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(comm, viewer));
8359566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(*viewer, PETSCVIEWERBINARY));
8369566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(*viewer, mode));
8379566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetName(*viewer, name));
8389566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions(*viewer));
8395c6c1daeSBarry Smith   PetscFunctionReturn(0);
8405c6c1daeSBarry Smith }
8415c6c1daeSBarry Smith 
8425c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
8439371c9d4SSatish Balay static PetscErrorCode PetscViewerBinaryWriteReadMPIIO(PetscViewer viewer, void *data, PetscInt num, PetscInt *count, PetscDataType dtype, PetscBool write) {
84422a8f86cSLisandro Dalcin   MPI_Comm            comm    = PetscObjectComm((PetscObject)viewer);
8455c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
84622a8f86cSLisandro Dalcin   MPI_File            mfdes   = vbinary->mfdes;
8475c6c1daeSBarry Smith   MPI_Datatype        mdtype;
84869e10bbaSLisandro Dalcin   PetscMPIInt         rank, cnt;
8495c6c1daeSBarry Smith   MPI_Status          status;
8505c6c1daeSBarry Smith   MPI_Aint            ul, dsize;
8515c6c1daeSBarry Smith 
8525c6c1daeSBarry Smith   PetscFunctionBegin;
8539566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
8549566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(num, &cnt));
8559566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToMPIDataType(dtype, &mdtype));
8565c6c1daeSBarry Smith   if (write) {
85748a46eb9SPierre Jolivet     if (rank == 0) PetscCall(MPIU_File_write_at(mfdes, vbinary->moff, data, cnt, mdtype, &status));
8585c6c1daeSBarry Smith   } else {
859dd400576SPatrick Sanan     if (rank == 0) {
8609566063dSJacob Faibussowitsch       PetscCall(MPIU_File_read_at(mfdes, vbinary->moff, data, cnt, mdtype, &status));
8619566063dSJacob Faibussowitsch       if (cnt > 0) PetscCallMPI(MPI_Get_count(&status, mdtype, &cnt));
8625c6c1daeSBarry Smith     }
8639566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&cnt, 1, MPI_INT, 0, comm));
8649566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(data, cnt, mdtype, 0, comm));
86569e10bbaSLisandro Dalcin   }
8669566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_get_extent(mdtype, &ul, &dsize));
8675c6c1daeSBarry Smith   vbinary->moff += dsize * cnt;
8689860990eSLisandro Dalcin   if (count) *count = cnt;
8695c6c1daeSBarry Smith   PetscFunctionReturn(0);
8705c6c1daeSBarry Smith }
8715c6c1daeSBarry Smith #endif
8725c6c1daeSBarry Smith 
8735c6c1daeSBarry Smith /*@C
8745c6c1daeSBarry Smith    PetscViewerBinaryRead - Reads from a binary file, all processors get the same result
8755c6c1daeSBarry Smith 
876d083f849SBarry Smith    Collective
8775c6c1daeSBarry Smith 
8785c6c1daeSBarry Smith    Input Parameters:
8795c6c1daeSBarry Smith +  viewer - the binary viewer
880907376e6SBarry Smith .  data - location of the data to be written
881060da220SMatthew G. Knepley .  num - number of items of data to read
882907376e6SBarry Smith -  dtype - type of data to read
8835c6c1daeSBarry Smith 
884f8e4bde8SMatthew G. Knepley    Output Parameters:
8855972f5f3SLisandro Dalcin .  count - number of items of data actually read, or NULL.
886f8e4bde8SMatthew G. Knepley 
8875c6c1daeSBarry Smith    Level: beginner
8885c6c1daeSBarry Smith 
889811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
890db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
891db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
8925c6c1daeSBarry Smith @*/
8939371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryRead(PetscViewer viewer, void *data, PetscInt num, PetscInt *count, PetscDataType dtype) {
89422a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
8955c6c1daeSBarry Smith 
89622a8f86cSLisandro Dalcin   PetscFunctionBegin;
89722a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
89822a8f86cSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, num, 3);
8999566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
90022a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
9015c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
902bc196f7cSDave May   if (vbinary->usempiio) {
9039566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWriteReadMPIIO(viewer, data, num, count, dtype, PETSC_FALSE));
9045c6c1daeSBarry Smith   } else {
9055c6c1daeSBarry Smith #endif
9069566063dSJacob Faibussowitsch     PetscCall(PetscBinarySynchronizedRead(PetscObjectComm((PetscObject)viewer), vbinary->fdes, data, num, count, dtype));
9075c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
9085c6c1daeSBarry Smith   }
9095c6c1daeSBarry Smith #endif
9105c6c1daeSBarry Smith   PetscFunctionReturn(0);
9115c6c1daeSBarry Smith }
9125c6c1daeSBarry Smith 
9135c6c1daeSBarry Smith /*@C
914811af0c4SBarry Smith    PetscViewerBinaryWrite - writes to a binary file, only from the first MPI rank
9155c6c1daeSBarry Smith 
916d083f849SBarry Smith    Collective
9175c6c1daeSBarry Smith 
9185c6c1daeSBarry Smith    Input Parameters:
9195c6c1daeSBarry Smith +  viewer - the binary viewer
9205c6c1daeSBarry Smith .  data - location of data
9215824c9d0SPatrick Sanan .  count - number of items of data to write
922f253e43cSLisandro Dalcin -  dtype - type of data to write
9235c6c1daeSBarry Smith 
9245c6c1daeSBarry Smith    Level: beginner
9255c6c1daeSBarry Smith 
926811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
927db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`, `PetscDataType`
928db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
9295c6c1daeSBarry Smith @*/
9309371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryWrite(PetscViewer viewer, const void *data, PetscInt count, PetscDataType dtype) {
93122a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
9325c6c1daeSBarry Smith 
9335c6c1daeSBarry Smith   PetscFunctionBegin;
93422a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
93522a8f86cSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, count, 3);
9369566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
93722a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
9385c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
939bc196f7cSDave May   if (vbinary->usempiio) {
9409566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWriteReadMPIIO(viewer, (void *)data, count, NULL, dtype, PETSC_TRUE));
9415c6c1daeSBarry Smith   } else {
9425c6c1daeSBarry Smith #endif
9439566063dSJacob Faibussowitsch     PetscCall(PetscBinarySynchronizedWrite(PetscObjectComm((PetscObject)viewer), vbinary->fdes, data, count, dtype));
9445c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
9455c6c1daeSBarry Smith   }
9465c6c1daeSBarry Smith #endif
9475c6c1daeSBarry Smith   PetscFunctionReturn(0);
9485c6c1daeSBarry Smith }
9495c6c1daeSBarry Smith 
9509371c9d4SSatish Balay static PetscErrorCode PetscViewerBinaryWriteReadAll(PetscViewer viewer, PetscBool write, void *data, PetscInt count, PetscInt start, PetscInt total, PetscDataType dtype) {
9515972f5f3SLisandro Dalcin   MPI_Comm              comm = PetscObjectComm((PetscObject)viewer);
9525972f5f3SLisandro Dalcin   PetscMPIInt           size, rank;
9535972f5f3SLisandro Dalcin   MPI_Datatype          mdtype;
954ec4bef21SJose E. Roman   PETSC_UNUSED MPI_Aint lb;
955ec4bef21SJose E. Roman   MPI_Aint              dsize;
9565972f5f3SLisandro Dalcin   PetscBool             useMPIIO;
9575972f5f3SLisandro Dalcin 
9585972f5f3SLisandro Dalcin   PetscFunctionBegin;
9595972f5f3SLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
960064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveBool(viewer, ((start >= 0) || (start == PETSC_DETERMINE)), 5);
961064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveBool(viewer, ((total >= 0) || (total == PETSC_DETERMINE)), 6);
962064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveInt(viewer, total, 6);
9639566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
9645972f5f3SLisandro Dalcin 
9659566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToMPIDataType(dtype, &mdtype));
9669566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_get_extent(mdtype, &lb, &dsize));
9679566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
9689566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
9695972f5f3SLisandro Dalcin 
9709566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetUseMPIIO(viewer, &useMPIIO));
9715972f5f3SLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
9725972f5f3SLisandro Dalcin   if (useMPIIO) {
9735972f5f3SLisandro Dalcin     MPI_File    mfdes;
9745972f5f3SLisandro Dalcin     MPI_Offset  off;
9755972f5f3SLisandro Dalcin     PetscMPIInt cnt;
9765972f5f3SLisandro Dalcin 
9775972f5f3SLisandro Dalcin     if (start == PETSC_DETERMINE) {
9789566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Scan(&count, &start, 1, MPIU_INT, MPI_SUM, comm));
9795972f5f3SLisandro Dalcin       start -= count;
9805972f5f3SLisandro Dalcin     }
9815972f5f3SLisandro Dalcin     if (total == PETSC_DETERMINE) {
9825972f5f3SLisandro Dalcin       total = start + count;
9839566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Bcast(&total, 1, MPIU_INT, size - 1, comm));
9845972f5f3SLisandro Dalcin     }
9859566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(count, &cnt));
9869566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryGetMPIIODescriptor(viewer, &mfdes));
9879566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryGetMPIIOOffset(viewer, &off));
9885972f5f3SLisandro Dalcin     off += (MPI_Offset)(start * dsize);
9895972f5f3SLisandro Dalcin     if (write) {
9909566063dSJacob Faibussowitsch       PetscCall(MPIU_File_write_at_all(mfdes, off, data, cnt, mdtype, MPI_STATUS_IGNORE));
9915972f5f3SLisandro Dalcin     } else {
9929566063dSJacob Faibussowitsch       PetscCall(MPIU_File_read_at_all(mfdes, off, data, cnt, mdtype, MPI_STATUS_IGNORE));
9935972f5f3SLisandro Dalcin     }
9945972f5f3SLisandro Dalcin     off = (MPI_Offset)(total * dsize);
9959566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryAddMPIIOOffset(viewer, off));
9965972f5f3SLisandro Dalcin     PetscFunctionReturn(0);
9975972f5f3SLisandro Dalcin   }
9985972f5f3SLisandro Dalcin #endif
9995972f5f3SLisandro Dalcin   {
10005972f5f3SLisandro Dalcin     int         fdes;
10015972f5f3SLisandro Dalcin     char       *workbuf = NULL;
1002dd400576SPatrick Sanan     PetscInt    tcount = rank == 0 ? 0 : count, maxcount = 0, message_count, flowcontrolcount;
10035972f5f3SLisandro Dalcin     PetscMPIInt tag, cnt, maxcnt, scnt = 0, rcnt = 0, j;
10045972f5f3SLisandro Dalcin     MPI_Status  status;
10055972f5f3SLisandro Dalcin 
10069566063dSJacob Faibussowitsch     PetscCall(PetscCommGetNewTag(comm, &tag));
10079566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Reduce(&tcount, &maxcount, 1, MPIU_INT, MPI_MAX, 0, comm));
10089566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(maxcount, &maxcnt));
10099566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(count, &cnt));
10105972f5f3SLisandro Dalcin 
10119566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryGetDescriptor(viewer, &fdes));
10129566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlowControlStart(viewer, &message_count, &flowcontrolcount));
1013dd400576SPatrick Sanan     if (rank == 0) {
10149566063dSJacob Faibussowitsch       PetscCall(PetscMalloc(maxcnt * dsize, &workbuf));
10155972f5f3SLisandro Dalcin       if (write) {
10169566063dSJacob Faibussowitsch         PetscCall(PetscBinaryWrite(fdes, data, cnt, dtype));
10175972f5f3SLisandro Dalcin       } else {
10189566063dSJacob Faibussowitsch         PetscCall(PetscBinaryRead(fdes, data, cnt, NULL, dtype));
10195972f5f3SLisandro Dalcin       }
10205972f5f3SLisandro Dalcin       for (j = 1; j < size; j++) {
10219566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlowControlStepMain(viewer, j, &message_count, flowcontrolcount));
10225972f5f3SLisandro Dalcin         if (write) {
10239566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Recv(workbuf, maxcnt, mdtype, j, tag, comm, &status));
10249566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Get_count(&status, mdtype, &rcnt));
10259566063dSJacob Faibussowitsch           PetscCall(PetscBinaryWrite(fdes, workbuf, rcnt, dtype));
10265972f5f3SLisandro Dalcin         } else {
10279566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Recv(&scnt, 1, MPI_INT, j, tag, comm, MPI_STATUS_IGNORE));
10289566063dSJacob Faibussowitsch           PetscCall(PetscBinaryRead(fdes, workbuf, scnt, NULL, dtype));
10299566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Send(workbuf, scnt, mdtype, j, tag, comm));
10305972f5f3SLisandro Dalcin         }
10315972f5f3SLisandro Dalcin       }
10329566063dSJacob Faibussowitsch       PetscCall(PetscFree(workbuf));
10339566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlowControlEndMain(viewer, &message_count));
10345972f5f3SLisandro Dalcin     } else {
10359566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlowControlStepWorker(viewer, rank, &message_count));
10365972f5f3SLisandro Dalcin       if (write) {
10379566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(data, cnt, mdtype, 0, tag, comm));
10385972f5f3SLisandro Dalcin       } else {
10399566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(&cnt, 1, MPI_INT, 0, tag, comm));
10409566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Recv(data, cnt, mdtype, 0, tag, comm, MPI_STATUS_IGNORE));
10415972f5f3SLisandro Dalcin       }
10429566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlowControlEndWorker(viewer, &message_count));
10435972f5f3SLisandro Dalcin     }
10445972f5f3SLisandro Dalcin   }
10455972f5f3SLisandro Dalcin   PetscFunctionReturn(0);
10465972f5f3SLisandro Dalcin }
10475972f5f3SLisandro Dalcin 
10485972f5f3SLisandro Dalcin /*@C
1049811af0c4SBarry Smith    PetscViewerBinaryReadAll - reads from a binary file from all MPI ranks, each rank receives its own portion of the data
10505972f5f3SLisandro Dalcin 
10515972f5f3SLisandro Dalcin    Collective
10525972f5f3SLisandro Dalcin 
10535972f5f3SLisandro Dalcin    Input Parameters:
10545972f5f3SLisandro Dalcin +  viewer - the binary viewer
10555972f5f3SLisandro Dalcin .  data - location of data
10565972f5f3SLisandro Dalcin .  count - local number of items of data to read
1057811af0c4SBarry Smith .  start - local start, can be `PETSC_DETERMINE`
1058811af0c4SBarry Smith .  total - global number of items of data to read, can be `PETSC_DETERMINE`
10595972f5f3SLisandro Dalcin -  dtype - type of data to read
10605972f5f3SLisandro Dalcin 
10615972f5f3SLisandro Dalcin    Level: advanced
10625972f5f3SLisandro Dalcin 
1063811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryRead()`, `PetscViewerBinaryWriteAll()`
10645972f5f3SLisandro Dalcin @*/
10659371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryReadAll(PetscViewer viewer, void *data, PetscInt count, PetscInt start, PetscInt total, PetscDataType dtype) {
10665972f5f3SLisandro Dalcin   PetscFunctionBegin;
10679566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWriteReadAll(viewer, PETSC_FALSE, data, count, start, total, dtype));
10685972f5f3SLisandro Dalcin   PetscFunctionReturn(0);
10695972f5f3SLisandro Dalcin }
10705972f5f3SLisandro Dalcin 
10715972f5f3SLisandro Dalcin /*@C
1072811af0c4SBarry Smith    PetscViewerBinaryWriteAll - writes to a binary file from all MPI ranks, each rank writes its own portion of the data
10735972f5f3SLisandro Dalcin 
10745972f5f3SLisandro Dalcin    Collective
10755972f5f3SLisandro Dalcin 
10765972f5f3SLisandro Dalcin    Input Parameters:
10775972f5f3SLisandro Dalcin +  viewer - the binary viewer
10785972f5f3SLisandro Dalcin .  data - location of data
10795972f5f3SLisandro Dalcin .  count - local number of items of data to write
1080811af0c4SBarry Smith .  start - local start, can be `PETSC_DETERMINE`
1081811af0c4SBarry Smith .  total - global number of items of data to write, can be `PETSC_DETERMINE`
10825972f5f3SLisandro Dalcin -  dtype - type of data to write
10835972f5f3SLisandro Dalcin 
10845972f5f3SLisandro Dalcin    Level: advanced
10855972f5f3SLisandro Dalcin 
1086811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryWriteAll()`, `PetscViewerBinaryReadAll()`
10875972f5f3SLisandro Dalcin @*/
10889371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryWriteAll(PetscViewer viewer, const void *data, PetscInt count, PetscInt start, PetscInt total, PetscDataType dtype) {
10895972f5f3SLisandro Dalcin   PetscFunctionBegin;
10909566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWriteReadAll(viewer, PETSC_TRUE, (void *)data, count, start, total, dtype));
10915972f5f3SLisandro Dalcin   PetscFunctionReturn(0);
10925972f5f3SLisandro Dalcin }
10935972f5f3SLisandro Dalcin 
10945c6c1daeSBarry Smith /*@C
1095811af0c4SBarry Smith    PetscViewerBinaryWriteStringArray - writes to a binary file, only from the first MPI rank, an array of strings
10965c6c1daeSBarry Smith 
1097d083f849SBarry Smith    Collective
10985c6c1daeSBarry Smith 
10995c6c1daeSBarry Smith    Input Parameters:
11005c6c1daeSBarry Smith +  viewer - the binary viewer
11015c6c1daeSBarry Smith -  data - location of the array of strings
11025c6c1daeSBarry Smith 
11035c6c1daeSBarry Smith    Level: intermediate
11045c6c1daeSBarry Smith 
1105811af0c4SBarry Smith     Note:
1106811af0c4SBarry Smith     The array of strings must be null terminated
11075c6c1daeSBarry Smith 
1108811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
1109db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
1110db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
11115c6c1daeSBarry Smith @*/
11129371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryWriteStringArray(PetscViewer viewer, const char *const *data) {
11135c6c1daeSBarry Smith   PetscInt i, n = 0, *sizes;
1114cc843e7aSLisandro Dalcin   size_t   len;
11155c6c1daeSBarry Smith 
111622a8f86cSLisandro Dalcin   PetscFunctionBegin;
11179566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
11185c6c1daeSBarry Smith   /* count number of strings */
11199371c9d4SSatish Balay   while (data[n++])
11209371c9d4SSatish Balay     ;
11215c6c1daeSBarry Smith   n--;
11229566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n + 1, &sizes));
11235c6c1daeSBarry Smith   sizes[0] = n;
11245c6c1daeSBarry Smith   for (i = 0; i < n; i++) {
11259566063dSJacob Faibussowitsch     PetscCall(PetscStrlen(data[i], &len));
1126cc843e7aSLisandro Dalcin     sizes[i + 1] = (PetscInt)len + 1; /* size includes space for the null terminator */
11275c6c1daeSBarry Smith   }
11289566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWrite(viewer, sizes, n + 1, PETSC_INT));
112948a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall(PetscViewerBinaryWrite(viewer, (void *)data[i], sizes[i + 1], PETSC_CHAR));
11309566063dSJacob Faibussowitsch   PetscCall(PetscFree(sizes));
11315c6c1daeSBarry Smith   PetscFunctionReturn(0);
11325c6c1daeSBarry Smith }
11335c6c1daeSBarry Smith 
11345c6c1daeSBarry Smith /*@C
1135811af0c4SBarry Smith    PetscViewerBinaryReadStringArray - reads a binary file an array of strings to all MPI ranks
11365c6c1daeSBarry Smith 
1137d083f849SBarry Smith    Collective
11385c6c1daeSBarry Smith 
11395c6c1daeSBarry Smith    Input Parameter:
11405c6c1daeSBarry Smith .  viewer - the binary viewer
11415c6c1daeSBarry Smith 
11425c6c1daeSBarry Smith    Output Parameter:
11435c6c1daeSBarry Smith .  data - location of the array of strings
11445c6c1daeSBarry Smith 
11455c6c1daeSBarry Smith    Level: intermediate
11465c6c1daeSBarry Smith 
1147811af0c4SBarry Smith     Note:
1148811af0c4SBarry Smith     The array of strings must null terminated
11495c6c1daeSBarry Smith 
1150811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
1151db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
1152db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
11535c6c1daeSBarry Smith @*/
11549371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryReadStringArray(PetscViewer viewer, char ***data) {
1155060da220SMatthew G. Knepley   PetscInt i, n, *sizes, N = 0;
11565c6c1daeSBarry Smith 
115722a8f86cSLisandro Dalcin   PetscFunctionBegin;
11589566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
11595c6c1daeSBarry Smith   /* count number of strings */
11609566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &n, 1, NULL, PETSC_INT));
11619566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &sizes));
11629566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, sizes, n, NULL, PETSC_INT));
1163a297a907SKarl Rupp   for (i = 0; i < n; i++) N += sizes[i];
11649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc((n + 1) * sizeof(char *) + N * sizeof(char), data));
11655c6c1daeSBarry Smith   (*data)[0] = (char *)((*data) + n + 1);
1166a297a907SKarl Rupp   for (i = 1; i < n; i++) (*data)[i] = (*data)[i - 1] + sizes[i - 1];
11679566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, (*data)[0], N, NULL, PETSC_CHAR));
116802c9f0b5SLisandro Dalcin   (*data)[n] = NULL;
11699566063dSJacob Faibussowitsch   PetscCall(PetscFree(sizes));
11705c6c1daeSBarry Smith   PetscFunctionReturn(0);
11715c6c1daeSBarry Smith }
11725c6c1daeSBarry Smith 
1173cc843e7aSLisandro Dalcin /*@C
1174cc843e7aSLisandro Dalcin      PetscViewerFileSetMode - Sets the open mode of file
1175cc843e7aSLisandro Dalcin 
1176811af0c4SBarry Smith     Logically Collective on viewer
1177cc843e7aSLisandro Dalcin 
1178cc843e7aSLisandro Dalcin   Input Parameters:
1179811af0c4SBarry Smith +  viewer - the `PetscViewer`; must be a a `PETSCVIEWERBINARY`, `PETSCVIEWERMATLAB`, `PETSCVIEWERHDF5`, or `PETSCVIEWERASCII`  `PetscViewer`
1180cc843e7aSLisandro Dalcin -  mode - open mode of file
1181cc843e7aSLisandro Dalcin $    FILE_MODE_WRITE - create new file for output
1182cc843e7aSLisandro Dalcin $    FILE_MODE_READ - open existing file for input
1183cc843e7aSLisandro Dalcin $    FILE_MODE_APPEND - open existing file for output
1184cc843e7aSLisandro Dalcin 
1185cc843e7aSLisandro Dalcin   Level: advanced
1186cc843e7aSLisandro Dalcin 
1187db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
1188cc843e7aSLisandro Dalcin @*/
11899371c9d4SSatish Balay PetscErrorCode PetscViewerFileSetMode(PetscViewer viewer, PetscFileMode mode) {
1190cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1191cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1192cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveEnum(viewer, mode, 2);
119308401ef6SPierre Jolivet   PetscCheck(mode != FILE_MODE_UNDEFINED, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Cannot set FILE_MODE_UNDEFINED");
1194f7d195e4SLawrence Mitchell   PetscCheck(mode >= FILE_MODE_UNDEFINED && mode <= FILE_MODE_APPEND_UPDATE, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_OUTOFRANGE, "Invalid file mode %d", (int)mode);
1195cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerFileSetMode_C", (PetscViewer, PetscFileMode), (viewer, mode));
1196cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1197cc843e7aSLisandro Dalcin }
1198cc843e7aSLisandro Dalcin 
11999371c9d4SSatish Balay static PetscErrorCode PetscViewerFileSetMode_Binary(PetscViewer viewer, PetscFileMode mode) {
1200cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1201cc843e7aSLisandro Dalcin 
1202cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1203cc73adaaSBarry Smith   PetscCheck(!viewer->setupcalled || vbinary->filemode == mode, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Cannot change mode to %s after setup", PetscFileModes[mode]);
1204cc843e7aSLisandro Dalcin   vbinary->filemode = mode;
1205cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1206cc843e7aSLisandro Dalcin }
1207cc843e7aSLisandro Dalcin 
1208cc843e7aSLisandro Dalcin /*@C
1209cc843e7aSLisandro Dalcin      PetscViewerFileGetMode - Gets the open mode of file
1210cc843e7aSLisandro Dalcin 
1211cc843e7aSLisandro Dalcin     Not Collective
1212cc843e7aSLisandro Dalcin 
1213cc843e7aSLisandro Dalcin   Input Parameter:
1214811af0c4SBarry Smith .  viewer - the `PetscViewer`; must be a `PETSCVIEWERBINARY`, `PETSCVIEWERMATLAB`, `PETSCVIEWERHDF5`, or `PETSCVIEWERASCII`  `PetscViewer`
1215cc843e7aSLisandro Dalcin 
1216cc843e7aSLisandro Dalcin   Output Parameter:
1217cc843e7aSLisandro Dalcin .  mode - open mode of file
1218cc843e7aSLisandro Dalcin $    FILE_MODE_WRITE - create new file for binary output
1219cc843e7aSLisandro Dalcin $    FILE_MODE_READ - open existing file for binary input
1220cc843e7aSLisandro Dalcin $    FILE_MODE_APPEND - open existing file for binary output
1221cc843e7aSLisandro Dalcin 
1222cc843e7aSLisandro Dalcin   Level: advanced
1223cc843e7aSLisandro Dalcin 
1224db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
1225cc843e7aSLisandro Dalcin @*/
12269371c9d4SSatish Balay PetscErrorCode PetscViewerFileGetMode(PetscViewer viewer, PetscFileMode *mode) {
1227cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1228cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1229cc843e7aSLisandro Dalcin   PetscValidPointer(mode, 2);
1230cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerFileGetMode_C", (PetscViewer, PetscFileMode *), (viewer, mode));
1231cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1232cc843e7aSLisandro Dalcin }
1233cc843e7aSLisandro Dalcin 
12349371c9d4SSatish Balay static PetscErrorCode PetscViewerFileGetMode_Binary(PetscViewer viewer, PetscFileMode *mode) {
1235cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1236cc843e7aSLisandro Dalcin 
1237cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1238cc843e7aSLisandro Dalcin   *mode = vbinary->filemode;
1239cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1240cc843e7aSLisandro Dalcin }
1241cc843e7aSLisandro Dalcin 
12429371c9d4SSatish Balay static PetscErrorCode PetscViewerFileSetName_Binary(PetscViewer viewer, const char name[]) {
1243cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1244cc843e7aSLisandro Dalcin 
1245cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1246cc843e7aSLisandro Dalcin   if (viewer->setupcalled && vbinary->filename) {
1247cc843e7aSLisandro Dalcin     /* gzip can be run after the file with the previous filename has been closed */
12489566063dSJacob Faibussowitsch     PetscCall(PetscFree(vbinary->ogzfilename));
12499566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(vbinary->filename, &vbinary->ogzfilename));
1250cc843e7aSLisandro Dalcin   }
12519566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary->filename));
12529566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &vbinary->filename));
1253cc843e7aSLisandro Dalcin   viewer->setupcalled = PETSC_FALSE;
1254cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1255cc843e7aSLisandro Dalcin }
1256cc843e7aSLisandro Dalcin 
12579371c9d4SSatish Balay static PetscErrorCode PetscViewerFileGetName_Binary(PetscViewer viewer, const char **name) {
12585c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
12595c6c1daeSBarry Smith 
12605c6c1daeSBarry Smith   PetscFunctionBegin;
12615c6c1daeSBarry Smith   *name = vbinary->filename;
12625c6c1daeSBarry Smith   PetscFunctionReturn(0);
12635c6c1daeSBarry Smith }
12645c6c1daeSBarry Smith 
12655c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
12669371c9d4SSatish Balay static PetscErrorCode PetscViewerFileSetUp_BinaryMPIIO(PetscViewer viewer) {
12675c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1268e8a65b7dSLisandro Dalcin   int                 amode;
12695c6c1daeSBarry Smith 
12705c6c1daeSBarry Smith   PetscFunctionBegin;
1271a297a907SKarl Rupp   vbinary->storecompressed = PETSC_FALSE;
12725c6c1daeSBarry Smith 
1273cc843e7aSLisandro Dalcin   vbinary->moff = 0;
1274cc843e7aSLisandro Dalcin   switch (vbinary->filemode) {
1275e8a65b7dSLisandro Dalcin   case FILE_MODE_READ: amode = MPI_MODE_RDONLY; break;
1276e8a65b7dSLisandro Dalcin   case FILE_MODE_WRITE: amode = MPI_MODE_WRONLY | MPI_MODE_CREATE; break;
127722a8f86cSLisandro Dalcin   case FILE_MODE_APPEND: amode = MPI_MODE_WRONLY | MPI_MODE_CREATE | MPI_MODE_APPEND; break;
12787e4fd573SVaclav Hapla   case FILE_MODE_UNDEFINED: SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerSetUp()");
127998921bdaSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[vbinary->filemode]);
12805c6c1daeSBarry Smith   }
12819566063dSJacob Faibussowitsch   PetscCallMPI(MPI_File_open(PetscObjectComm((PetscObject)viewer), vbinary->filename, amode, MPI_INFO_NULL, &vbinary->mfdes));
128222a8f86cSLisandro Dalcin   /*
128322a8f86cSLisandro Dalcin       The MPI standard does not have MPI_MODE_TRUNCATE. We emulate this behavior by setting the file size to zero.
128422a8f86cSLisandro Dalcin   */
12859566063dSJacob Faibussowitsch   if (vbinary->filemode == FILE_MODE_WRITE) PetscCallMPI(MPI_File_set_size(vbinary->mfdes, 0));
128622a8f86cSLisandro Dalcin   /*
128722a8f86cSLisandro Dalcin       Initially, all processes view the file as a linear byte stream. Therefore, for files opened with MPI_MODE_APPEND,
128822a8f86cSLisandro Dalcin       MPI_File_get_position[_shared](fh, &offset) returns the absolute byte position at the end of file.
128922a8f86cSLisandro Dalcin       Otherwise, we would need to call MPI_File_get_byte_offset(fh, offset, &byte_offset) to convert
129022a8f86cSLisandro Dalcin       the offset in etype units to an absolute byte position.
129122a8f86cSLisandro Dalcin    */
12929566063dSJacob Faibussowitsch   if (vbinary->filemode == FILE_MODE_APPEND) PetscCallMPI(MPI_File_get_position(vbinary->mfdes, &vbinary->moff));
1293cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1294cc843e7aSLisandro Dalcin }
1295cc843e7aSLisandro Dalcin #endif
12965c6c1daeSBarry Smith 
12979371c9d4SSatish Balay static PetscErrorCode PetscViewerFileSetUp_BinarySTDIO(PetscViewer viewer) {
1298cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1299cc843e7aSLisandro Dalcin   const char         *fname;
1300cc843e7aSLisandro Dalcin   char                bname[PETSC_MAX_PATH_LEN], *gz;
1301cc843e7aSLisandro Dalcin   PetscBool           found;
1302cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
13035c6c1daeSBarry Smith 
1304cc843e7aSLisandro Dalcin   PetscFunctionBegin;
13059566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
13065c6c1daeSBarry Smith 
1307cc843e7aSLisandro Dalcin   /* if file name ends in .gz strip that off and note user wants file compressed */
1308cc843e7aSLisandro Dalcin   vbinary->storecompressed = PETSC_FALSE;
1309cc843e7aSLisandro Dalcin   if (vbinary->filemode == FILE_MODE_WRITE) {
13109566063dSJacob Faibussowitsch     PetscCall(PetscStrstr(vbinary->filename, ".gz", &gz));
13119371c9d4SSatish Balay     if (gz && gz[3] == 0) {
13129371c9d4SSatish Balay       *gz                      = 0;
13139371c9d4SSatish Balay       vbinary->storecompressed = PETSC_TRUE;
13149371c9d4SSatish Balay     }
1315cc843e7aSLisandro Dalcin   }
1316cc843e7aSLisandro Dalcin #if !defined(PETSC_HAVE_POPEN)
131728b400f6SJacob Faibussowitsch   PetscCheck(!vbinary->storecompressed, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP_SYS, "Cannot run gzip on this machine");
1318cc843e7aSLisandro Dalcin #endif
1319cc843e7aSLisandro Dalcin 
1320cc843e7aSLisandro Dalcin   fname = vbinary->filename;
1321cc843e7aSLisandro Dalcin   if (vbinary->filemode == FILE_MODE_READ) { /* possibly get the file from remote site or compressed file */
13229566063dSJacob Faibussowitsch     PetscCall(PetscFileRetrieve(PetscObjectComm((PetscObject)viewer), fname, bname, PETSC_MAX_PATH_LEN, &found));
132328b400f6SJacob Faibussowitsch     PetscCheck(found, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_OPEN, "Cannot locate file: %s", fname);
1324cc843e7aSLisandro Dalcin     fname = bname;
13255c6c1daeSBarry Smith   }
13265c6c1daeSBarry Smith 
1327cc843e7aSLisandro Dalcin   vbinary->fdes = -1;
1328dd400576SPatrick Sanan   if (rank == 0) { /* only first processor opens file*/
1329cc843e7aSLisandro Dalcin     PetscFileMode mode = vbinary->filemode;
1330cc843e7aSLisandro Dalcin     if (mode == FILE_MODE_APPEND) {
1331cc843e7aSLisandro Dalcin       /* check if asked to append to a non-existing file */
13329566063dSJacob Faibussowitsch       PetscCall(PetscTestFile(fname, '\0', &found));
1333cc843e7aSLisandro Dalcin       if (!found) mode = FILE_MODE_WRITE;
1334cc843e7aSLisandro Dalcin     }
13359566063dSJacob Faibussowitsch     PetscCall(PetscBinaryOpen(fname, mode, &vbinary->fdes));
1336cc843e7aSLisandro Dalcin   }
1337cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1338cc843e7aSLisandro Dalcin }
1339cc843e7aSLisandro Dalcin 
13409371c9d4SSatish Balay static PetscErrorCode PetscViewerFileSetUp_BinaryInfo(PetscViewer viewer) {
1341cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1342cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
1343cc843e7aSLisandro Dalcin   PetscBool           found;
1344cc843e7aSLisandro Dalcin 
1345cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1346cc843e7aSLisandro Dalcin   vbinary->fdes_info = NULL;
13479566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
1348dd400576SPatrick Sanan   if (!vbinary->skipinfo && (vbinary->filemode == FILE_MODE_READ || rank == 0)) {
1349cc843e7aSLisandro Dalcin     char infoname[PETSC_MAX_PATH_LEN], iname[PETSC_MAX_PATH_LEN], *gz;
1350cc843e7aSLisandro Dalcin 
13519566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(infoname, vbinary->filename, sizeof(infoname)));
1352cc843e7aSLisandro Dalcin     /* remove .gz if it ends file name */
13539566063dSJacob Faibussowitsch     PetscCall(PetscStrstr(infoname, ".gz", &gz));
1354cc843e7aSLisandro Dalcin     if (gz && gz[3] == 0) *gz = 0;
1355cc843e7aSLisandro Dalcin 
13569566063dSJacob Faibussowitsch     PetscCall(PetscStrlcat(infoname, ".info", sizeof(infoname)));
1357cc843e7aSLisandro Dalcin     if (vbinary->filemode == FILE_MODE_READ) {
13589566063dSJacob Faibussowitsch       PetscCall(PetscFixFilename(infoname, iname));
13599566063dSJacob Faibussowitsch       PetscCall(PetscFileRetrieve(PetscObjectComm((PetscObject)viewer), iname, infoname, PETSC_MAX_PATH_LEN, &found));
13609566063dSJacob Faibussowitsch       if (found) PetscCall(PetscOptionsInsertFile(PetscObjectComm((PetscObject)viewer), ((PetscObject)viewer)->options, infoname, PETSC_FALSE));
1361dd400576SPatrick Sanan     } else if (rank == 0) { /* write or append */
1362cc843e7aSLisandro Dalcin       const char *omode  = (vbinary->filemode == FILE_MODE_APPEND) ? "a" : "w";
1363cc843e7aSLisandro Dalcin       vbinary->fdes_info = fopen(infoname, omode);
136428b400f6SJacob Faibussowitsch       PetscCheck(vbinary->fdes_info, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open .info file %s for writing", infoname);
13655c6c1daeSBarry Smith     }
13665c6c1daeSBarry Smith   }
13675c6c1daeSBarry Smith   PetscFunctionReturn(0);
13685c6c1daeSBarry Smith }
13695c6c1daeSBarry Smith 
13709371c9d4SSatish Balay static PetscErrorCode PetscViewerSetUp_Binary(PetscViewer viewer) {
13715c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1372cc843e7aSLisandro Dalcin   PetscBool           usempiio;
1373cc843e7aSLisandro Dalcin 
13745c6c1daeSBarry Smith   PetscFunctionBegin;
13759566063dSJacob Faibussowitsch   if (!vbinary->setfromoptionscalled) PetscCall(PetscViewerSetFromOptions(viewer));
137628b400f6SJacob Faibussowitsch   PetscCheck(vbinary->filename, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetName()");
137708401ef6SPierre Jolivet   PetscCheck(vbinary->filemode != (PetscFileMode)-1, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode()");
13789566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_Binary(viewer));
1379cc843e7aSLisandro Dalcin 
13809566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetUseMPIIO(viewer, &usempiio));
1381cc843e7aSLisandro Dalcin   if (usempiio) {
1382cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
13839566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetUp_BinaryMPIIO(viewer));
1384cc843e7aSLisandro Dalcin #endif
1385cc843e7aSLisandro Dalcin   } else {
13869566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetUp_BinarySTDIO(viewer));
1387cc843e7aSLisandro Dalcin   }
13889566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetUp_BinaryInfo(viewer));
1389cc843e7aSLisandro Dalcin 
13909566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectState((PetscObject)viewer, "File: %s", vbinary->filename));
13915c6c1daeSBarry Smith   PetscFunctionReturn(0);
13925c6c1daeSBarry Smith }
13935c6c1daeSBarry Smith 
13949371c9d4SSatish Balay static PetscErrorCode PetscViewerView_Binary(PetscViewer v, PetscViewer viewer) {
1395cb6ad94fSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
1396cb6ad94fSLisandro Dalcin   const char         *fname   = vbinary->filename ? vbinary->filename : "not yet set";
1397cc843e7aSLisandro Dalcin   const char         *fmode   = vbinary->filemode != (PetscFileMode)-1 ? PetscFileModes[vbinary->filemode] : "not yet set";
1398cb6ad94fSLisandro Dalcin   PetscBool           usempiio;
13992bf49c77SBarry Smith 
14002bf49c77SBarry Smith   PetscFunctionBegin;
14019566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetUseMPIIO(v, &usempiio));
14029566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Filename: %s\n", fname));
14039566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Mode: %s (%s)\n", fmode, usempiio ? "mpiio" : "stdio"));
14042bf49c77SBarry Smith   PetscFunctionReturn(0);
14052bf49c77SBarry Smith }
14062bf49c77SBarry Smith 
14079371c9d4SSatish Balay static PetscErrorCode PetscViewerSetFromOptions_Binary(PetscViewer viewer, PetscOptionItems *PetscOptionsObject) {
140822a8f86cSLisandro Dalcin   PetscViewer_Binary *binary = (PetscViewer_Binary *)viewer->data;
1409d22fd6bcSDave May   char                defaultname[PETSC_MAX_PATH_LEN];
141003a1d59fSDave May   PetscBool           flg;
1411e0385b85SDave May 
141203a1d59fSDave May   PetscFunctionBegin;
141322a8f86cSLisandro Dalcin   if (viewer->setupcalled) PetscFunctionReturn(0);
1414d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Binary PetscViewer Options");
14159566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(defaultname, PETSC_MAX_PATH_LEN - 1, "binaryoutput"));
14169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-viewer_binary_filename", "Specify filename", "PetscViewerFileSetName", defaultname, defaultname, sizeof(defaultname), &flg));
14179566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscViewerFileSetName_Binary(viewer, defaultname));
14189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_skip_info", "Skip writing/reading .info file", "PetscViewerBinarySetSkipInfo", binary->skipinfo, &binary->skipinfo, NULL));
14199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_skip_options", "Skip parsing Vec/Mat load options", "PetscViewerBinarySetSkipOptions", binary->skipoptions, &binary->skipoptions, NULL));
14209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_skip_header", "Skip writing/reading header information", "PetscViewerBinarySetSkipHeader", binary->skipheader, &binary->skipheader, NULL));
142103a1d59fSDave May #if defined(PETSC_HAVE_MPIIO)
14229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_mpiio", "Use MPI-IO functionality to write/read binary file", "PetscViewerBinarySetUseMPIIO", binary->usempiio, &binary->usempiio, NULL));
14235972f5f3SLisandro Dalcin #else
14249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_mpiio", "Use MPI-IO functionality to write/read binary file (NOT AVAILABLE)", "PetscViewerBinarySetUseMPIIO", PETSC_FALSE, NULL, NULL));
142503a1d59fSDave May #endif
1426d0609cedSBarry Smith   PetscOptionsHeadEnd();
1427bc196f7cSDave May   binary->setfromoptionscalled = PETSC_TRUE;
142803a1d59fSDave May   PetscFunctionReturn(0);
142903a1d59fSDave May }
143003a1d59fSDave May 
14318556b5ebSBarry Smith /*MC
14328556b5ebSBarry Smith    PETSCVIEWERBINARY - A viewer that saves to binary files
14338556b5ebSBarry Smith 
14341b266c99SBarry Smith   Level: beginner
14351b266c99SBarry Smith 
1436811af0c4SBarry Smith .seealso: `PetscViewerBinaryOpen()`, `PETSC_VIEWER_STDOUT_()`, `PETSC_VIEWER_STDOUT_SELF`, `PETSC_VIEWER_STDOUT_WORLD`, `PetscViewerCreate()`, `PetscViewerASCIIOpen()`,
1437811af0c4SBarry Smith           `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`, `PETSCVIEWERDRAW`, `PETSCVIEWERSOCKET`
1438811af0c4SBarry Smith           `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()`,
1439811af0c4SBarry Smith           `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`
14408556b5ebSBarry Smith M*/
14418556b5ebSBarry Smith 
14429371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscViewerCreate_Binary(PetscViewer v) {
14435c6c1daeSBarry Smith   PetscViewer_Binary *vbinary;
14445c6c1daeSBarry Smith 
14455c6c1daeSBarry Smith   PetscFunctionBegin;
1446*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&vbinary));
14475c6c1daeSBarry Smith   v->data = (void *)vbinary;
1448cc843e7aSLisandro Dalcin 
144903a1d59fSDave May   v->ops->setfromoptions   = PetscViewerSetFromOptions_Binary;
14505c6c1daeSBarry Smith   v->ops->destroy          = PetscViewerDestroy_Binary;
14512bf49c77SBarry Smith   v->ops->view             = PetscViewerView_Binary;
1452c98fd787SBarry Smith   v->ops->setup            = PetscViewerSetUp_Binary;
1453cc843e7aSLisandro Dalcin   v->ops->flush            = NULL; /* Should we support Flush() ? */
1454cc843e7aSLisandro Dalcin   v->ops->getsubviewer     = PetscViewerGetSubViewer_Binary;
1455cc843e7aSLisandro Dalcin   v->ops->restoresubviewer = PetscViewerRestoreSubViewer_Binary;
1456cc843e7aSLisandro Dalcin   v->ops->read             = PetscViewerBinaryRead;
1457cc843e7aSLisandro Dalcin 
1458cc843e7aSLisandro Dalcin   vbinary->fdes = -1;
1459e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
1460cc843e7aSLisandro Dalcin   vbinary->usempiio = PETSC_FALSE;
1461e8a65b7dSLisandro Dalcin   vbinary->mfdes    = MPI_FILE_NULL;
1462e8a65b7dSLisandro Dalcin   vbinary->mfsub    = MPI_FILE_NULL;
1463e8a65b7dSLisandro Dalcin #endif
1464cc843e7aSLisandro Dalcin   vbinary->filename        = NULL;
14657e4fd573SVaclav Hapla   vbinary->filemode        = FILE_MODE_UNDEFINED;
146602c9f0b5SLisandro Dalcin   vbinary->fdes_info       = NULL;
14675c6c1daeSBarry Smith   vbinary->skipinfo        = PETSC_FALSE;
14685c6c1daeSBarry Smith   vbinary->skipoptions     = PETSC_TRUE;
14695c6c1daeSBarry Smith   vbinary->skipheader      = PETSC_FALSE;
14705c6c1daeSBarry Smith   vbinary->storecompressed = PETSC_FALSE;
1471f90597f1SStefano Zampini   vbinary->ogzfilename     = NULL;
14725c6c1daeSBarry Smith   vbinary->flowcontrol     = 256; /* seems a good number for Cray XT-5 */
14735c6c1daeSBarry Smith 
1474cc843e7aSLisandro Dalcin   vbinary->setfromoptionscalled = PETSC_FALSE;
1475cc843e7aSLisandro Dalcin 
14769566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetFlowControl_C", PetscViewerBinaryGetFlowControl_Binary));
14779566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetFlowControl_C", PetscViewerBinarySetFlowControl_Binary));
14789566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipHeader_C", PetscViewerBinaryGetSkipHeader_Binary));
14799566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipHeader_C", PetscViewerBinarySetSkipHeader_Binary));
14809566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipOptions_C", PetscViewerBinaryGetSkipOptions_Binary));
14819566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipOptions_C", PetscViewerBinarySetSkipOptions_Binary));
14829566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipInfo_C", PetscViewerBinaryGetSkipInfo_Binary));
14839566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipInfo_C", PetscViewerBinarySetSkipInfo_Binary));
14849566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetInfoPointer_C", PetscViewerBinaryGetInfoPointer_Binary));
14859566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_Binary));
14869566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_Binary));
14879566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_Binary));
14889566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_Binary));
14895c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
14909566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetUseMPIIO_C", PetscViewerBinaryGetUseMPIIO_Binary));
14919566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetUseMPIIO_C", PetscViewerBinarySetUseMPIIO_Binary));
14925c6c1daeSBarry Smith #endif
14935c6c1daeSBarry Smith   PetscFunctionReturn(0);
14945c6c1daeSBarry Smith }
14955c6c1daeSBarry Smith 
14965c6c1daeSBarry Smith /* ---------------------------------------------------------------------*/
14975c6c1daeSBarry Smith /*
14985c6c1daeSBarry Smith     The variable Petsc_Viewer_Binary_keyval is used to indicate an MPI attribute that
14995c6c1daeSBarry Smith   is attached to a communicator, in this case the attribute is a PetscViewer.
15005c6c1daeSBarry Smith */
1501d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_Binary_keyval = MPI_KEYVAL_INVALID;
15025c6c1daeSBarry Smith 
15035c6c1daeSBarry Smith /*@C
1504811af0c4SBarry Smith      PETSC_VIEWER_BINARY_ - Creates a `PETSCVIEWERBINARY` `PetscViewer` shared by all processors
15055c6c1daeSBarry Smith                      in a communicator.
15065c6c1daeSBarry Smith 
1507d083f849SBarry Smith      Collective
15085c6c1daeSBarry Smith 
15095c6c1daeSBarry Smith      Input Parameter:
1510811af0c4SBarry Smith .    comm - the MPI communicator to share the `PETSCVIEWERBINARY`
15115c6c1daeSBarry Smith 
15125c6c1daeSBarry Smith      Level: intermediate
15135c6c1daeSBarry Smith 
15145c6c1daeSBarry Smith    Options Database Keys:
151510699b91SBarry Smith +    -viewer_binary_filename <name> - filename in which to store the binary data, defaults to binaryoutput
151610699b91SBarry Smith .    -viewer_binary_skip_info - true means do not create .info file for this viewer
151710699b91SBarry Smith .    -viewer_binary_skip_options - true means do not use the options database for this viewer
151810699b91SBarry Smith .    -viewer_binary_skip_header - true means do not store the usual header information in the binary file
151910699b91SBarry Smith -    -viewer_binary_mpiio - true means use the file via MPI-IO, maybe faster for large files and many MPI ranks
15205c6c1daeSBarry Smith 
1521811af0c4SBarry Smith    Environmental variable:
152210699b91SBarry Smith -   PETSC_VIEWER_BINARY_FILENAME - filename in which to store the binary data, defaults to binaryoutput
15235c6c1daeSBarry Smith 
1524811af0c4SBarry Smith      Note:
1525811af0c4SBarry Smith      Unlike almost all other PETSc routines, `PETSC_VIEWER_BINARY_` does not return
15265c6c1daeSBarry Smith      an error code.  The binary PetscViewer is usually used in the form
15275c6c1daeSBarry Smith $       XXXView(XXX object,PETSC_VIEWER_BINARY_(comm));
15285c6c1daeSBarry Smith 
1529811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PETSC_VIEWER_BINARY_WORLD`, `PETSC_VIEWER_BINARY_SELF`, `PetscViewerBinaryOpen()`, `PetscViewerCreate()`,
1530db781477SPatrick Sanan           `PetscViewerDestroy()`
15315c6c1daeSBarry Smith @*/
15329371c9d4SSatish Balay PetscViewer PETSC_VIEWER_BINARY_(MPI_Comm comm) {
15335c6c1daeSBarry Smith   PetscErrorCode ierr;
15345c6c1daeSBarry Smith   PetscBool      flg;
15355c6c1daeSBarry Smith   PetscViewer    viewer;
15365c6c1daeSBarry Smith   char           fname[PETSC_MAX_PATH_LEN];
15375c6c1daeSBarry Smith   MPI_Comm       ncomm;
15385c6c1daeSBarry Smith 
15395c6c1daeSBarry Smith   PetscFunctionBegin;
15409371c9d4SSatish Balay   ierr = PetscCommDuplicate(comm, &ncomm, NULL);
15419371c9d4SSatish Balay   if (ierr) {
15429371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
15439371c9d4SSatish Balay     PetscFunctionReturn(NULL);
15449371c9d4SSatish Balay   }
15455c6c1daeSBarry Smith   if (Petsc_Viewer_Binary_keyval == MPI_KEYVAL_INVALID) {
154602c9f0b5SLisandro Dalcin     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Viewer_Binary_keyval, NULL);
15479371c9d4SSatish Balay     if (ierr) {
15489371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
15499371c9d4SSatish Balay       PetscFunctionReturn(NULL);
15509371c9d4SSatish Balay     }
15515c6c1daeSBarry Smith   }
155247435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm, Petsc_Viewer_Binary_keyval, (void **)&viewer, (int *)&flg);
15539371c9d4SSatish Balay   if (ierr) {
15549371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
15559371c9d4SSatish Balay     PetscFunctionReturn(NULL);
15569371c9d4SSatish Balay   }
15575c6c1daeSBarry Smith   if (!flg) { /* PetscViewer not yet created */
15585c6c1daeSBarry Smith     ierr = PetscOptionsGetenv(ncomm, "PETSC_VIEWER_BINARY_FILENAME", fname, PETSC_MAX_PATH_LEN, &flg);
15599371c9d4SSatish Balay     if (ierr) {
15609371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
15619371c9d4SSatish Balay       PetscFunctionReturn(NULL);
15629371c9d4SSatish Balay     }
15635c6c1daeSBarry Smith     if (!flg) {
15645c6c1daeSBarry Smith       ierr = PetscStrcpy(fname, "binaryoutput");
15659371c9d4SSatish Balay       if (ierr) {
15669371c9d4SSatish Balay         PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
15679371c9d4SSatish Balay         PetscFunctionReturn(NULL);
15689371c9d4SSatish Balay       }
15695c6c1daeSBarry Smith     }
15705c6c1daeSBarry Smith     ierr = PetscViewerBinaryOpen(ncomm, fname, FILE_MODE_WRITE, &viewer);
15719371c9d4SSatish Balay     if (ierr) {
15729371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
15739371c9d4SSatish Balay       PetscFunctionReturn(NULL);
15749371c9d4SSatish Balay     }
15755c6c1daeSBarry Smith     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
15769371c9d4SSatish Balay     if (ierr) {
15779371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
15789371c9d4SSatish Balay       PetscFunctionReturn(NULL);
15799371c9d4SSatish Balay     }
158047435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm, Petsc_Viewer_Binary_keyval, (void *)viewer);
15819371c9d4SSatish Balay     if (ierr) {
15829371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
15839371c9d4SSatish Balay       PetscFunctionReturn(NULL);
15849371c9d4SSatish Balay     }
15855c6c1daeSBarry Smith   }
15865c6c1daeSBarry Smith   ierr = PetscCommDestroy(&ncomm);
15879371c9d4SSatish Balay   if (ierr) {
15889371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
15899371c9d4SSatish Balay     PetscFunctionReturn(NULL);
15909371c9d4SSatish Balay   }
15915c6c1daeSBarry Smith   PetscFunctionReturn(viewer);
15925c6c1daeSBarry Smith }
1593