xref: /petsc/src/sys/classes/viewer/impls/binary/binv.c (revision c410d8cc8908fa7552d5d2cadc5c41ff54401973)
1af0996ceSBarry Smith #include <petsc/private/viewerimpl.h> /*I   "petscviewer.h"   I*/
25c6c1daeSBarry Smith 
376667918SBarry Smith /*
476667918SBarry Smith    This needs to start the same as PetscViewer_Socket.
576667918SBarry Smith */
65c6c1daeSBarry Smith typedef struct {
75c6c1daeSBarry Smith   int       fdes;        /* file descriptor, ignored if using MPI IO */
876667918SBarry Smith   PetscInt  flowcontrol; /* allow only <flowcontrol> messages outstanding at a time while doing IO */
976667918SBarry Smith   PetscBool skipheader;  /* don't write header, only raw data */
105c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
11bc196f7cSDave May   PetscBool  usempiio;
125c6c1daeSBarry Smith   MPI_File   mfdes; /* ignored unless using MPI IO */
13e8a65b7dSLisandro Dalcin   MPI_File   mfsub; /* subviewer support */
145c6c1daeSBarry Smith   MPI_Offset moff;
155c6c1daeSBarry Smith #endif
16cc843e7aSLisandro Dalcin   char         *filename;            /* file name */
17cc843e7aSLisandro Dalcin   PetscFileMode filemode;            /* read/write/append mode */
185c6c1daeSBarry Smith   FILE         *fdes_info;           /* optional file containing info on binary file*/
195c6c1daeSBarry Smith   PetscBool     storecompressed;     /* gzip the write binary file when closing it*/
20f90597f1SStefano Zampini   char         *ogzfilename;         /* gzip can be run after the filename has been updated */
215c6c1daeSBarry Smith   PetscBool     skipinfo;            /* Don't create info file for writing; don't use for reading */
225c6c1daeSBarry Smith   PetscBool     skipoptions;         /* don't use PETSc options database when loading */
23a261c58fSBarry Smith   PetscBool     matlabheaderwritten; /* if format is PETSC_VIEWER_BINARY_MATLAB has the MATLAB .info header been written yet */
24c98fd787SBarry Smith   PetscBool     setfromoptionscalled;
255c6c1daeSBarry Smith } PetscViewer_Binary;
265c6c1daeSBarry Smith 
27d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryClearFunctionList(PetscViewer v)
28d71ae5a4SJacob Faibussowitsch {
292e956fe4SStefano Zampini   PetscFunctionBegin;
302e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetFlowControl_C", NULL));
312e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetFlowControl_C", NULL));
322e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipHeader_C", NULL));
332e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipHeader_C", NULL));
342e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipOptions_C", NULL));
352e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipOptions_C", NULL));
362e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipInfo_C", NULL));
372e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipInfo_C", NULL));
382e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetInfoPointer_C", NULL));
392e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", NULL));
402e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", NULL));
412e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", NULL));
422e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", NULL));
432e956fe4SStefano Zampini #if defined(PETSC_HAVE_MPIIO)
442e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetUseMPIIO_C", NULL));
452e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetUseMPIIO_C", NULL));
462e956fe4SStefano Zampini #endif
473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
482e956fe4SStefano Zampini }
492e956fe4SStefano Zampini 
50cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
51d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySyncMPIIO(PetscViewer viewer)
52d71ae5a4SJacob Faibussowitsch {
53cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
54cc843e7aSLisandro Dalcin 
55cc843e7aSLisandro Dalcin   PetscFunctionBegin;
563ba16761SJacob Faibussowitsch   if (vbinary->filemode == FILE_MODE_READ) PetscFunctionReturn(PETSC_SUCCESS);
5748a46eb9SPierre Jolivet   if (vbinary->mfsub != MPI_FILE_NULL) PetscCallMPI(MPI_File_sync(vbinary->mfsub));
58cc843e7aSLisandro Dalcin   if (vbinary->mfdes != MPI_FILE_NULL) {
599566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer)));
609566063dSJacob Faibussowitsch     PetscCallMPI(MPI_File_sync(vbinary->mfdes));
61cc843e7aSLisandro Dalcin   }
623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
63cc843e7aSLisandro Dalcin }
64cc843e7aSLisandro Dalcin #endif
65cc843e7aSLisandro Dalcin 
66d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerGetSubViewer_Binary(PetscViewer viewer, MPI_Comm comm, PetscViewer *outviewer)
67d71ae5a4SJacob Faibussowitsch {
68e8a65b7dSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
69cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
705c6c1daeSBarry Smith 
715c6c1daeSBarry Smith   PetscFunctionBegin;
729566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
73e8a65b7dSLisandro Dalcin 
74e8a65b7dSLisandro Dalcin   /* Return subviewer in process zero */
759566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
76dd400576SPatrick Sanan   if (rank == 0) {
77e5afcf28SBarry Smith     PetscMPIInt flg;
78e5afcf28SBarry Smith 
799566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_compare(PETSC_COMM_SELF, comm, &flg));
80cc73adaaSBarry Smith     PetscCheck(flg == MPI_IDENT || flg == MPI_CONGRUENT, PETSC_COMM_SELF, PETSC_ERR_SUP, "PetscViewerGetSubViewer() for PETSCVIEWERBINARY requires a singleton MPI_Comm");
819566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(comm, outviewer));
829566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(*outviewer, PETSCVIEWERBINARY));
839566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy((*outviewer)->data, vbinary, sizeof(PetscViewer_Binary)));
8412f4c3a9SLisandro Dalcin     (*outviewer)->setupcalled = PETSC_TRUE;
8596509d17SLisandro Dalcin   } else {
8696509d17SLisandro Dalcin     *outviewer = NULL;
8796509d17SLisandro Dalcin   }
88e8a65b7dSLisandro Dalcin 
89e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
90e8a65b7dSLisandro Dalcin   if (vbinary->usempiio && *outviewer) {
91e8a65b7dSLisandro Dalcin     PetscViewer_Binary *obinary = (PetscViewer_Binary *)(*outviewer)->data;
92e8a65b7dSLisandro Dalcin     /* Parent viewer opens a new MPI file handle on PETSC_COMM_SELF and keeps track of it for future reuse */
93e8a65b7dSLisandro Dalcin     if (vbinary->mfsub == MPI_FILE_NULL) {
94e8a65b7dSLisandro Dalcin       int amode;
95cc843e7aSLisandro Dalcin       switch (vbinary->filemode) {
96d71ae5a4SJacob Faibussowitsch       case FILE_MODE_READ:
97d71ae5a4SJacob Faibussowitsch         amode = MPI_MODE_RDONLY;
98d71ae5a4SJacob Faibussowitsch         break;
99d71ae5a4SJacob Faibussowitsch       case FILE_MODE_WRITE:
100d71ae5a4SJacob Faibussowitsch         amode = MPI_MODE_WRONLY;
101d71ae5a4SJacob Faibussowitsch         break;
102d71ae5a4SJacob Faibussowitsch       case FILE_MODE_APPEND:
103d71ae5a4SJacob Faibussowitsch         amode = MPI_MODE_WRONLY;
104d71ae5a4SJacob Faibussowitsch         break;
105d71ae5a4SJacob Faibussowitsch       default:
106d71ae5a4SJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[vbinary->filemode]);
107e8a65b7dSLisandro Dalcin       }
1089566063dSJacob Faibussowitsch       PetscCallMPI(MPI_File_open(PETSC_COMM_SELF, vbinary->filename, amode, MPI_INFO_NULL, &vbinary->mfsub));
109e8a65b7dSLisandro Dalcin     }
110e8a65b7dSLisandro Dalcin     /* Subviewer gets the MPI file handle on PETSC_COMM_SELF */
111e8a65b7dSLisandro Dalcin     obinary->mfdes = vbinary->mfsub;
112e8a65b7dSLisandro Dalcin     obinary->mfsub = MPI_FILE_NULL;
113e8a65b7dSLisandro Dalcin     obinary->moff  = vbinary->moff;
114e8a65b7dSLisandro Dalcin   }
115e8a65b7dSLisandro Dalcin #endif
116cc843e7aSLisandro Dalcin 
117cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
1189566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinarySyncMPIIO(viewer));
119cc843e7aSLisandro Dalcin #endif
1203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1215c6c1daeSBarry Smith }
1225c6c1daeSBarry Smith 
123d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerRestoreSubViewer_Binary(PetscViewer viewer, MPI_Comm comm, PetscViewer *outviewer)
124d71ae5a4SJacob Faibussowitsch {
125e8a65b7dSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
126cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
127e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
128e8a65b7dSLisandro Dalcin   MPI_Offset moff = 0;
129e8a65b7dSLisandro Dalcin #endif
1305c6c1daeSBarry Smith 
1315c6c1daeSBarry Smith   PetscFunctionBegin;
1329566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
133c5853193SPierre Jolivet   PetscCheck(rank == 0 || !*outviewer, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Subviewer not obtained from viewer");
134e8a65b7dSLisandro Dalcin 
135e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
136e8a65b7dSLisandro Dalcin   if (vbinary->usempiio && *outviewer) {
137e8a65b7dSLisandro Dalcin     PetscViewer_Binary *obinary = (PetscViewer_Binary *)(*outviewer)->data;
13808401ef6SPierre Jolivet     PetscCheck(obinary->mfdes == vbinary->mfsub, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Subviewer not obtained from viewer");
1399566063dSJacob Faibussowitsch     if (obinary->mfsub != MPI_FILE_NULL) PetscCallMPI(MPI_File_close(&obinary->mfsub));
140e8a65b7dSLisandro Dalcin     moff = obinary->moff;
141e8a65b7dSLisandro Dalcin   }
142e8a65b7dSLisandro Dalcin #endif
143e8a65b7dSLisandro Dalcin 
144e8a65b7dSLisandro Dalcin   if (*outviewer) {
145e8a65b7dSLisandro Dalcin     PetscViewer_Binary *obinary = (PetscViewer_Binary *)(*outviewer)->data;
14608401ef6SPierre Jolivet     PetscCheck(obinary->fdes == vbinary->fdes, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Subviewer not obtained from viewer");
1479566063dSJacob Faibussowitsch     PetscCall(PetscFree((*outviewer)->data));
1482e956fe4SStefano Zampini     PetscCall(PetscViewerBinaryClearFunctionList(*outviewer));
1499566063dSJacob Faibussowitsch     PetscCall(PetscHeaderDestroy(outviewer));
1505c6c1daeSBarry Smith   }
151e8a65b7dSLisandro Dalcin 
152e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
153e8a65b7dSLisandro Dalcin   if (vbinary->usempiio) {
154e8a65b7dSLisandro Dalcin     PetscInt64 ioff = (PetscInt64)moff; /* We could use MPI_OFFSET datatype (requires MPI 2.2) */
1559566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&ioff, 1, MPIU_INT64, 0, PetscObjectComm((PetscObject)viewer)));
156e8a65b7dSLisandro Dalcin     vbinary->moff = (MPI_Offset)ioff;
157e8a65b7dSLisandro Dalcin   }
158e8a65b7dSLisandro Dalcin #endif
159cc843e7aSLisandro Dalcin 
160cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
1619566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinarySyncMPIIO(viewer));
162cc843e7aSLisandro Dalcin #endif
1633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1645c6c1daeSBarry Smith }
1655c6c1daeSBarry Smith 
1665c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
1675c6c1daeSBarry Smith /*@C
16885b8072bSPatrick Sanan     PetscViewerBinaryGetMPIIOOffset - Gets the current global offset that should be passed to `MPI_File_set_view()` or `MPI_File_{write|read}_at[_all]()`
1695c6c1daeSBarry Smith 
170cf53795eSBarry Smith     Not Collective; No Fortran Support
1715c6c1daeSBarry Smith 
1725c6c1daeSBarry Smith     Input Parameter:
173*c410d8ccSBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
1745c6c1daeSBarry Smith 
1755c6c1daeSBarry Smith     Output Parameter:
17622a8f86cSLisandro Dalcin .   off - the current global offset
1775c6c1daeSBarry Smith 
1785c6c1daeSBarry Smith     Level: advanced
1795c6c1daeSBarry Smith 
180811af0c4SBarry Smith     Note:
181811af0c4SBarry Smith     Use `PetscViewerBinaryAddMPIIOOffset()` to increase this value after you have written a view.
182811af0c4SBarry Smith 
183d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryAddMPIIOOffset()`
1845c6c1daeSBarry Smith @*/
185d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetMPIIOOffset(PetscViewer viewer, MPI_Offset *off)
186d71ae5a4SJacob Faibussowitsch {
18722a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
1885c6c1daeSBarry Smith 
1895c6c1daeSBarry Smith   PetscFunctionBegin;
19022a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
19122a8f86cSLisandro Dalcin   PetscValidPointer(off, 2);
19222a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
1935c6c1daeSBarry Smith   *off    = vbinary->moff;
1943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1955c6c1daeSBarry Smith }
1965c6c1daeSBarry Smith 
1975c6c1daeSBarry Smith /*@C
19822a8f86cSLisandro Dalcin     PetscViewerBinaryAddMPIIOOffset - Adds to the current global offset
1995c6c1daeSBarry Smith 
200cf53795eSBarry Smith     Logically Collective; No Fortran Support
2015c6c1daeSBarry Smith 
2025c6c1daeSBarry Smith     Input Parameters:
203811af0c4SBarry Smith +   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
20422a8f86cSLisandro Dalcin -   off - the addition to the global offset
2055c6c1daeSBarry Smith 
2065c6c1daeSBarry Smith     Level: advanced
2075c6c1daeSBarry Smith 
208811af0c4SBarry Smith     Note:
209811af0c4SBarry Smith     Use `PetscViewerBinaryGetMPIIOOffset()` to get the value that you should pass to `MPI_File_set_view()` or `MPI_File_{write|read}_at[_all]()`
210811af0c4SBarry Smith 
211d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
2125c6c1daeSBarry Smith @*/
213d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryAddMPIIOOffset(PetscViewer viewer, MPI_Offset off)
214d71ae5a4SJacob Faibussowitsch {
21522a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
2165c6c1daeSBarry Smith 
2175c6c1daeSBarry Smith   PetscFunctionBegin;
21822a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
21922a8f86cSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, (PetscInt)off, 2);
22022a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
2215c6c1daeSBarry Smith   vbinary->moff += off;
2223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2235c6c1daeSBarry Smith }
2245c6c1daeSBarry Smith 
2255c6c1daeSBarry Smith /*@C
226811af0c4SBarry Smith     PetscViewerBinaryGetMPIIODescriptor - Extracts the MPI IO file descriptor from a `PetscViewer`.
2275c6c1daeSBarry Smith 
228cf53795eSBarry Smith     Not Collective; No Fortran Support
2295c6c1daeSBarry Smith 
2305c6c1daeSBarry Smith     Input Parameter:
231811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
2325c6c1daeSBarry Smith 
2335c6c1daeSBarry Smith     Output Parameter:
2345c6c1daeSBarry Smith .   fdes - file descriptor
2355c6c1daeSBarry Smith 
2365c6c1daeSBarry Smith     Level: advanced
2375c6c1daeSBarry Smith 
238d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
2395c6c1daeSBarry Smith @*/
240d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetMPIIODescriptor(PetscViewer viewer, MPI_File *fdes)
241d71ae5a4SJacob Faibussowitsch {
24222a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
2435c6c1daeSBarry Smith 
2445c6c1daeSBarry Smith   PetscFunctionBegin;
24522a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
24622a8f86cSLisandro Dalcin   PetscValidPointer(fdes, 2);
2479566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
24822a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
2495c6c1daeSBarry Smith   *fdes   = vbinary->mfdes;
2503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2515c6c1daeSBarry Smith }
252cc843e7aSLisandro Dalcin #endif
2535c6c1daeSBarry Smith 
254cc843e7aSLisandro Dalcin /*@
255cc843e7aSLisandro Dalcin     PetscViewerBinarySetUseMPIIO - Sets a binary viewer to use MPI-IO for reading/writing. Must be called
256811af0c4SBarry Smith         before `PetscViewerFileSetName()`
257cc843e7aSLisandro Dalcin 
258c3339decSBarry Smith     Logically Collective
259cc843e7aSLisandro Dalcin 
260cc843e7aSLisandro Dalcin     Input Parameters:
261811af0c4SBarry Smith +   viewer - the `PetscViewer`; must be a `PETSCVIEWERBINARY`
262811af0c4SBarry Smith -   use - `PETSC_TRUE` means MPI-IO will be used
263cc843e7aSLisandro Dalcin 
264811af0c4SBarry Smith     Options Database Key:
265cc843e7aSLisandro Dalcin     -viewer_binary_mpiio : Flag for using MPI-IO
266cc843e7aSLisandro Dalcin 
267cc843e7aSLisandro Dalcin     Level: advanced
268cc843e7aSLisandro Dalcin 
269d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`,
270db781477SPatrick Sanan           `PetscViewerBinaryGetUseMPIIO()`
271cc843e7aSLisandro Dalcin @*/
272d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetUseMPIIO(PetscViewer viewer, PetscBool use)
273d71ae5a4SJacob Faibussowitsch {
274a8aae444SDave May   PetscFunctionBegin;
275cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
276cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, use, 2);
277cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetUseMPIIO_C", (PetscViewer, PetscBool), (viewer, use));
2783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
279cc843e7aSLisandro Dalcin }
280cc843e7aSLisandro Dalcin 
281cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
282d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetUseMPIIO_Binary(PetscViewer viewer, PetscBool use)
283d71ae5a4SJacob Faibussowitsch {
284cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
285cc843e7aSLisandro Dalcin   PetscFunctionBegin;
286cc73adaaSBarry Smith   PetscCheck(!viewer->setupcalled || vbinary->usempiio == use, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Cannot change MPIIO to %s after setup", PetscBools[use]);
287cc843e7aSLisandro Dalcin   vbinary->usempiio = use;
2883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
289a8aae444SDave May }
290a8aae444SDave May #endif
291a8aae444SDave May 
292cc843e7aSLisandro Dalcin /*@
2933f423023SBarry Smith     PetscViewerBinaryGetUseMPIIO - Returns `PETSC_TRUE` if the binary viewer uses MPI-IO.
2945c6c1daeSBarry Smith 
2955c6c1daeSBarry Smith     Not Collective
2965c6c1daeSBarry Smith 
2975c6c1daeSBarry Smith     Input Parameter:
298811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`; must be a `PETSCVIEWERBINARY`
2995c6c1daeSBarry Smith 
3005c6c1daeSBarry Smith     Output Parameter:
301811af0c4SBarry Smith .   use - `PETSC_TRUE` if MPI-IO is being used
3025c6c1daeSBarry Smith 
3035c6c1daeSBarry Smith     Level: advanced
3045c6c1daeSBarry Smith 
305bc196f7cSDave May     Note:
3063f423023SBarry Smith     If MPI-IO is not available, this function will always return `PETSC_FALSE`
307bc196f7cSDave May 
308d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
3095c6c1daeSBarry Smith @*/
310d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetUseMPIIO(PetscViewer viewer, PetscBool *use)
311d71ae5a4SJacob Faibussowitsch {
3125c6c1daeSBarry Smith   PetscFunctionBegin;
313cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
314cc843e7aSLisandro Dalcin   PetscValidBoolPointer(use, 2);
315cc843e7aSLisandro Dalcin   *use = PETSC_FALSE;
316cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinaryGetUseMPIIO_C", (PetscViewer, PetscBool *), (viewer, use));
3173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3185c6c1daeSBarry Smith }
3195c6c1daeSBarry Smith 
320cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
321d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetUseMPIIO_Binary(PetscViewer viewer, PetscBool *use)
322d71ae5a4SJacob Faibussowitsch {
3235c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
3245c6c1daeSBarry Smith 
3255c6c1daeSBarry Smith   PetscFunctionBegin;
326cc843e7aSLisandro Dalcin   *use = vbinary->usempiio;
3273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
328cc843e7aSLisandro Dalcin }
329cc843e7aSLisandro Dalcin #endif
330cc843e7aSLisandro Dalcin 
331cc843e7aSLisandro Dalcin /*@
332*c410d8ccSBarry Smith     PetscViewerBinarySetFlowControl - Sets how many messages are allowed to be outstanding at the same time during parallel IO reads/writes
333cc843e7aSLisandro Dalcin 
334cc843e7aSLisandro Dalcin     Not Collective
335cc843e7aSLisandro Dalcin 
336d8d19677SJose E. Roman     Input Parameters:
337*c410d8ccSBarry Smith +   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
338cc843e7aSLisandro Dalcin -   fc - the number of messages, defaults to 256 if this function was not called
339cc843e7aSLisandro Dalcin 
340cc843e7aSLisandro Dalcin     Level: advanced
341cc843e7aSLisandro Dalcin 
342d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetFlowControl()`
343cc843e7aSLisandro Dalcin @*/
344d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetFlowControl(PetscViewer viewer, PetscInt fc)
345d71ae5a4SJacob Faibussowitsch {
346cc843e7aSLisandro Dalcin   PetscFunctionBegin;
347cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
348cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, fc, 2);
349cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetFlowControl_C", (PetscViewer, PetscInt), (viewer, fc));
3503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3515c6c1daeSBarry Smith }
3525c6c1daeSBarry Smith 
353d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetFlowControl_Binary(PetscViewer viewer, PetscInt fc)
354d71ae5a4SJacob Faibussowitsch {
355cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
356cc843e7aSLisandro Dalcin 
357cc843e7aSLisandro Dalcin   PetscFunctionBegin;
35808401ef6SPierre Jolivet   PetscCheck(fc > 1, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_OUTOFRANGE, "Flow control count must be greater than 1, %" PetscInt_FMT " was set", fc);
359cc843e7aSLisandro Dalcin   vbinary->flowcontrol = fc;
3603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
361cc843e7aSLisandro Dalcin }
362cc843e7aSLisandro Dalcin 
363cc843e7aSLisandro Dalcin /*@
364*c410d8ccSBarry Smith     PetscViewerBinaryGetFlowControl - Returns how many messages are allowed to be outstanding at the same time during parallel IO reads/writes
3655c6c1daeSBarry Smith 
3665c6c1daeSBarry Smith     Not Collective
3675c6c1daeSBarry Smith 
3685c6c1daeSBarry Smith     Input Parameter:
369811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
3705c6c1daeSBarry Smith 
3715c6c1daeSBarry Smith     Output Parameter:
3725c6c1daeSBarry Smith .   fc - the number of messages
3735c6c1daeSBarry Smith 
3745c6c1daeSBarry Smith     Level: advanced
3755c6c1daeSBarry Smith 
376d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinarySetFlowControl()`
3775c6c1daeSBarry Smith @*/
378d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetFlowControl(PetscViewer viewer, PetscInt *fc)
379d71ae5a4SJacob Faibussowitsch {
3805c6c1daeSBarry Smith   PetscFunctionBegin;
381cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
382cc843e7aSLisandro Dalcin   PetscValidIntPointer(fc, 2);
383cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetFlowControl_C", (PetscViewer, PetscInt *), (viewer, fc));
3843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3855c6c1daeSBarry Smith }
3865c6c1daeSBarry Smith 
38776667918SBarry Smith PETSC_INTERN PetscErrorCode PetscViewerBinaryGetFlowControl_Binary(PetscViewer viewer, PetscInt *fc)
388d71ae5a4SJacob Faibussowitsch {
3895c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
3905c6c1daeSBarry Smith 
3915c6c1daeSBarry Smith   PetscFunctionBegin;
392cc843e7aSLisandro Dalcin   *fc = vbinary->flowcontrol;
3933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3945c6c1daeSBarry Smith }
3955c6c1daeSBarry Smith 
3965c6c1daeSBarry Smith /*@C
397811af0c4SBarry Smith     PetscViewerBinaryGetDescriptor - Extracts the file descriptor from a `PetscViewer` of `PetscViewerType` `PETSCVIEWERBINARY`.
3985c6c1daeSBarry Smith 
399cd05f99aSJacob Faibussowitsch     Collective because it may trigger a `PetscViewerSetUp()` call; No Fortran Support
4005c6c1daeSBarry Smith 
4015c6c1daeSBarry Smith     Input Parameter:
402811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
4035c6c1daeSBarry Smith 
4045c6c1daeSBarry Smith     Output Parameter:
4055c6c1daeSBarry Smith .   fdes - file descriptor
4065c6c1daeSBarry Smith 
4075c6c1daeSBarry Smith     Level: advanced
4085c6c1daeSBarry Smith 
409811af0c4SBarry Smith     Note:
410811af0c4SBarry Smith       For writable binary `PetscViewer`s, the descriptor will only be valid for the
411811af0c4SBarry Smith     first processor in the communicator that shares the `PetscViewer`. For readable
412*c410d8ccSBarry Smith     files it will only be valid on processes that have the file. If MPI rank 0 does not
413*c410d8ccSBarry Smith     have the file it generates an error even if another MPI process does have the file.
4145c6c1daeSBarry Smith 
415d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`
4165c6c1daeSBarry Smith @*/
417d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetDescriptor(PetscViewer viewer, int *fdes)
418d71ae5a4SJacob Faibussowitsch {
41922a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
4205c6c1daeSBarry Smith 
4215c6c1daeSBarry Smith   PetscFunctionBegin;
42222a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
42322a8f86cSLisandro Dalcin   PetscValidPointer(fdes, 2);
4249566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
42522a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
4265c6c1daeSBarry Smith   *fdes   = vbinary->fdes;
4273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4285c6c1daeSBarry Smith }
4295c6c1daeSBarry Smith 
4305c6c1daeSBarry Smith /*@
431*c410d8ccSBarry Smith     PetscViewerBinarySkipInfo - Binary file will not have `.info` file created with it
4325c6c1daeSBarry Smith 
4335c6c1daeSBarry Smith     Not Collective
4345c6c1daeSBarry Smith 
435fd292e60Sprj-     Input Parameter:
436811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerCreate()`
4375c6c1daeSBarry Smith 
4385c6c1daeSBarry Smith     Options Database Key:
439*c410d8ccSBarry Smith .   -viewer_binary_skip_info - true indicates do not generate `.info` file
4405c6c1daeSBarry Smith 
4415c6c1daeSBarry Smith     Level: advanced
4425c6c1daeSBarry Smith 
44395452b02SPatrick Sanan     Notes:
444811af0c4SBarry Smith     This must be called after `PetscViewerSetType()`. If you use `PetscViewerBinaryOpen()` then
4453f423023SBarry Smith     you can only skip the info file with the `-viewer_binary_skip_info` flag. To use the function you must open the
446811af0c4SBarry Smith     viewer with `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinarySkipInfo()`.
4475c6c1daeSBarry Smith 
448*c410d8ccSBarry Smith     The `.info` files contains meta information about the data in the binary file, for example the block size if it was
4495c6c1daeSBarry Smith     set for a vector or matrix.
4505c6c1daeSBarry Smith 
451811af0c4SBarry Smith     This routine is deprecated, use `PetscViewerBinarySetSkipInfo()`
452811af0c4SBarry Smith 
453d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySetSkipOptions()`,
454db781477SPatrick Sanan           `PetscViewerBinaryGetSkipOptions()`, `PetscViewerBinaryGetSkipInfo()`
4555c6c1daeSBarry Smith @*/
456d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySkipInfo(PetscViewer viewer)
457d71ae5a4SJacob Faibussowitsch {
4585c6c1daeSBarry Smith   PetscFunctionBegin;
4599566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinarySetSkipInfo(viewer, PETSC_TRUE));
4603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
461807ea322SDave May }
462807ea322SDave May 
463807ea322SDave May /*@
464*c410d8ccSBarry Smith     PetscViewerBinarySetSkipInfo - Binary file will not have `.info` file created with it
465807ea322SDave May 
466807ea322SDave May     Not Collective
467807ea322SDave May 
468d8d19677SJose E. Roman     Input Parameters:
4693f423023SBarry Smith +   viewer - PetscViewer context, obtained from `PetscViewerCreate()`
470*c410d8ccSBarry Smith -   skip - `PETSC_TRUE` implies the `.info` file will not be generated
471807ea322SDave May 
472807ea322SDave May     Options Database Key:
473*c410d8ccSBarry Smith .   -viewer_binary_skip_info - true indicates do not generate `.info` file
474807ea322SDave May 
475807ea322SDave May     Level: advanced
476807ea322SDave May 
477*c410d8ccSBarry Smith     Notes:
478*c410d8ccSBarry Smith     This must be called after `PetscViewerSetType()`. If you use `PetscViewerBinaryOpen()` then
479*c410d8ccSBarry Smith     you can only skip the info file with the `-viewer_binary_skip_info` flag. To use the function you must open the
480*c410d8ccSBarry Smith     viewer with `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinarySkipInfo()`.
481*c410d8ccSBarry Smith 
482*c410d8ccSBarry Smith     The `.info` file contains meta information about the data in the binary file, for example the block size if it was
483*c410d8ccSBarry Smith     set for a vector or matrix.
484*c410d8ccSBarry Smith 
485d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySetSkipOptions()`,
486*c410d8ccSBarry Smith           `PetscViewerBinaryGetSkipOptions()`, `PetscViewerBinaryGetSkipInfo()`, `PetscViewerBinaryGetInfoPointer()`
487807ea322SDave May @*/
488d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetSkipInfo(PetscViewer viewer, PetscBool skip)
489d71ae5a4SJacob Faibussowitsch {
490807ea322SDave May   PetscFunctionBegin;
491cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
492cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, skip, 2);
493cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetSkipInfo_C", (PetscViewer, PetscBool), (viewer, skip));
4943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
495807ea322SDave May }
496807ea322SDave May 
497d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetSkipInfo_Binary(PetscViewer viewer, PetscBool skip)
498d71ae5a4SJacob Faibussowitsch {
499807ea322SDave May   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
500807ea322SDave May 
501807ea322SDave May   PetscFunctionBegin;
502cc843e7aSLisandro Dalcin   vbinary->skipinfo = skip;
5033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
504807ea322SDave May }
505807ea322SDave May 
506807ea322SDave May /*@
507*c410d8ccSBarry Smith     PetscViewerBinaryGetSkipInfo - check if viewer wrote a `.info` file
508807ea322SDave May 
509807ea322SDave May     Not Collective
510807ea322SDave May 
511807ea322SDave May     Input Parameter:
512811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
513807ea322SDave May 
514807ea322SDave May     Output Parameter:
515*c410d8ccSBarry Smith .   skip - `PETSC_TRUE` implies the `.info` file was not generated
516807ea322SDave May 
517807ea322SDave May     Level: advanced
518807ea322SDave May 
519811af0c4SBarry Smith     Note:
520811af0c4SBarry Smith     This must be called after `PetscViewerSetType()`
521807ea322SDave May 
522d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`,
523*c410d8ccSBarry Smith           `PetscViewerBinarySetSkipOptions()`, `PetscViewerBinarySetSkipInfo()`, `PetscViewerBinaryGetInfoPointer()`
524807ea322SDave May @*/
525d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetSkipInfo(PetscViewer viewer, PetscBool *skip)
526d71ae5a4SJacob Faibussowitsch {
527807ea322SDave May   PetscFunctionBegin;
528cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
529cc843e7aSLisandro Dalcin   PetscValidBoolPointer(skip, 2);
530cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetSkipInfo_C", (PetscViewer, PetscBool *), (viewer, skip));
5313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
532807ea322SDave May }
533807ea322SDave May 
534d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetSkipInfo_Binary(PetscViewer viewer, PetscBool *skip)
535d71ae5a4SJacob Faibussowitsch {
536807ea322SDave May   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
537807ea322SDave May 
538807ea322SDave May   PetscFunctionBegin;
539cc843e7aSLisandro Dalcin   *skip = vbinary->skipinfo;
5403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
541807ea322SDave May }
542807ea322SDave May 
5435c6c1daeSBarry Smith /*@
544*c410d8ccSBarry Smith     PetscViewerBinarySetSkipOptions - do not use values in the PETSc options database when loading objects
5455c6c1daeSBarry Smith 
5465c6c1daeSBarry Smith     Not Collective
5475c6c1daeSBarry Smith 
5485c6c1daeSBarry Smith     Input Parameters:
549811af0c4SBarry Smith +   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
550811af0c4SBarry Smith -   skip - `PETSC_TRUE` means do not use the options from the options database
5515c6c1daeSBarry Smith 
5525c6c1daeSBarry Smith     Options Database Key:
553811af0c4SBarry Smith .   -viewer_binary_skip_options <true or false> - true means do not use the options from the options database
5545c6c1daeSBarry Smith 
5555c6c1daeSBarry Smith     Level: advanced
5565c6c1daeSBarry Smith 
557811af0c4SBarry Smith     Note:
558811af0c4SBarry Smith     This must be called after `PetscViewerSetType()`
5595c6c1daeSBarry Smith 
560d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
561db781477SPatrick Sanan           `PetscViewerBinaryGetSkipOptions()`
5625c6c1daeSBarry Smith @*/
563d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetSkipOptions(PetscViewer viewer, PetscBool skip)
564d71ae5a4SJacob Faibussowitsch {
5655c6c1daeSBarry Smith   PetscFunctionBegin;
566cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
567cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, skip, 2);
568cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetSkipOptions_C", (PetscViewer, PetscBool), (viewer, skip));
5693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
570807ea322SDave May }
571807ea322SDave May 
572d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetSkipOptions_Binary(PetscViewer viewer, PetscBool skip)
573d71ae5a4SJacob Faibussowitsch {
574cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
575807ea322SDave May 
576807ea322SDave May   PetscFunctionBegin;
577cc843e7aSLisandro Dalcin   vbinary->skipoptions = skip;
5783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5795c6c1daeSBarry Smith }
5805c6c1daeSBarry Smith 
5815c6c1daeSBarry Smith /*@
5825c6c1daeSBarry Smith     PetscViewerBinaryGetSkipOptions - checks if viewer uses the PETSc options database when loading objects
5835c6c1daeSBarry Smith 
5845c6c1daeSBarry Smith     Not Collective
5855c6c1daeSBarry Smith 
5865c6c1daeSBarry Smith     Input Parameter:
587811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
5885c6c1daeSBarry Smith 
5895c6c1daeSBarry Smith     Output Parameter:
590811af0c4SBarry Smith .   skip - `PETSC_TRUE` means do not use
5915c6c1daeSBarry Smith 
5925c6c1daeSBarry Smith     Level: advanced
5935c6c1daeSBarry Smith 
594811af0c4SBarry Smith     Note:
595811af0c4SBarry Smith     This must be called after `PetscViewerSetType()`
5965c6c1daeSBarry Smith 
597d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
598db781477SPatrick Sanan           `PetscViewerBinarySetSkipOptions()`
5995c6c1daeSBarry Smith @*/
600d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetSkipOptions(PetscViewer viewer, PetscBool *skip)
601d71ae5a4SJacob Faibussowitsch {
6025c6c1daeSBarry Smith   PetscFunctionBegin;
603cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
604cc843e7aSLisandro Dalcin   PetscValidBoolPointer(skip, 2);
605cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetSkipOptions_C", (PetscViewer, PetscBool *), (viewer, skip));
6063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6075c6c1daeSBarry Smith }
6085c6c1daeSBarry Smith 
609d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetSkipOptions_Binary(PetscViewer viewer, PetscBool *skip)
610d71ae5a4SJacob Faibussowitsch {
611d21b9a37SPierre Jolivet   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
6125c6c1daeSBarry Smith 
6135c6c1daeSBarry Smith   PetscFunctionBegin;
614cc843e7aSLisandro Dalcin   *skip = vbinary->skipoptions;
6153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6165c6c1daeSBarry Smith }
6175c6c1daeSBarry Smith 
6185c6c1daeSBarry Smith /*@
6195c6c1daeSBarry Smith     PetscViewerBinarySetSkipHeader - do not write a header with size information on output, just raw data
6205c6c1daeSBarry Smith 
6215c6c1daeSBarry Smith     Not Collective
6225c6c1daeSBarry Smith 
6235c6c1daeSBarry Smith     Input Parameters:
624811af0c4SBarry Smith +   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
625811af0c4SBarry Smith -   skip - `PETSC_TRUE` means do not write header
6265c6c1daeSBarry Smith 
6275c6c1daeSBarry Smith     Options Database Key:
628811af0c4SBarry Smith .   -viewer_binary_skip_header <true or false> - true means do not write header
6295c6c1daeSBarry Smith 
6305c6c1daeSBarry Smith     Level: advanced
6315c6c1daeSBarry Smith 
632*c410d8ccSBarry Smith     Notes:
633811af0c4SBarry Smith       This must be called after `PetscViewerSetType()`
6345c6c1daeSBarry Smith 
635*c410d8ccSBarry Smith     If this option is selected, the output file cannot be read with the `XXXLoad()` such as `VecLoad()`
636*c410d8ccSBarry Smith 
637d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
638db781477SPatrick Sanan           `PetscViewerBinaryGetSkipHeader()`
6395c6c1daeSBarry Smith @*/
640d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetSkipHeader(PetscViewer viewer, PetscBool skip)
641d71ae5a4SJacob Faibussowitsch {
6425c6c1daeSBarry Smith   PetscFunctionBegin;
643cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
644cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, skip, 2);
645cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetSkipHeader_C", (PetscViewer, PetscBool), (viewer, skip));
6463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6475c6c1daeSBarry Smith }
6485c6c1daeSBarry Smith 
649d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetSkipHeader_Binary(PetscViewer viewer, PetscBool skip)
650d71ae5a4SJacob Faibussowitsch {
6515c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
6525c6c1daeSBarry Smith 
6535c6c1daeSBarry Smith   PetscFunctionBegin;
654cc843e7aSLisandro Dalcin   vbinary->skipheader = skip;
6553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6565c6c1daeSBarry Smith }
6575c6c1daeSBarry Smith 
6585c6c1daeSBarry Smith /*@
6595c6c1daeSBarry Smith     PetscViewerBinaryGetSkipHeader - checks whether to write a header with size information on output, or just raw data
6605c6c1daeSBarry Smith 
6615c6c1daeSBarry Smith     Not Collective
6625c6c1daeSBarry Smith 
6635c6c1daeSBarry Smith     Input Parameter:
664811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
6655c6c1daeSBarry Smith 
6665c6c1daeSBarry Smith     Output Parameter:
667811af0c4SBarry Smith .   skip - `PETSC_TRUE` means do not write header
6685c6c1daeSBarry Smith 
6695c6c1daeSBarry Smith     Level: advanced
6705c6c1daeSBarry Smith 
67195452b02SPatrick Sanan     Notes:
67295452b02SPatrick Sanan     This must be called after PetscViewerSetType()
6735c6c1daeSBarry Smith 
674811af0c4SBarry Smith     Returns `PETSC_FALSE` for `PETSCSOCKETVIEWER`, you cannot skip the header for it.
6755c6c1daeSBarry Smith 
676d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
677db781477SPatrick Sanan           `PetscViewerBinarySetSkipHeader()`
6785c6c1daeSBarry Smith @*/
679d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetSkipHeader(PetscViewer viewer, PetscBool *skip)
680d71ae5a4SJacob Faibussowitsch {
6815c6c1daeSBarry Smith   PetscFunctionBegin;
682cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
683cc843e7aSLisandro Dalcin   PetscValidBoolPointer(skip, 2);
684cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetSkipHeader_C", (PetscViewer, PetscBool *), (viewer, skip));
6853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6865c6c1daeSBarry Smith }
6875c6c1daeSBarry Smith 
688d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetSkipHeader_Binary(PetscViewer viewer, PetscBool *skip)
689d71ae5a4SJacob Faibussowitsch {
690f3b211e4SSatish Balay   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
6915c6c1daeSBarry Smith 
6925c6c1daeSBarry Smith   PetscFunctionBegin;
693cc843e7aSLisandro Dalcin   *skip = vbinary->skipheader;
6943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6955c6c1daeSBarry Smith }
6965c6c1daeSBarry Smith 
6975c6c1daeSBarry Smith /*@C
6985c6c1daeSBarry Smith     PetscViewerBinaryGetInfoPointer - Extracts the file pointer for the ASCII
699*c410d8ccSBarry Smith           `.info` file associated with a binary file.
7005c6c1daeSBarry Smith 
701cf53795eSBarry Smith     Not Collective; No Fortran Support
7025c6c1daeSBarry Smith 
7035c6c1daeSBarry Smith     Input Parameter:
704811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
7055c6c1daeSBarry Smith 
7065c6c1daeSBarry Smith     Output Parameter:
707*c410d8ccSBarry Smith .   file - file pointer  Always returns `NULL` if not a binary viewer
7085c6c1daeSBarry Smith 
7095c6c1daeSBarry Smith     Level: advanced
7105c6c1daeSBarry Smith 
711811af0c4SBarry Smith     Note:
712811af0c4SBarry Smith       For writable binary `PetscViewer`s, the file pointer will only be valid for the
713*c410d8ccSBarry Smith     first processor in the MPI communicator that shares the `PetscViewer`.
7145c6c1daeSBarry Smith 
715*c410d8ccSBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinaryGetSkipInfo()`,
716*c410d8ccSBarry Smith           `PetscViewerBinarySetSkipInfo()`
7175c6c1daeSBarry Smith @*/
718d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetInfoPointer(PetscViewer viewer, FILE **file)
719d71ae5a4SJacob Faibussowitsch {
7205c6c1daeSBarry Smith   PetscFunctionBegin;
721cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
722cc843e7aSLisandro Dalcin   PetscValidPointer(file, 2);
7230298fd71SBarry Smith   *file = NULL;
724cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinaryGetInfoPointer_C", (PetscViewer, FILE **), (viewer, file));
7253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7265c6c1daeSBarry Smith }
7275c6c1daeSBarry Smith 
728d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetInfoPointer_Binary(PetscViewer viewer, FILE **file)
729d71ae5a4SJacob Faibussowitsch {
730cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
7315c6c1daeSBarry Smith 
7325c6c1daeSBarry Smith   PetscFunctionBegin;
7339566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
734cc843e7aSLisandro Dalcin   *file = vbinary->fdes_info;
735cc843e7aSLisandro Dalcin   if (viewer->format == PETSC_VIEWER_BINARY_MATLAB && !vbinary->matlabheaderwritten) {
7365c6c1daeSBarry Smith     if (vbinary->fdes_info) {
737cc843e7aSLisandro Dalcin       FILE *info = vbinary->fdes_info;
7389566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- begin code written by PetscViewerBinary for MATLAB format ---#\n"));
7399566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#$$ Set.filename = '%s';\n", vbinary->filename));
7409566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#$$ fd = PetscOpenFile(Set.filename);\n"));
7419566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n"));
742cc843e7aSLisandro Dalcin     }
743cc843e7aSLisandro Dalcin     vbinary->matlabheaderwritten = PETSC_TRUE;
7445c6c1daeSBarry Smith   }
7453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7465c6c1daeSBarry Smith }
7475c6c1daeSBarry Smith 
748e0385b85SDave May #if defined(PETSC_HAVE_MPIIO)
749d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_BinaryMPIIO(PetscViewer v)
750d71ae5a4SJacob Faibussowitsch {
751e0385b85SDave May   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
752e0385b85SDave May 
753e0385b85SDave May   PetscFunctionBegin;
75448a46eb9SPierre Jolivet   if (vbinary->mfdes != MPI_FILE_NULL) PetscCallMPI(MPI_File_close(&vbinary->mfdes));
75548a46eb9SPierre Jolivet   if (vbinary->mfsub != MPI_FILE_NULL) PetscCallMPI(MPI_File_close(&vbinary->mfsub));
756cc843e7aSLisandro Dalcin   vbinary->moff = 0;
7573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
758e0385b85SDave May }
759e0385b85SDave May #endif
760e0385b85SDave May 
761d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_BinarySTDIO(PetscViewer v)
762d71ae5a4SJacob Faibussowitsch {
763cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
764cc843e7aSLisandro Dalcin 
765cc843e7aSLisandro Dalcin   PetscFunctionBegin;
766cc843e7aSLisandro Dalcin   if (vbinary->fdes != -1) {
7679566063dSJacob Faibussowitsch     PetscCall(PetscBinaryClose(vbinary->fdes));
768cc843e7aSLisandro Dalcin     vbinary->fdes = -1;
769cc843e7aSLisandro Dalcin     if (vbinary->storecompressed) {
770cc843e7aSLisandro Dalcin       char        cmd[8 + PETSC_MAX_PATH_LEN], out[64 + PETSC_MAX_PATH_LEN] = "";
771cc843e7aSLisandro Dalcin       const char *gzfilename = vbinary->ogzfilename ? vbinary->ogzfilename : vbinary->filename;
772cc843e7aSLisandro Dalcin       /* compress the file */
7739566063dSJacob Faibussowitsch       PetscCall(PetscStrncpy(cmd, "gzip -f ", sizeof(cmd)));
7749566063dSJacob Faibussowitsch       PetscCall(PetscStrlcat(cmd, gzfilename, sizeof(cmd)));
775cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_POPEN)
776cc843e7aSLisandro Dalcin       {
777cc843e7aSLisandro Dalcin         FILE *fp;
7789566063dSJacob Faibussowitsch         PetscCall(PetscPOpen(PETSC_COMM_SELF, NULL, cmd, "r", &fp));
779cc73adaaSBarry Smith         PetscCheck(!fgets(out, (int)(sizeof(out) - 1), fp), PETSC_COMM_SELF, PETSC_ERR_LIB, "Error from command %s\n%s", cmd, out);
7809566063dSJacob Faibussowitsch         PetscCall(PetscPClose(PETSC_COMM_SELF, fp));
781cc843e7aSLisandro Dalcin       }
782cc843e7aSLisandro Dalcin #endif
783cc843e7aSLisandro Dalcin     }
784cc843e7aSLisandro Dalcin   }
7859566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary->ogzfilename));
7863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
787cc843e7aSLisandro Dalcin }
788cc843e7aSLisandro Dalcin 
789d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_BinaryInfo(PetscViewer v)
790d71ae5a4SJacob Faibussowitsch {
791cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
792cc843e7aSLisandro Dalcin 
793cc843e7aSLisandro Dalcin   PetscFunctionBegin;
794cc843e7aSLisandro Dalcin   if (v->format == PETSC_VIEWER_BINARY_MATLAB && vbinary->matlabheaderwritten) {
795cc843e7aSLisandro Dalcin     if (vbinary->fdes_info) {
796cc843e7aSLisandro Dalcin       FILE *info = vbinary->fdes_info;
7979566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- begin code written by PetscViewerBinary for MATLAB format ---#\n"));
7989566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#$$ close(fd);\n"));
7999566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n"));
800cc843e7aSLisandro Dalcin     }
801cc843e7aSLisandro Dalcin   }
802cc843e7aSLisandro Dalcin   if (vbinary->fdes_info) {
803cc843e7aSLisandro Dalcin     FILE *info         = vbinary->fdes_info;
804cc843e7aSLisandro Dalcin     vbinary->fdes_info = NULL;
805cc73adaaSBarry Smith     PetscCheck(!fclose(info), PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file");
806cc843e7aSLisandro Dalcin   }
8073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
808cc843e7aSLisandro Dalcin }
809cc843e7aSLisandro Dalcin 
810d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_Binary(PetscViewer v)
811d71ae5a4SJacob Faibussowitsch {
812cc843e7aSLisandro Dalcin   PetscFunctionBegin;
813cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
8149566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_BinaryMPIIO(v));
815cc843e7aSLisandro Dalcin #endif
8169566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_BinarySTDIO(v));
8179566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_BinaryInfo(v));
8183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
819cc843e7aSLisandro Dalcin }
820cc843e7aSLisandro Dalcin 
821d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerDestroy_Binary(PetscViewer v)
822d71ae5a4SJacob Faibussowitsch {
8235c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
8245c6c1daeSBarry Smith 
8255c6c1daeSBarry Smith   PetscFunctionBegin;
8269566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_Binary(v));
8279566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary->filename));
8289566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary));
8292e956fe4SStefano Zampini   PetscCall(PetscViewerBinaryClearFunctionList(v));
8303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
831e0385b85SDave May }
8325c6c1daeSBarry Smith 
8335c6c1daeSBarry Smith /*@C
8345c6c1daeSBarry Smith    PetscViewerBinaryOpen - Opens a file for binary input/output.
8355c6c1daeSBarry Smith 
836d083f849SBarry Smith    Collective
8375c6c1daeSBarry Smith 
8385c6c1daeSBarry Smith    Input Parameters:
8395c6c1daeSBarry Smith +  comm - MPI communicator
8405c6c1daeSBarry Smith .  name - name of file
841cc843e7aSLisandro Dalcin -  mode - open mode of file
8423f423023SBarry Smith .vb
8433f423023SBarry Smith     FILE_MODE_WRITE - create new file for binary output
8443f423023SBarry Smith     FILE_MODE_READ - open existing file for binary input
8453f423023SBarry Smith     FILE_MODE_APPEND - open existing file for binary output
8463f423023SBarry Smith .ve
8475c6c1daeSBarry Smith 
8485c6c1daeSBarry Smith    Output Parameter:
849cc843e7aSLisandro Dalcin .  viewer - PetscViewer for binary input/output to use with the specified file
8505c6c1daeSBarry Smith 
8515c6c1daeSBarry Smith     Options Database Keys:
852811af0c4SBarry Smith +    -viewer_binary_filename <name> - name of file to use
853811af0c4SBarry Smith .    -viewer_binary_skip_info - true to skip opening an info file
854811af0c4SBarry Smith .    -viewer_binary_skip_options - true to not use options database while creating viewer
855811af0c4SBarry Smith .    -viewer_binary_skip_header - true to skip output object headers to the file
856811af0c4SBarry Smith -    -viewer_binary_mpiio - true to use MPI-IO for input and output to the file (more scalable for large problems)
8575c6c1daeSBarry Smith 
8585c6c1daeSBarry Smith    Level: beginner
8595c6c1daeSBarry Smith 
8605c6c1daeSBarry Smith    Note:
861811af0c4SBarry Smith    This `PetscViewer` should be destroyed with `PetscViewerDestroy()`.
8625c6c1daeSBarry Smith 
8635c6c1daeSBarry Smith     For reading files, the filename may begin with ftp:// or http:// and/or
8645c6c1daeSBarry Smith     end with .gz; in this case file is brought over and uncompressed.
8655c6c1daeSBarry Smith 
8665c6c1daeSBarry Smith     For creating files, if the file name ends with .gz it is automatically
8675c6c1daeSBarry Smith     compressed when closed.
8685c6c1daeSBarry Smith 
869d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
870db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
871db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`, `PetscViewerBinarySetUseMPIIO()`,
872db781477SPatrick Sanan           `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
8735c6c1daeSBarry Smith @*/
874d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryOpen(MPI_Comm comm, const char name[], PetscFileMode mode, PetscViewer *viewer)
875d71ae5a4SJacob Faibussowitsch {
8765c6c1daeSBarry Smith   PetscFunctionBegin;
8779566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(comm, viewer));
8789566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(*viewer, PETSCVIEWERBINARY));
8799566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(*viewer, mode));
8809566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetName(*viewer, name));
8819566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions(*viewer));
8823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8835c6c1daeSBarry Smith }
8845c6c1daeSBarry Smith 
8855c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
886d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryWriteReadMPIIO(PetscViewer viewer, void *data, PetscInt num, PetscInt *count, PetscDataType dtype, PetscBool write)
887d71ae5a4SJacob Faibussowitsch {
88822a8f86cSLisandro Dalcin   MPI_Comm            comm    = PetscObjectComm((PetscObject)viewer);
8895c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
89022a8f86cSLisandro Dalcin   MPI_File            mfdes   = vbinary->mfdes;
8915c6c1daeSBarry Smith   MPI_Datatype        mdtype;
89269e10bbaSLisandro Dalcin   PetscMPIInt         rank, cnt;
8935c6c1daeSBarry Smith   MPI_Status          status;
8945c6c1daeSBarry Smith   MPI_Aint            ul, dsize;
8955c6c1daeSBarry Smith 
8965c6c1daeSBarry Smith   PetscFunctionBegin;
8979566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
8989566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(num, &cnt));
8999566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToMPIDataType(dtype, &mdtype));
9005c6c1daeSBarry Smith   if (write) {
90148a46eb9SPierre Jolivet     if (rank == 0) PetscCall(MPIU_File_write_at(mfdes, vbinary->moff, data, cnt, mdtype, &status));
9025c6c1daeSBarry Smith   } else {
903dd400576SPatrick Sanan     if (rank == 0) {
9049566063dSJacob Faibussowitsch       PetscCall(MPIU_File_read_at(mfdes, vbinary->moff, data, cnt, mdtype, &status));
9059566063dSJacob Faibussowitsch       if (cnt > 0) PetscCallMPI(MPI_Get_count(&status, mdtype, &cnt));
9065c6c1daeSBarry Smith     }
9079566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&cnt, 1, MPI_INT, 0, comm));
9089566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(data, cnt, mdtype, 0, comm));
90969e10bbaSLisandro Dalcin   }
9109566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_get_extent(mdtype, &ul, &dsize));
9115c6c1daeSBarry Smith   vbinary->moff += dsize * cnt;
9129860990eSLisandro Dalcin   if (count) *count = cnt;
9133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9145c6c1daeSBarry Smith }
9155c6c1daeSBarry Smith #endif
9165c6c1daeSBarry Smith 
9175c6c1daeSBarry Smith /*@C
9185c6c1daeSBarry Smith    PetscViewerBinaryRead - Reads from a binary file, all processors get the same result
9195c6c1daeSBarry Smith 
920d083f849SBarry Smith    Collective
9215c6c1daeSBarry Smith 
9225c6c1daeSBarry Smith    Input Parameters:
923*c410d8ccSBarry Smith +  viewer - the `PETSCVIEWERBINARY` viewer
924060da220SMatthew G. Knepley .  num - number of items of data to read
925907376e6SBarry Smith -  dtype - type of data to read
9265c6c1daeSBarry Smith 
927*c410d8ccSBarry Smith    Output Parameters:
928*c410d8ccSBarry Smith +  data - location of the read data, treated as an array of the type indicated by `dtype`
929*c410d8ccSBarry Smith -  count - number of items of data actually read, or `NULL`.
930f8e4bde8SMatthew G. Knepley 
9315c6c1daeSBarry Smith    Level: beginner
9325c6c1daeSBarry Smith 
933d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
934db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
935db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
9365c6c1daeSBarry Smith @*/
937d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryRead(PetscViewer viewer, void *data, PetscInt num, PetscInt *count, PetscDataType dtype)
938d71ae5a4SJacob Faibussowitsch {
93922a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
9405c6c1daeSBarry Smith 
94122a8f86cSLisandro Dalcin   PetscFunctionBegin;
94222a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
94322a8f86cSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, num, 3);
9449566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
94522a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
9465c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
947bc196f7cSDave May   if (vbinary->usempiio) {
9489566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWriteReadMPIIO(viewer, data, num, count, dtype, PETSC_FALSE));
9495c6c1daeSBarry Smith   } else {
9505c6c1daeSBarry Smith #endif
9519566063dSJacob Faibussowitsch     PetscCall(PetscBinarySynchronizedRead(PetscObjectComm((PetscObject)viewer), vbinary->fdes, data, num, count, dtype));
9525c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
9535c6c1daeSBarry Smith   }
9545c6c1daeSBarry Smith #endif
9553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9565c6c1daeSBarry Smith }
9575c6c1daeSBarry Smith 
9585c6c1daeSBarry Smith /*@C
959811af0c4SBarry Smith    PetscViewerBinaryWrite - writes to a binary file, only from the first MPI rank
9605c6c1daeSBarry Smith 
961d083f849SBarry Smith    Collective
9625c6c1daeSBarry Smith 
9635c6c1daeSBarry Smith    Input Parameters:
964*c410d8ccSBarry Smith +  viewer - the `PETSCVIEWERBINARY` viewer
965*c410d8ccSBarry Smith .  data - location of data, treated as an array of the type indicated by `dtype`
9665824c9d0SPatrick Sanan .  count - number of items of data to write
967f253e43cSLisandro Dalcin -  dtype - type of data to write
9685c6c1daeSBarry Smith 
9695c6c1daeSBarry Smith    Level: beginner
9705c6c1daeSBarry Smith 
971d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
972db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`, `PetscDataType`
973db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
9745c6c1daeSBarry Smith @*/
975d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryWrite(PetscViewer viewer, const void *data, PetscInt count, PetscDataType dtype)
976d71ae5a4SJacob Faibussowitsch {
97722a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
9785c6c1daeSBarry Smith 
9795c6c1daeSBarry Smith   PetscFunctionBegin;
98022a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
98122a8f86cSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, count, 3);
9829566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
98322a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
9845c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
985bc196f7cSDave May   if (vbinary->usempiio) {
9869566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWriteReadMPIIO(viewer, (void *)data, count, NULL, dtype, PETSC_TRUE));
9875c6c1daeSBarry Smith   } else {
9885c6c1daeSBarry Smith #endif
9899566063dSJacob Faibussowitsch     PetscCall(PetscBinarySynchronizedWrite(PetscObjectComm((PetscObject)viewer), vbinary->fdes, data, count, dtype));
9905c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
9915c6c1daeSBarry Smith   }
9925c6c1daeSBarry Smith #endif
9933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9945c6c1daeSBarry Smith }
9955c6c1daeSBarry Smith 
996d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryWriteReadAll(PetscViewer viewer, PetscBool write, void *data, PetscInt count, PetscInt start, PetscInt total, PetscDataType dtype)
997d71ae5a4SJacob Faibussowitsch {
9985972f5f3SLisandro Dalcin   MPI_Comm              comm = PetscObjectComm((PetscObject)viewer);
9995972f5f3SLisandro Dalcin   PetscMPIInt           size, rank;
10005972f5f3SLisandro Dalcin   MPI_Datatype          mdtype;
1001ec4bef21SJose E. Roman   PETSC_UNUSED MPI_Aint lb;
1002ec4bef21SJose E. Roman   MPI_Aint              dsize;
10035972f5f3SLisandro Dalcin   PetscBool             useMPIIO;
10045972f5f3SLisandro Dalcin 
10055972f5f3SLisandro Dalcin   PetscFunctionBegin;
10065972f5f3SLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
1007064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveBool(viewer, ((start >= 0) || (start == PETSC_DETERMINE)), 5);
1008064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveBool(viewer, ((total >= 0) || (total == PETSC_DETERMINE)), 6);
1009064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveInt(viewer, total, 6);
10109566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
10115972f5f3SLisandro Dalcin 
10129566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToMPIDataType(dtype, &mdtype));
10139566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_get_extent(mdtype, &lb, &dsize));
10149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
10159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
10165972f5f3SLisandro Dalcin 
10179566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetUseMPIIO(viewer, &useMPIIO));
10185972f5f3SLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
10195972f5f3SLisandro Dalcin   if (useMPIIO) {
10205972f5f3SLisandro Dalcin     MPI_File    mfdes;
10215972f5f3SLisandro Dalcin     MPI_Offset  off;
10225972f5f3SLisandro Dalcin     PetscMPIInt cnt;
10235972f5f3SLisandro Dalcin 
10245972f5f3SLisandro Dalcin     if (start == PETSC_DETERMINE) {
10259566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Scan(&count, &start, 1, MPIU_INT, MPI_SUM, comm));
10265972f5f3SLisandro Dalcin       start -= count;
10275972f5f3SLisandro Dalcin     }
10285972f5f3SLisandro Dalcin     if (total == PETSC_DETERMINE) {
10295972f5f3SLisandro Dalcin       total = start + count;
10309566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Bcast(&total, 1, MPIU_INT, size - 1, comm));
10315972f5f3SLisandro Dalcin     }
10329566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(count, &cnt));
10339566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryGetMPIIODescriptor(viewer, &mfdes));
10349566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryGetMPIIOOffset(viewer, &off));
10355972f5f3SLisandro Dalcin     off += (MPI_Offset)(start * dsize);
10365972f5f3SLisandro Dalcin     if (write) {
10379566063dSJacob Faibussowitsch       PetscCall(MPIU_File_write_at_all(mfdes, off, data, cnt, mdtype, MPI_STATUS_IGNORE));
10385972f5f3SLisandro Dalcin     } else {
10399566063dSJacob Faibussowitsch       PetscCall(MPIU_File_read_at_all(mfdes, off, data, cnt, mdtype, MPI_STATUS_IGNORE));
10405972f5f3SLisandro Dalcin     }
10415972f5f3SLisandro Dalcin     off = (MPI_Offset)(total * dsize);
10429566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryAddMPIIOOffset(viewer, off));
10433ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
10445972f5f3SLisandro Dalcin   }
10455972f5f3SLisandro Dalcin #endif
10465972f5f3SLisandro Dalcin   {
10475972f5f3SLisandro Dalcin     int         fdes;
10485972f5f3SLisandro Dalcin     char       *workbuf = NULL;
1049dd400576SPatrick Sanan     PetscInt    tcount = rank == 0 ? 0 : count, maxcount = 0, message_count, flowcontrolcount;
10505972f5f3SLisandro Dalcin     PetscMPIInt tag, cnt, maxcnt, scnt = 0, rcnt = 0, j;
10515972f5f3SLisandro Dalcin     MPI_Status  status;
10525972f5f3SLisandro Dalcin 
10539566063dSJacob Faibussowitsch     PetscCall(PetscCommGetNewTag(comm, &tag));
10549566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Reduce(&tcount, &maxcount, 1, MPIU_INT, MPI_MAX, 0, comm));
10559566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(maxcount, &maxcnt));
10569566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(count, &cnt));
10575972f5f3SLisandro Dalcin 
10589566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryGetDescriptor(viewer, &fdes));
10599566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlowControlStart(viewer, &message_count, &flowcontrolcount));
1060dd400576SPatrick Sanan     if (rank == 0) {
10619566063dSJacob Faibussowitsch       PetscCall(PetscMalloc(maxcnt * dsize, &workbuf));
10625972f5f3SLisandro Dalcin       if (write) {
10639566063dSJacob Faibussowitsch         PetscCall(PetscBinaryWrite(fdes, data, cnt, dtype));
10645972f5f3SLisandro Dalcin       } else {
10659566063dSJacob Faibussowitsch         PetscCall(PetscBinaryRead(fdes, data, cnt, NULL, dtype));
10665972f5f3SLisandro Dalcin       }
10675972f5f3SLisandro Dalcin       for (j = 1; j < size; j++) {
10689566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlowControlStepMain(viewer, j, &message_count, flowcontrolcount));
10695972f5f3SLisandro Dalcin         if (write) {
10709566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Recv(workbuf, maxcnt, mdtype, j, tag, comm, &status));
10719566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Get_count(&status, mdtype, &rcnt));
10729566063dSJacob Faibussowitsch           PetscCall(PetscBinaryWrite(fdes, workbuf, rcnt, dtype));
10735972f5f3SLisandro Dalcin         } else {
10749566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Recv(&scnt, 1, MPI_INT, j, tag, comm, MPI_STATUS_IGNORE));
10759566063dSJacob Faibussowitsch           PetscCall(PetscBinaryRead(fdes, workbuf, scnt, NULL, dtype));
10769566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Send(workbuf, scnt, mdtype, j, tag, comm));
10775972f5f3SLisandro Dalcin         }
10785972f5f3SLisandro Dalcin       }
10799566063dSJacob Faibussowitsch       PetscCall(PetscFree(workbuf));
10809566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlowControlEndMain(viewer, &message_count));
10815972f5f3SLisandro Dalcin     } else {
10829566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlowControlStepWorker(viewer, rank, &message_count));
10835972f5f3SLisandro Dalcin       if (write) {
10849566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(data, cnt, mdtype, 0, tag, comm));
10855972f5f3SLisandro Dalcin       } else {
10869566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(&cnt, 1, MPI_INT, 0, tag, comm));
10879566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Recv(data, cnt, mdtype, 0, tag, comm, MPI_STATUS_IGNORE));
10885972f5f3SLisandro Dalcin       }
10899566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlowControlEndWorker(viewer, &message_count));
10905972f5f3SLisandro Dalcin     }
10915972f5f3SLisandro Dalcin   }
10923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10935972f5f3SLisandro Dalcin }
10945972f5f3SLisandro Dalcin 
10955972f5f3SLisandro Dalcin /*@C
1096*c410d8ccSBarry Smith    PetscViewerBinaryReadAll - reads from a binary file from all MPI processes, each rank receives its own portion of the data
10975972f5f3SLisandro Dalcin 
10985972f5f3SLisandro Dalcin    Collective
10995972f5f3SLisandro Dalcin 
11005972f5f3SLisandro Dalcin    Input Parameters:
1101*c410d8ccSBarry Smith +  viewer - the `PETSCVIEWERBINARY` viewer
11025972f5f3SLisandro Dalcin .  count - local number of items of data to read
1103811af0c4SBarry Smith .  start - local start, can be `PETSC_DETERMINE`
1104811af0c4SBarry Smith .  total - global number of items of data to read, can be `PETSC_DETERMINE`
11055972f5f3SLisandro Dalcin -  dtype - type of data to read
11065972f5f3SLisandro Dalcin 
1107*c410d8ccSBarry Smith    Output Parameter:
1108*c410d8ccSBarry Smith .  data - location of data, treated as an array of type indicated by `dtype`
1109*c410d8ccSBarry Smith 
11105972f5f3SLisandro Dalcin    Level: advanced
11115972f5f3SLisandro Dalcin 
1112d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryRead()`, `PetscViewerBinaryWriteAll()`
11135972f5f3SLisandro Dalcin @*/
1114d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryReadAll(PetscViewer viewer, void *data, PetscInt count, PetscInt start, PetscInt total, PetscDataType dtype)
1115d71ae5a4SJacob Faibussowitsch {
11165972f5f3SLisandro Dalcin   PetscFunctionBegin;
11179566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWriteReadAll(viewer, PETSC_FALSE, data, count, start, total, dtype));
11183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11195972f5f3SLisandro Dalcin }
11205972f5f3SLisandro Dalcin 
11215972f5f3SLisandro Dalcin /*@C
1122*c410d8ccSBarry Smith    PetscViewerBinaryWriteAll - writes to a binary file from all MPI processes, each rank writes its own portion of the data
11235972f5f3SLisandro Dalcin 
11245972f5f3SLisandro Dalcin    Collective
11255972f5f3SLisandro Dalcin 
11265972f5f3SLisandro Dalcin    Input Parameters:
1127*c410d8ccSBarry Smith +  viewer - the `PETSCVIEWERBINARY` viewer
11285972f5f3SLisandro Dalcin .  data - location of data
1129*c410d8ccSBarry Smith .  count - local number of items of data to write, treated as an array of type indicated by `dtype`
1130811af0c4SBarry Smith .  start - local start, can be `PETSC_DETERMINE`
1131811af0c4SBarry Smith .  total - global number of items of data to write, can be `PETSC_DETERMINE`
11325972f5f3SLisandro Dalcin -  dtype - type of data to write
11335972f5f3SLisandro Dalcin 
11345972f5f3SLisandro Dalcin    Level: advanced
11355972f5f3SLisandro Dalcin 
1136d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryWriteAll()`, `PetscViewerBinaryReadAll()`
11375972f5f3SLisandro Dalcin @*/
1138d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryWriteAll(PetscViewer viewer, const void *data, PetscInt count, PetscInt start, PetscInt total, PetscDataType dtype)
1139d71ae5a4SJacob Faibussowitsch {
11405972f5f3SLisandro Dalcin   PetscFunctionBegin;
11419566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWriteReadAll(viewer, PETSC_TRUE, (void *)data, count, start, total, dtype));
11423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11435972f5f3SLisandro Dalcin }
11445972f5f3SLisandro Dalcin 
11455c6c1daeSBarry Smith /*@C
1146811af0c4SBarry Smith    PetscViewerBinaryWriteStringArray - writes to a binary file, only from the first MPI rank, an array of strings
11475c6c1daeSBarry Smith 
1148d083f849SBarry Smith    Collective
11495c6c1daeSBarry Smith 
11505c6c1daeSBarry Smith    Input Parameters:
1151*c410d8ccSBarry Smith +  viewer - the `PETSCVIEWERBINARY` viewer
11525c6c1daeSBarry Smith -  data - location of the array of strings
11535c6c1daeSBarry Smith 
11545c6c1daeSBarry Smith    Level: intermediate
11555c6c1daeSBarry Smith 
1156811af0c4SBarry Smith     Note:
11573f423023SBarry Smith     The array of strings must be `NULL` terminated
11585c6c1daeSBarry Smith 
1159d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
1160db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
1161db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
11625c6c1daeSBarry Smith @*/
1163d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryWriteStringArray(PetscViewer viewer, const char *const *data)
1164d71ae5a4SJacob Faibussowitsch {
11655c6c1daeSBarry Smith   PetscInt i, n = 0, *sizes;
1166cc843e7aSLisandro Dalcin   size_t   len;
11675c6c1daeSBarry Smith 
116822a8f86cSLisandro Dalcin   PetscFunctionBegin;
11699566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
11705c6c1daeSBarry Smith   /* count number of strings */
11719371c9d4SSatish Balay   while (data[n++])
11729371c9d4SSatish Balay     ;
11735c6c1daeSBarry Smith   n--;
11749566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n + 1, &sizes));
11755c6c1daeSBarry Smith   sizes[0] = n;
11765c6c1daeSBarry Smith   for (i = 0; i < n; i++) {
11779566063dSJacob Faibussowitsch     PetscCall(PetscStrlen(data[i], &len));
1178cc843e7aSLisandro Dalcin     sizes[i + 1] = (PetscInt)len + 1; /* size includes space for the null terminator */
11795c6c1daeSBarry Smith   }
11809566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWrite(viewer, sizes, n + 1, PETSC_INT));
118148a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall(PetscViewerBinaryWrite(viewer, (void *)data[i], sizes[i + 1], PETSC_CHAR));
11829566063dSJacob Faibussowitsch   PetscCall(PetscFree(sizes));
11833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11845c6c1daeSBarry Smith }
11855c6c1daeSBarry Smith 
11865c6c1daeSBarry Smith /*@C
1187*c410d8ccSBarry Smith    PetscViewerBinaryReadStringArray - reads a binary file an array of strings to all MPI processes
11885c6c1daeSBarry Smith 
1189d083f849SBarry Smith    Collective
11905c6c1daeSBarry Smith 
11915c6c1daeSBarry Smith    Input Parameter:
1192*c410d8ccSBarry Smith .  viewer - the `PETSCVIEWERBINARY` viewer
11935c6c1daeSBarry Smith 
11945c6c1daeSBarry Smith    Output Parameter:
11955c6c1daeSBarry Smith .  data - location of the array of strings
11965c6c1daeSBarry Smith 
11975c6c1daeSBarry Smith    Level: intermediate
11985c6c1daeSBarry Smith 
1199811af0c4SBarry Smith     Note:
12003f423023SBarry Smith     The array of strings must `NULL` terminated
12015c6c1daeSBarry Smith 
1202d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
1203db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
1204db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
12055c6c1daeSBarry Smith @*/
1206d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryReadStringArray(PetscViewer viewer, char ***data)
1207d71ae5a4SJacob Faibussowitsch {
1208060da220SMatthew G. Knepley   PetscInt i, n, *sizes, N = 0;
12095c6c1daeSBarry Smith 
121022a8f86cSLisandro Dalcin   PetscFunctionBegin;
12119566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
12125c6c1daeSBarry Smith   /* count number of strings */
12139566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &n, 1, NULL, PETSC_INT));
12149566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &sizes));
12159566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, sizes, n, NULL, PETSC_INT));
1216a297a907SKarl Rupp   for (i = 0; i < n; i++) N += sizes[i];
12179566063dSJacob Faibussowitsch   PetscCall(PetscMalloc((n + 1) * sizeof(char *) + N * sizeof(char), data));
12185c6c1daeSBarry Smith   (*data)[0] = (char *)((*data) + n + 1);
1219a297a907SKarl Rupp   for (i = 1; i < n; i++) (*data)[i] = (*data)[i - 1] + sizes[i - 1];
12209566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, (*data)[0], N, NULL, PETSC_CHAR));
122102c9f0b5SLisandro Dalcin   (*data)[n] = NULL;
12229566063dSJacob Faibussowitsch   PetscCall(PetscFree(sizes));
12233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12245c6c1daeSBarry Smith }
12255c6c1daeSBarry Smith 
1226cc843e7aSLisandro Dalcin /*@C
1227cc843e7aSLisandro Dalcin      PetscViewerFileSetMode - Sets the open mode of file
1228cc843e7aSLisandro Dalcin 
1229c3339decSBarry Smith     Logically Collective
1230cc843e7aSLisandro Dalcin 
1231cc843e7aSLisandro Dalcin   Input Parameters:
1232811af0c4SBarry Smith +  viewer - the `PetscViewer`; must be a a `PETSCVIEWERBINARY`, `PETSCVIEWERMATLAB`, `PETSCVIEWERHDF5`, or `PETSCVIEWERASCII`  `PetscViewer`
1233cc843e7aSLisandro Dalcin -  mode - open mode of file
1234cf53795eSBarry Smith .vb
1235cf53795eSBarry Smith     FILE_MODE_WRITE - create new file for output
1236cf53795eSBarry Smith     FILE_MODE_READ - open existing file for input
1237cf53795eSBarry Smith     FILE_MODE_APPEND - open existing file for output
1238cf53795eSBarry Smith .ve
1239cc843e7aSLisandro Dalcin 
1240cc843e7aSLisandro Dalcin   Level: advanced
1241cc843e7aSLisandro Dalcin 
1242d1f92df0SBarry Smith .seealso: [](sec_viewers), `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
1243cc843e7aSLisandro Dalcin @*/
1244d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerFileSetMode(PetscViewer viewer, PetscFileMode mode)
1245d71ae5a4SJacob Faibussowitsch {
1246cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1247cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1248cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveEnum(viewer, mode, 2);
124908401ef6SPierre Jolivet   PetscCheck(mode != FILE_MODE_UNDEFINED, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Cannot set FILE_MODE_UNDEFINED");
1250f7d195e4SLawrence Mitchell   PetscCheck(mode >= FILE_MODE_UNDEFINED && mode <= FILE_MODE_APPEND_UPDATE, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_OUTOFRANGE, "Invalid file mode %d", (int)mode);
1251cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerFileSetMode_C", (PetscViewer, PetscFileMode), (viewer, mode));
12523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1253cc843e7aSLisandro Dalcin }
1254cc843e7aSLisandro Dalcin 
1255d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetMode_Binary(PetscViewer viewer, PetscFileMode mode)
1256d71ae5a4SJacob Faibussowitsch {
1257cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1258cc843e7aSLisandro Dalcin 
1259cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1260cc73adaaSBarry Smith   PetscCheck(!viewer->setupcalled || vbinary->filemode == mode, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Cannot change mode to %s after setup", PetscFileModes[mode]);
1261cc843e7aSLisandro Dalcin   vbinary->filemode = mode;
12623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1263cc843e7aSLisandro Dalcin }
1264cc843e7aSLisandro Dalcin 
1265cc843e7aSLisandro Dalcin /*@C
1266*c410d8ccSBarry Smith      PetscViewerFileGetMode - Gets the open mode of a file associated with a `PetscViewer`
1267cc843e7aSLisandro Dalcin 
1268cc843e7aSLisandro Dalcin     Not Collective
1269cc843e7aSLisandro Dalcin 
1270cc843e7aSLisandro Dalcin   Input Parameter:
1271811af0c4SBarry Smith .  viewer - the `PetscViewer`; must be a `PETSCVIEWERBINARY`, `PETSCVIEWERMATLAB`, `PETSCVIEWERHDF5`, or `PETSCVIEWERASCII`  `PetscViewer`
1272cc843e7aSLisandro Dalcin 
1273cc843e7aSLisandro Dalcin   Output Parameter:
1274cc843e7aSLisandro Dalcin .  mode - open mode of file
1275cf53795eSBarry Smith .vb
1276cf53795eSBarry Smith     FILE_MODE_WRITE - create new file for binary output
1277cf53795eSBarry Smith     FILE_MODE_READ - open existing file for binary input
1278cf53795eSBarry Smith     FILE_MODE_APPEND - open existing file for binary output
1279cf53795eSBarry Smith .ve
1280cc843e7aSLisandro Dalcin 
1281cc843e7aSLisandro Dalcin   Level: advanced
1282cc843e7aSLisandro Dalcin 
1283d1f92df0SBarry Smith .seealso: [](sec_viewers), `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
1284cc843e7aSLisandro Dalcin @*/
1285d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerFileGetMode(PetscViewer viewer, PetscFileMode *mode)
1286d71ae5a4SJacob Faibussowitsch {
1287cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1288cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1289cc843e7aSLisandro Dalcin   PetscValidPointer(mode, 2);
1290cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerFileGetMode_C", (PetscViewer, PetscFileMode *), (viewer, mode));
12913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1292cc843e7aSLisandro Dalcin }
1293cc843e7aSLisandro Dalcin 
1294d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetMode_Binary(PetscViewer viewer, PetscFileMode *mode)
1295d71ae5a4SJacob Faibussowitsch {
1296cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1297cc843e7aSLisandro Dalcin 
1298cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1299cc843e7aSLisandro Dalcin   *mode = vbinary->filemode;
13003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1301cc843e7aSLisandro Dalcin }
1302cc843e7aSLisandro Dalcin 
1303d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetName_Binary(PetscViewer viewer, const char name[])
1304d71ae5a4SJacob Faibussowitsch {
1305cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1306cc843e7aSLisandro Dalcin 
1307cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1308cc843e7aSLisandro Dalcin   if (viewer->setupcalled && vbinary->filename) {
1309cc843e7aSLisandro Dalcin     /* gzip can be run after the file with the previous filename has been closed */
13109566063dSJacob Faibussowitsch     PetscCall(PetscFree(vbinary->ogzfilename));
13119566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(vbinary->filename, &vbinary->ogzfilename));
1312cc843e7aSLisandro Dalcin   }
13139566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary->filename));
13149566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &vbinary->filename));
1315cc843e7aSLisandro Dalcin   viewer->setupcalled = PETSC_FALSE;
13163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1317cc843e7aSLisandro Dalcin }
1318cc843e7aSLisandro Dalcin 
1319d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetName_Binary(PetscViewer viewer, const char **name)
1320d71ae5a4SJacob Faibussowitsch {
13215c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
13225c6c1daeSBarry Smith 
13235c6c1daeSBarry Smith   PetscFunctionBegin;
13245c6c1daeSBarry Smith   *name = vbinary->filename;
13253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13265c6c1daeSBarry Smith }
13275c6c1daeSBarry Smith 
13285c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
1329d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetUp_BinaryMPIIO(PetscViewer viewer)
1330d71ae5a4SJacob Faibussowitsch {
13315c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1332e8a65b7dSLisandro Dalcin   int                 amode;
13335c6c1daeSBarry Smith 
13345c6c1daeSBarry Smith   PetscFunctionBegin;
1335a297a907SKarl Rupp   vbinary->storecompressed = PETSC_FALSE;
13365c6c1daeSBarry Smith 
1337cc843e7aSLisandro Dalcin   vbinary->moff = 0;
1338cc843e7aSLisandro Dalcin   switch (vbinary->filemode) {
1339d71ae5a4SJacob Faibussowitsch   case FILE_MODE_READ:
1340d71ae5a4SJacob Faibussowitsch     amode = MPI_MODE_RDONLY;
1341d71ae5a4SJacob Faibussowitsch     break;
1342d71ae5a4SJacob Faibussowitsch   case FILE_MODE_WRITE:
1343d71ae5a4SJacob Faibussowitsch     amode = MPI_MODE_WRONLY | MPI_MODE_CREATE;
1344d71ae5a4SJacob Faibussowitsch     break;
1345d71ae5a4SJacob Faibussowitsch   case FILE_MODE_APPEND:
1346d71ae5a4SJacob Faibussowitsch     amode = MPI_MODE_WRONLY | MPI_MODE_CREATE | MPI_MODE_APPEND;
1347d71ae5a4SJacob Faibussowitsch     break;
1348d71ae5a4SJacob Faibussowitsch   case FILE_MODE_UNDEFINED:
1349d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerSetUp()");
1350d71ae5a4SJacob Faibussowitsch   default:
1351d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[vbinary->filemode]);
13525c6c1daeSBarry Smith   }
13539566063dSJacob Faibussowitsch   PetscCallMPI(MPI_File_open(PetscObjectComm((PetscObject)viewer), vbinary->filename, amode, MPI_INFO_NULL, &vbinary->mfdes));
135422a8f86cSLisandro Dalcin   /*
135522a8f86cSLisandro Dalcin       The MPI standard does not have MPI_MODE_TRUNCATE. We emulate this behavior by setting the file size to zero.
135622a8f86cSLisandro Dalcin   */
13579566063dSJacob Faibussowitsch   if (vbinary->filemode == FILE_MODE_WRITE) PetscCallMPI(MPI_File_set_size(vbinary->mfdes, 0));
135822a8f86cSLisandro Dalcin   /*
135922a8f86cSLisandro Dalcin       Initially, all processes view the file as a linear byte stream. Therefore, for files opened with MPI_MODE_APPEND,
136022a8f86cSLisandro Dalcin       MPI_File_get_position[_shared](fh, &offset) returns the absolute byte position at the end of file.
136122a8f86cSLisandro Dalcin       Otherwise, we would need to call MPI_File_get_byte_offset(fh, offset, &byte_offset) to convert
136222a8f86cSLisandro Dalcin       the offset in etype units to an absolute byte position.
136322a8f86cSLisandro Dalcin    */
13649566063dSJacob Faibussowitsch   if (vbinary->filemode == FILE_MODE_APPEND) PetscCallMPI(MPI_File_get_position(vbinary->mfdes, &vbinary->moff));
13653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1366cc843e7aSLisandro Dalcin }
1367cc843e7aSLisandro Dalcin #endif
13685c6c1daeSBarry Smith 
1369d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetUp_BinarySTDIO(PetscViewer viewer)
1370d71ae5a4SJacob Faibussowitsch {
1371cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1372cc843e7aSLisandro Dalcin   const char         *fname;
1373bbcf679cSJacob Faibussowitsch   char                bname[PETSC_MAX_PATH_LEN], *gz = NULL;
1374cc843e7aSLisandro Dalcin   PetscBool           found;
1375cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
13765c6c1daeSBarry Smith 
1377cc843e7aSLisandro Dalcin   PetscFunctionBegin;
13789566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
13795c6c1daeSBarry Smith 
1380cc843e7aSLisandro Dalcin   /* if file name ends in .gz strip that off and note user wants file compressed */
1381cc843e7aSLisandro Dalcin   vbinary->storecompressed = PETSC_FALSE;
1382cc843e7aSLisandro Dalcin   if (vbinary->filemode == FILE_MODE_WRITE) {
13839566063dSJacob Faibussowitsch     PetscCall(PetscStrstr(vbinary->filename, ".gz", &gz));
13849371c9d4SSatish Balay     if (gz && gz[3] == 0) {
13859371c9d4SSatish Balay       *gz                      = 0;
13869371c9d4SSatish Balay       vbinary->storecompressed = PETSC_TRUE;
13879371c9d4SSatish Balay     }
1388cc843e7aSLisandro Dalcin   }
1389cc843e7aSLisandro Dalcin #if !defined(PETSC_HAVE_POPEN)
139028b400f6SJacob Faibussowitsch   PetscCheck(!vbinary->storecompressed, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP_SYS, "Cannot run gzip on this machine");
1391cc843e7aSLisandro Dalcin #endif
1392cc843e7aSLisandro Dalcin 
1393cc843e7aSLisandro Dalcin   fname = vbinary->filename;
1394cc843e7aSLisandro Dalcin   if (vbinary->filemode == FILE_MODE_READ) { /* possibly get the file from remote site or compressed file */
13959566063dSJacob Faibussowitsch     PetscCall(PetscFileRetrieve(PetscObjectComm((PetscObject)viewer), fname, bname, PETSC_MAX_PATH_LEN, &found));
139628b400f6SJacob Faibussowitsch     PetscCheck(found, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_OPEN, "Cannot locate file: %s", fname);
1397cc843e7aSLisandro Dalcin     fname = bname;
13985c6c1daeSBarry Smith   }
13995c6c1daeSBarry Smith 
1400cc843e7aSLisandro Dalcin   vbinary->fdes = -1;
1401dd400576SPatrick Sanan   if (rank == 0) { /* only first processor opens file*/
1402cc843e7aSLisandro Dalcin     PetscFileMode mode = vbinary->filemode;
1403cc843e7aSLisandro Dalcin     if (mode == FILE_MODE_APPEND) {
1404cc843e7aSLisandro Dalcin       /* check if asked to append to a non-existing file */
14059566063dSJacob Faibussowitsch       PetscCall(PetscTestFile(fname, '\0', &found));
1406cc843e7aSLisandro Dalcin       if (!found) mode = FILE_MODE_WRITE;
1407cc843e7aSLisandro Dalcin     }
14089566063dSJacob Faibussowitsch     PetscCall(PetscBinaryOpen(fname, mode, &vbinary->fdes));
1409cc843e7aSLisandro Dalcin   }
14103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1411cc843e7aSLisandro Dalcin }
1412cc843e7aSLisandro Dalcin 
1413d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetUp_BinaryInfo(PetscViewer viewer)
1414d71ae5a4SJacob Faibussowitsch {
1415cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1416cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
1417cc843e7aSLisandro Dalcin   PetscBool           found;
1418cc843e7aSLisandro Dalcin 
1419cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1420cc843e7aSLisandro Dalcin   vbinary->fdes_info = NULL;
14219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
1422dd400576SPatrick Sanan   if (!vbinary->skipinfo && (vbinary->filemode == FILE_MODE_READ || rank == 0)) {
1423cc843e7aSLisandro Dalcin     char infoname[PETSC_MAX_PATH_LEN], iname[PETSC_MAX_PATH_LEN], *gz;
1424cc843e7aSLisandro Dalcin 
14259566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(infoname, vbinary->filename, sizeof(infoname)));
1426cc843e7aSLisandro Dalcin     /* remove .gz if it ends file name */
14279566063dSJacob Faibussowitsch     PetscCall(PetscStrstr(infoname, ".gz", &gz));
1428cc843e7aSLisandro Dalcin     if (gz && gz[3] == 0) *gz = 0;
1429cc843e7aSLisandro Dalcin 
14309566063dSJacob Faibussowitsch     PetscCall(PetscStrlcat(infoname, ".info", sizeof(infoname)));
1431cc843e7aSLisandro Dalcin     if (vbinary->filemode == FILE_MODE_READ) {
14329566063dSJacob Faibussowitsch       PetscCall(PetscFixFilename(infoname, iname));
14339566063dSJacob Faibussowitsch       PetscCall(PetscFileRetrieve(PetscObjectComm((PetscObject)viewer), iname, infoname, PETSC_MAX_PATH_LEN, &found));
14349566063dSJacob Faibussowitsch       if (found) PetscCall(PetscOptionsInsertFile(PetscObjectComm((PetscObject)viewer), ((PetscObject)viewer)->options, infoname, PETSC_FALSE));
1435dd400576SPatrick Sanan     } else if (rank == 0) { /* write or append */
1436cc843e7aSLisandro Dalcin       const char *omode  = (vbinary->filemode == FILE_MODE_APPEND) ? "a" : "w";
1437cc843e7aSLisandro Dalcin       vbinary->fdes_info = fopen(infoname, omode);
143828b400f6SJacob Faibussowitsch       PetscCheck(vbinary->fdes_info, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open .info file %s for writing", infoname);
14395c6c1daeSBarry Smith     }
14405c6c1daeSBarry Smith   }
14413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14425c6c1daeSBarry Smith }
14435c6c1daeSBarry Smith 
1444d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerSetUp_Binary(PetscViewer viewer)
1445d71ae5a4SJacob Faibussowitsch {
14465c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1447cc843e7aSLisandro Dalcin   PetscBool           usempiio;
1448cc843e7aSLisandro Dalcin 
14495c6c1daeSBarry Smith   PetscFunctionBegin;
14509566063dSJacob Faibussowitsch   if (!vbinary->setfromoptionscalled) PetscCall(PetscViewerSetFromOptions(viewer));
145128b400f6SJacob Faibussowitsch   PetscCheck(vbinary->filename, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetName()");
145208401ef6SPierre Jolivet   PetscCheck(vbinary->filemode != (PetscFileMode)-1, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode()");
14539566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_Binary(viewer));
1454cc843e7aSLisandro Dalcin 
14559566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetUseMPIIO(viewer, &usempiio));
1456cc843e7aSLisandro Dalcin   if (usempiio) {
1457cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
14589566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetUp_BinaryMPIIO(viewer));
1459cc843e7aSLisandro Dalcin #endif
1460cc843e7aSLisandro Dalcin   } else {
14619566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetUp_BinarySTDIO(viewer));
1462cc843e7aSLisandro Dalcin   }
14639566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetUp_BinaryInfo(viewer));
1464cc843e7aSLisandro Dalcin 
14659566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectState((PetscObject)viewer, "File: %s", vbinary->filename));
14663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14675c6c1daeSBarry Smith }
14685c6c1daeSBarry Smith 
1469d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerView_Binary(PetscViewer v, PetscViewer viewer)
1470d71ae5a4SJacob Faibussowitsch {
1471cb6ad94fSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
1472cb6ad94fSLisandro Dalcin   const char         *fname   = vbinary->filename ? vbinary->filename : "not yet set";
1473cc843e7aSLisandro Dalcin   const char         *fmode   = vbinary->filemode != (PetscFileMode)-1 ? PetscFileModes[vbinary->filemode] : "not yet set";
1474cb6ad94fSLisandro Dalcin   PetscBool           usempiio;
14752bf49c77SBarry Smith 
14762bf49c77SBarry Smith   PetscFunctionBegin;
14779566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetUseMPIIO(v, &usempiio));
14789566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Filename: %s\n", fname));
14799566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Mode: %s (%s)\n", fmode, usempiio ? "mpiio" : "stdio"));
14803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14812bf49c77SBarry Smith }
14822bf49c77SBarry Smith 
1483d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerSetFromOptions_Binary(PetscViewer viewer, PetscOptionItems *PetscOptionsObject)
1484d71ae5a4SJacob Faibussowitsch {
148522a8f86cSLisandro Dalcin   PetscViewer_Binary *binary = (PetscViewer_Binary *)viewer->data;
1486d22fd6bcSDave May   char                defaultname[PETSC_MAX_PATH_LEN];
148703a1d59fSDave May   PetscBool           flg;
1488e0385b85SDave May 
148903a1d59fSDave May   PetscFunctionBegin;
14903ba16761SJacob Faibussowitsch   if (viewer->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
1491d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Binary PetscViewer Options");
14929566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(defaultname, PETSC_MAX_PATH_LEN - 1, "binaryoutput"));
14939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-viewer_binary_filename", "Specify filename", "PetscViewerFileSetName", defaultname, defaultname, sizeof(defaultname), &flg));
14949566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscViewerFileSetName_Binary(viewer, defaultname));
14959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_skip_info", "Skip writing/reading .info file", "PetscViewerBinarySetSkipInfo", binary->skipinfo, &binary->skipinfo, NULL));
14969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_skip_options", "Skip parsing Vec/Mat load options", "PetscViewerBinarySetSkipOptions", binary->skipoptions, &binary->skipoptions, NULL));
14979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_skip_header", "Skip writing/reading header information", "PetscViewerBinarySetSkipHeader", binary->skipheader, &binary->skipheader, NULL));
149803a1d59fSDave May #if defined(PETSC_HAVE_MPIIO)
14999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_mpiio", "Use MPI-IO functionality to write/read binary file", "PetscViewerBinarySetUseMPIIO", binary->usempiio, &binary->usempiio, NULL));
15005972f5f3SLisandro Dalcin #else
15019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_mpiio", "Use MPI-IO functionality to write/read binary file (NOT AVAILABLE)", "PetscViewerBinarySetUseMPIIO", PETSC_FALSE, NULL, NULL));
150203a1d59fSDave May #endif
1503d0609cedSBarry Smith   PetscOptionsHeadEnd();
1504bc196f7cSDave May   binary->setfromoptionscalled = PETSC_TRUE;
15053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
150603a1d59fSDave May }
150703a1d59fSDave May 
15088556b5ebSBarry Smith /*MC
15098556b5ebSBarry Smith    PETSCVIEWERBINARY - A viewer that saves to binary files
15108556b5ebSBarry Smith 
15111b266c99SBarry Smith   Level: beginner
15121b266c99SBarry Smith 
1513d1f92df0SBarry Smith .seealso: [](sec_viewers), `PetscViewerBinaryOpen()`, `PETSC_VIEWER_STDOUT_()`, `PETSC_VIEWER_STDOUT_SELF`, `PETSC_VIEWER_STDOUT_WORLD`, `PetscViewerCreate()`, `PetscViewerASCIIOpen()`,
1514811af0c4SBarry Smith           `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`, `PETSCVIEWERDRAW`, `PETSCVIEWERSOCKET`
1515811af0c4SBarry Smith           `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()`,
1516811af0c4SBarry Smith           `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`
15178556b5ebSBarry Smith M*/
15188556b5ebSBarry Smith 
1519d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscViewerCreate_Binary(PetscViewer v)
1520d71ae5a4SJacob Faibussowitsch {
15215c6c1daeSBarry Smith   PetscViewer_Binary *vbinary;
15225c6c1daeSBarry Smith 
15235c6c1daeSBarry Smith   PetscFunctionBegin;
15244dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&vbinary));
15255c6c1daeSBarry Smith   v->data = (void *)vbinary;
1526cc843e7aSLisandro Dalcin 
152703a1d59fSDave May   v->ops->setfromoptions   = PetscViewerSetFromOptions_Binary;
15285c6c1daeSBarry Smith   v->ops->destroy          = PetscViewerDestroy_Binary;
15292bf49c77SBarry Smith   v->ops->view             = PetscViewerView_Binary;
1530c98fd787SBarry Smith   v->ops->setup            = PetscViewerSetUp_Binary;
1531cc843e7aSLisandro Dalcin   v->ops->flush            = NULL; /* Should we support Flush() ? */
1532cc843e7aSLisandro Dalcin   v->ops->getsubviewer     = PetscViewerGetSubViewer_Binary;
1533cc843e7aSLisandro Dalcin   v->ops->restoresubviewer = PetscViewerRestoreSubViewer_Binary;
1534cc843e7aSLisandro Dalcin   v->ops->read             = PetscViewerBinaryRead;
1535cc843e7aSLisandro Dalcin 
1536cc843e7aSLisandro Dalcin   vbinary->fdes = -1;
1537e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
1538cc843e7aSLisandro Dalcin   vbinary->usempiio = PETSC_FALSE;
1539e8a65b7dSLisandro Dalcin   vbinary->mfdes    = MPI_FILE_NULL;
1540e8a65b7dSLisandro Dalcin   vbinary->mfsub    = MPI_FILE_NULL;
1541e8a65b7dSLisandro Dalcin #endif
1542cc843e7aSLisandro Dalcin   vbinary->filename        = NULL;
15437e4fd573SVaclav Hapla   vbinary->filemode        = FILE_MODE_UNDEFINED;
154402c9f0b5SLisandro Dalcin   vbinary->fdes_info       = NULL;
15455c6c1daeSBarry Smith   vbinary->skipinfo        = PETSC_FALSE;
15465c6c1daeSBarry Smith   vbinary->skipoptions     = PETSC_TRUE;
15475c6c1daeSBarry Smith   vbinary->skipheader      = PETSC_FALSE;
15485c6c1daeSBarry Smith   vbinary->storecompressed = PETSC_FALSE;
1549f90597f1SStefano Zampini   vbinary->ogzfilename     = NULL;
15505c6c1daeSBarry Smith   vbinary->flowcontrol     = 256; /* seems a good number for Cray XT-5 */
15515c6c1daeSBarry Smith 
1552cc843e7aSLisandro Dalcin   vbinary->setfromoptionscalled = PETSC_FALSE;
1553cc843e7aSLisandro Dalcin 
15549566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetFlowControl_C", PetscViewerBinaryGetFlowControl_Binary));
15559566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetFlowControl_C", PetscViewerBinarySetFlowControl_Binary));
15569566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipHeader_C", PetscViewerBinaryGetSkipHeader_Binary));
15579566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipHeader_C", PetscViewerBinarySetSkipHeader_Binary));
15589566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipOptions_C", PetscViewerBinaryGetSkipOptions_Binary));
15599566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipOptions_C", PetscViewerBinarySetSkipOptions_Binary));
15609566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipInfo_C", PetscViewerBinaryGetSkipInfo_Binary));
15619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipInfo_C", PetscViewerBinarySetSkipInfo_Binary));
15629566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetInfoPointer_C", PetscViewerBinaryGetInfoPointer_Binary));
15639566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_Binary));
15649566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_Binary));
15659566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_Binary));
15669566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_Binary));
15675c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
15689566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetUseMPIIO_C", PetscViewerBinaryGetUseMPIIO_Binary));
15699566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetUseMPIIO_C", PetscViewerBinarySetUseMPIIO_Binary));
15705c6c1daeSBarry Smith #endif
15713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15725c6c1daeSBarry Smith }
15735c6c1daeSBarry Smith 
15745c6c1daeSBarry Smith /*
15755c6c1daeSBarry Smith     The variable Petsc_Viewer_Binary_keyval is used to indicate an MPI attribute that
15765c6c1daeSBarry Smith   is attached to a communicator, in this case the attribute is a PetscViewer.
15775c6c1daeSBarry Smith */
1578d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_Binary_keyval = MPI_KEYVAL_INVALID;
15795c6c1daeSBarry Smith 
15805c6c1daeSBarry Smith /*@C
1581811af0c4SBarry Smith      PETSC_VIEWER_BINARY_ - Creates a `PETSCVIEWERBINARY` `PetscViewer` shared by all processors
15825c6c1daeSBarry Smith                      in a communicator.
15835c6c1daeSBarry Smith 
1584d083f849SBarry Smith      Collective
15855c6c1daeSBarry Smith 
15865c6c1daeSBarry Smith      Input Parameter:
1587811af0c4SBarry Smith .    comm - the MPI communicator to share the `PETSCVIEWERBINARY`
15885c6c1daeSBarry Smith 
15895c6c1daeSBarry Smith      Level: intermediate
15905c6c1daeSBarry Smith 
15915c6c1daeSBarry Smith    Options Database Keys:
159210699b91SBarry Smith +    -viewer_binary_filename <name> - filename in which to store the binary data, defaults to binaryoutput
159310699b91SBarry Smith .    -viewer_binary_skip_info - true means do not create .info file for this viewer
159410699b91SBarry Smith .    -viewer_binary_skip_options - true means do not use the options database for this viewer
159510699b91SBarry Smith .    -viewer_binary_skip_header - true means do not store the usual header information in the binary file
159610699b91SBarry Smith -    -viewer_binary_mpiio - true means use the file via MPI-IO, maybe faster for large files and many MPI ranks
15975c6c1daeSBarry Smith 
1598811af0c4SBarry Smith    Environmental variable:
159910699b91SBarry Smith -   PETSC_VIEWER_BINARY_FILENAME - filename in which to store the binary data, defaults to binaryoutput
16005c6c1daeSBarry Smith 
1601811af0c4SBarry Smith      Note:
1602811af0c4SBarry Smith      Unlike almost all other PETSc routines, `PETSC_VIEWER_BINARY_` does not return
16035c6c1daeSBarry Smith      an error code.  The binary PetscViewer is usually used in the form
16045c6c1daeSBarry Smith $       XXXView(XXX object,PETSC_VIEWER_BINARY_(comm));
16055c6c1daeSBarry Smith 
1606d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PETSC_VIEWER_BINARY_WORLD`, `PETSC_VIEWER_BINARY_SELF`, `PetscViewerBinaryOpen()`, `PetscViewerCreate()`,
1607db781477SPatrick Sanan           `PetscViewerDestroy()`
16085c6c1daeSBarry Smith @*/
1609d71ae5a4SJacob Faibussowitsch PetscViewer PETSC_VIEWER_BINARY_(MPI_Comm comm)
1610d71ae5a4SJacob Faibussowitsch {
16115c6c1daeSBarry Smith   PetscErrorCode ierr;
16123ba16761SJacob Faibussowitsch   PetscMPIInt    mpi_ierr;
16135c6c1daeSBarry Smith   PetscBool      flg;
16145c6c1daeSBarry Smith   PetscViewer    viewer;
16155c6c1daeSBarry Smith   char           fname[PETSC_MAX_PATH_LEN];
16165c6c1daeSBarry Smith   MPI_Comm       ncomm;
16175c6c1daeSBarry Smith 
16185c6c1daeSBarry Smith   PetscFunctionBegin;
16199371c9d4SSatish Balay   ierr = PetscCommDuplicate(comm, &ncomm, NULL);
16209371c9d4SSatish Balay   if (ierr) {
16213ba16761SJacob Faibussowitsch     ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16229371c9d4SSatish Balay     PetscFunctionReturn(NULL);
16239371c9d4SSatish Balay   }
16245c6c1daeSBarry Smith   if (Petsc_Viewer_Binary_keyval == MPI_KEYVAL_INVALID) {
16253ba16761SJacob Faibussowitsch     mpi_ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Viewer_Binary_keyval, NULL);
16263ba16761SJacob Faibussowitsch     if (mpi_ierr) {
16273ba16761SJacob Faibussowitsch       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16289371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16299371c9d4SSatish Balay     }
16305c6c1daeSBarry Smith   }
16313ba16761SJacob Faibussowitsch   mpi_ierr = MPI_Comm_get_attr(ncomm, Petsc_Viewer_Binary_keyval, (void **)&viewer, (int *)&flg);
16323ba16761SJacob Faibussowitsch   if (mpi_ierr) {
16333ba16761SJacob Faibussowitsch     ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16349371c9d4SSatish Balay     PetscFunctionReturn(NULL);
16359371c9d4SSatish Balay   }
16365c6c1daeSBarry Smith   if (!flg) { /* PetscViewer not yet created */
16375c6c1daeSBarry Smith     ierr = PetscOptionsGetenv(ncomm, "PETSC_VIEWER_BINARY_FILENAME", fname, PETSC_MAX_PATH_LEN, &flg);
16389371c9d4SSatish Balay     if (ierr) {
16393ba16761SJacob Faibussowitsch       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16409371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16419371c9d4SSatish Balay     }
16425c6c1daeSBarry Smith     if (!flg) {
1643c6a7a370SJeremy L Thompson       ierr = PetscStrncpy(fname, "binaryoutput", sizeof(fname));
16449371c9d4SSatish Balay       if (ierr) {
16453ba16761SJacob Faibussowitsch         ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16469371c9d4SSatish Balay         PetscFunctionReturn(NULL);
16479371c9d4SSatish Balay       }
16485c6c1daeSBarry Smith     }
16495c6c1daeSBarry Smith     ierr = PetscViewerBinaryOpen(ncomm, fname, FILE_MODE_WRITE, &viewer);
16509371c9d4SSatish Balay     if (ierr) {
16513ba16761SJacob Faibussowitsch       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16529371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16539371c9d4SSatish Balay     }
16545c6c1daeSBarry Smith     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
16559371c9d4SSatish Balay     if (ierr) {
16563ba16761SJacob Faibussowitsch       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16579371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16589371c9d4SSatish Balay     }
16593ba16761SJacob Faibussowitsch     mpi_ierr = MPI_Comm_set_attr(ncomm, Petsc_Viewer_Binary_keyval, (void *)viewer);
16603ba16761SJacob Faibussowitsch     if (mpi_ierr) {
16613ba16761SJacob Faibussowitsch       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16629371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16639371c9d4SSatish Balay     }
16645c6c1daeSBarry Smith   }
16655c6c1daeSBarry Smith   ierr = PetscCommDestroy(&ncomm);
16669371c9d4SSatish Balay   if (ierr) {
16673ba16761SJacob Faibussowitsch     ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16689371c9d4SSatish Balay     PetscFunctionReturn(NULL);
16699371c9d4SSatish Balay   }
16705c6c1daeSBarry Smith   PetscFunctionReturn(viewer);
16715c6c1daeSBarry Smith }
1672