xref: /petsc/src/sys/classes/viewer/impls/binary/binv.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1af0996ceSBarry Smith #include <petsc/private/viewerimpl.h> /*I   "petscviewer.h"   I*/
25c6c1daeSBarry Smith 
35c6c1daeSBarry Smith typedef struct {
45c6c1daeSBarry Smith   int fdes; /* file descriptor, ignored if using MPI IO */
55c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
6bc196f7cSDave May   PetscBool  usempiio;
75c6c1daeSBarry Smith   MPI_File   mfdes; /* ignored unless using MPI IO */
8e8a65b7dSLisandro Dalcin   MPI_File   mfsub; /* subviewer support */
95c6c1daeSBarry Smith   MPI_Offset moff;
105c6c1daeSBarry Smith #endif
11cc843e7aSLisandro Dalcin   char         *filename;            /* file name */
12cc843e7aSLisandro Dalcin   PetscFileMode filemode;            /* read/write/append mode */
135c6c1daeSBarry Smith   FILE         *fdes_info;           /* optional file containing info on binary file*/
145c6c1daeSBarry Smith   PetscBool     storecompressed;     /* gzip the write binary file when closing it*/
15f90597f1SStefano Zampini   char         *ogzfilename;         /* gzip can be run after the filename has been updated */
165c6c1daeSBarry Smith   PetscBool     skipinfo;            /* Don't create info file for writing; don't use for reading */
175c6c1daeSBarry Smith   PetscBool     skipoptions;         /* don't use PETSc options database when loading */
185c6c1daeSBarry Smith   PetscInt      flowcontrol;         /* allow only <flowcontrol> messages outstanding at a time while doing IO */
195c6c1daeSBarry Smith   PetscBool     skipheader;          /* don't write header, only raw data */
20a261c58fSBarry Smith   PetscBool     matlabheaderwritten; /* if format is PETSC_VIEWER_BINARY_MATLAB has the MATLAB .info header been written yet */
21c98fd787SBarry Smith   PetscBool     setfromoptionscalled;
225c6c1daeSBarry Smith } PetscViewer_Binary;
235c6c1daeSBarry Smith 
24*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryClearFunctionList(PetscViewer v)
25*d71ae5a4SJacob Faibussowitsch {
262e956fe4SStefano Zampini   PetscFunctionBegin;
272e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetFlowControl_C", NULL));
282e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetFlowControl_C", NULL));
292e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipHeader_C", NULL));
302e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipHeader_C", NULL));
312e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipOptions_C", NULL));
322e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipOptions_C", NULL));
332e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipInfo_C", NULL));
342e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipInfo_C", NULL));
352e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetInfoPointer_C", NULL));
362e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", NULL));
372e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", NULL));
382e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", NULL));
392e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", NULL));
402e956fe4SStefano Zampini #if defined(PETSC_HAVE_MPIIO)
412e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetUseMPIIO_C", NULL));
422e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetUseMPIIO_C", NULL));
432e956fe4SStefano Zampini #endif
442e956fe4SStefano Zampini   PetscFunctionReturn(0);
452e956fe4SStefano Zampini }
462e956fe4SStefano Zampini 
47cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
48*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySyncMPIIO(PetscViewer viewer)
49*d71ae5a4SJacob Faibussowitsch {
50cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
51cc843e7aSLisandro Dalcin 
52cc843e7aSLisandro Dalcin   PetscFunctionBegin;
53cc843e7aSLisandro Dalcin   if (vbinary->filemode == FILE_MODE_READ) PetscFunctionReturn(0);
5448a46eb9SPierre Jolivet   if (vbinary->mfsub != MPI_FILE_NULL) PetscCallMPI(MPI_File_sync(vbinary->mfsub));
55cc843e7aSLisandro Dalcin   if (vbinary->mfdes != MPI_FILE_NULL) {
569566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer)));
579566063dSJacob Faibussowitsch     PetscCallMPI(MPI_File_sync(vbinary->mfdes));
58cc843e7aSLisandro Dalcin   }
59cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
60cc843e7aSLisandro Dalcin }
61cc843e7aSLisandro Dalcin #endif
62cc843e7aSLisandro Dalcin 
63*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerGetSubViewer_Binary(PetscViewer viewer, MPI_Comm comm, PetscViewer *outviewer)
64*d71ae5a4SJacob Faibussowitsch {
65e8a65b7dSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
66cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
675c6c1daeSBarry Smith 
685c6c1daeSBarry Smith   PetscFunctionBegin;
699566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
70e8a65b7dSLisandro Dalcin 
71e8a65b7dSLisandro Dalcin   /* Return subviewer in process zero */
729566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
73dd400576SPatrick Sanan   if (rank == 0) {
74e5afcf28SBarry Smith     PetscMPIInt flg;
75e5afcf28SBarry Smith 
769566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_compare(PETSC_COMM_SELF, comm, &flg));
77cc73adaaSBarry Smith     PetscCheck(flg == MPI_IDENT || flg == MPI_CONGRUENT, PETSC_COMM_SELF, PETSC_ERR_SUP, "PetscViewerGetSubViewer() for PETSCVIEWERBINARY requires a singleton MPI_Comm");
789566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(comm, outviewer));
799566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(*outviewer, PETSCVIEWERBINARY));
809566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy((*outviewer)->data, vbinary, sizeof(PetscViewer_Binary)));
8112f4c3a9SLisandro Dalcin     (*outviewer)->setupcalled = PETSC_TRUE;
8296509d17SLisandro Dalcin   } else {
8396509d17SLisandro Dalcin     *outviewer = NULL;
8496509d17SLisandro Dalcin   }
85e8a65b7dSLisandro Dalcin 
86e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
87e8a65b7dSLisandro Dalcin   if (vbinary->usempiio && *outviewer) {
88e8a65b7dSLisandro Dalcin     PetscViewer_Binary *obinary = (PetscViewer_Binary *)(*outviewer)->data;
89e8a65b7dSLisandro Dalcin     /* Parent viewer opens a new MPI file handle on PETSC_COMM_SELF and keeps track of it for future reuse */
90e8a65b7dSLisandro Dalcin     if (vbinary->mfsub == MPI_FILE_NULL) {
91e8a65b7dSLisandro Dalcin       int amode;
92cc843e7aSLisandro Dalcin       switch (vbinary->filemode) {
93*d71ae5a4SJacob Faibussowitsch       case FILE_MODE_READ:
94*d71ae5a4SJacob Faibussowitsch         amode = MPI_MODE_RDONLY;
95*d71ae5a4SJacob Faibussowitsch         break;
96*d71ae5a4SJacob Faibussowitsch       case FILE_MODE_WRITE:
97*d71ae5a4SJacob Faibussowitsch         amode = MPI_MODE_WRONLY;
98*d71ae5a4SJacob Faibussowitsch         break;
99*d71ae5a4SJacob Faibussowitsch       case FILE_MODE_APPEND:
100*d71ae5a4SJacob Faibussowitsch         amode = MPI_MODE_WRONLY;
101*d71ae5a4SJacob Faibussowitsch         break;
102*d71ae5a4SJacob Faibussowitsch       default:
103*d71ae5a4SJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[vbinary->filemode]);
104e8a65b7dSLisandro Dalcin       }
1059566063dSJacob Faibussowitsch       PetscCallMPI(MPI_File_open(PETSC_COMM_SELF, vbinary->filename, amode, MPI_INFO_NULL, &vbinary->mfsub));
106e8a65b7dSLisandro Dalcin     }
107e8a65b7dSLisandro Dalcin     /* Subviewer gets the MPI file handle on PETSC_COMM_SELF */
108e8a65b7dSLisandro Dalcin     obinary->mfdes = vbinary->mfsub;
109e8a65b7dSLisandro Dalcin     obinary->mfsub = MPI_FILE_NULL;
110e8a65b7dSLisandro Dalcin     obinary->moff  = vbinary->moff;
111e8a65b7dSLisandro Dalcin   }
112e8a65b7dSLisandro Dalcin #endif
113cc843e7aSLisandro Dalcin 
114cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
1159566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinarySyncMPIIO(viewer));
116cc843e7aSLisandro Dalcin #endif
1175c6c1daeSBarry Smith   PetscFunctionReturn(0);
1185c6c1daeSBarry Smith }
1195c6c1daeSBarry Smith 
120*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerRestoreSubViewer_Binary(PetscViewer viewer, MPI_Comm comm, PetscViewer *outviewer)
121*d71ae5a4SJacob Faibussowitsch {
122e8a65b7dSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
123cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
124e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
125e8a65b7dSLisandro Dalcin   MPI_Offset moff = 0;
126e8a65b7dSLisandro Dalcin #endif
1275c6c1daeSBarry Smith 
1285c6c1daeSBarry Smith   PetscFunctionBegin;
1299566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
130c5853193SPierre Jolivet   PetscCheck(rank == 0 || !*outviewer, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Subviewer not obtained from viewer");
131e8a65b7dSLisandro Dalcin 
132e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
133e8a65b7dSLisandro Dalcin   if (vbinary->usempiio && *outviewer) {
134e8a65b7dSLisandro Dalcin     PetscViewer_Binary *obinary = (PetscViewer_Binary *)(*outviewer)->data;
13508401ef6SPierre Jolivet     PetscCheck(obinary->mfdes == vbinary->mfsub, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Subviewer not obtained from viewer");
1369566063dSJacob Faibussowitsch     if (obinary->mfsub != MPI_FILE_NULL) PetscCallMPI(MPI_File_close(&obinary->mfsub));
137e8a65b7dSLisandro Dalcin     moff = obinary->moff;
138e8a65b7dSLisandro Dalcin   }
139e8a65b7dSLisandro Dalcin #endif
140e8a65b7dSLisandro Dalcin 
141e8a65b7dSLisandro Dalcin   if (*outviewer) {
142e8a65b7dSLisandro Dalcin     PetscViewer_Binary *obinary = (PetscViewer_Binary *)(*outviewer)->data;
14308401ef6SPierre Jolivet     PetscCheck(obinary->fdes == vbinary->fdes, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Subviewer not obtained from viewer");
1449566063dSJacob Faibussowitsch     PetscCall(PetscFree((*outviewer)->data));
1452e956fe4SStefano Zampini     PetscCall(PetscViewerBinaryClearFunctionList(*outviewer));
1469566063dSJacob Faibussowitsch     PetscCall(PetscHeaderDestroy(outviewer));
1475c6c1daeSBarry Smith   }
148e8a65b7dSLisandro Dalcin 
149e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
150e8a65b7dSLisandro Dalcin   if (vbinary->usempiio) {
151e8a65b7dSLisandro Dalcin     PetscInt64 ioff = (PetscInt64)moff; /* We could use MPI_OFFSET datatype (requires MPI 2.2) */
1529566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&ioff, 1, MPIU_INT64, 0, PetscObjectComm((PetscObject)viewer)));
153e8a65b7dSLisandro Dalcin     vbinary->moff = (MPI_Offset)ioff;
154e8a65b7dSLisandro Dalcin   }
155e8a65b7dSLisandro Dalcin #endif
156cc843e7aSLisandro Dalcin 
157cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
1589566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinarySyncMPIIO(viewer));
159cc843e7aSLisandro Dalcin #endif
1605c6c1daeSBarry Smith   PetscFunctionReturn(0);
1615c6c1daeSBarry Smith }
1625c6c1daeSBarry Smith 
1635c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
1645c6c1daeSBarry Smith /*@C
16585b8072bSPatrick Sanan     PetscViewerBinaryGetMPIIOOffset - Gets the current global offset that should be passed to `MPI_File_set_view()` or `MPI_File_{write|read}_at[_all]()`
1665c6c1daeSBarry Smith 
1675c6c1daeSBarry Smith     Not Collective
1685c6c1daeSBarry Smith 
1695c6c1daeSBarry Smith     Input Parameter:
17085b8072bSPatrick Sanan .   viewer - PetscViewer context, obtained from `PetscViewerBinaryOpen()`
1715c6c1daeSBarry Smith 
1725c6c1daeSBarry Smith     Output Parameter:
17322a8f86cSLisandro Dalcin .   off - the current global offset
1745c6c1daeSBarry Smith 
1755c6c1daeSBarry Smith     Level: advanced
1765c6c1daeSBarry Smith 
177811af0c4SBarry Smith     Note:
178811af0c4SBarry Smith     Use `PetscViewerBinaryAddMPIIOOffset()` to increase this value after you have written a view.
179811af0c4SBarry Smith 
1805c6c1daeSBarry Smith     Fortran Note:
1815c6c1daeSBarry Smith     This routine is not supported in Fortran.
1825c6c1daeSBarry Smith 
183811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryAddMPIIOOffset()`
1845c6c1daeSBarry Smith @*/
185*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetMPIIOOffset(PetscViewer viewer, MPI_Offset *off)
186*d71ae5a4SJacob 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;
1945c6c1daeSBarry Smith   PetscFunctionReturn(0);
1955c6c1daeSBarry Smith }
1965c6c1daeSBarry Smith 
1975c6c1daeSBarry Smith /*@C
19822a8f86cSLisandro Dalcin     PetscViewerBinaryAddMPIIOOffset - Adds to the current global offset
1995c6c1daeSBarry Smith 
20022a8f86cSLisandro Dalcin     Logically Collective
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 
2115c6c1daeSBarry Smith     Fortran Note:
2125c6c1daeSBarry Smith     This routine is not supported in Fortran.
2135c6c1daeSBarry Smith 
214811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
2155c6c1daeSBarry Smith @*/
216*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryAddMPIIOOffset(PetscViewer viewer, MPI_Offset off)
217*d71ae5a4SJacob Faibussowitsch {
21822a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
2195c6c1daeSBarry Smith 
2205c6c1daeSBarry Smith   PetscFunctionBegin;
22122a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
22222a8f86cSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, (PetscInt)off, 2);
22322a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
2245c6c1daeSBarry Smith   vbinary->moff += off;
2255c6c1daeSBarry Smith   PetscFunctionReturn(0);
2265c6c1daeSBarry Smith }
2275c6c1daeSBarry Smith 
2285c6c1daeSBarry Smith /*@C
229811af0c4SBarry Smith     PetscViewerBinaryGetMPIIODescriptor - Extracts the MPI IO file descriptor from a `PetscViewer`.
2305c6c1daeSBarry Smith 
2315c6c1daeSBarry Smith     Not Collective
2325c6c1daeSBarry Smith 
2335c6c1daeSBarry Smith     Input Parameter:
234811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
2355c6c1daeSBarry Smith 
2365c6c1daeSBarry Smith     Output Parameter:
2375c6c1daeSBarry Smith .   fdes - file descriptor
2385c6c1daeSBarry Smith 
2395c6c1daeSBarry Smith     Level: advanced
2405c6c1daeSBarry Smith 
2415c6c1daeSBarry Smith     Fortran Note:
2425c6c1daeSBarry Smith     This routine is not supported in Fortran.
2435c6c1daeSBarry Smith 
244811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
2455c6c1daeSBarry Smith @*/
246*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetMPIIODescriptor(PetscViewer viewer, MPI_File *fdes)
247*d71ae5a4SJacob Faibussowitsch {
24822a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
2495c6c1daeSBarry Smith 
2505c6c1daeSBarry Smith   PetscFunctionBegin;
25122a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
25222a8f86cSLisandro Dalcin   PetscValidPointer(fdes, 2);
2539566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
25422a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
2555c6c1daeSBarry Smith   *fdes   = vbinary->mfdes;
2565c6c1daeSBarry Smith   PetscFunctionReturn(0);
2575c6c1daeSBarry Smith }
258cc843e7aSLisandro Dalcin #endif
2595c6c1daeSBarry Smith 
260cc843e7aSLisandro Dalcin /*@
261cc843e7aSLisandro Dalcin     PetscViewerBinarySetUseMPIIO - Sets a binary viewer to use MPI-IO for reading/writing. Must be called
262811af0c4SBarry Smith         before `PetscViewerFileSetName()`
263cc843e7aSLisandro Dalcin 
264811af0c4SBarry Smith     Logically Collective on viewer
265cc843e7aSLisandro Dalcin 
266cc843e7aSLisandro Dalcin     Input Parameters:
267811af0c4SBarry Smith +   viewer - the `PetscViewer`; must be a `PETSCVIEWERBINARY`
268811af0c4SBarry Smith -   use - `PETSC_TRUE` means MPI-IO will be used
269cc843e7aSLisandro Dalcin 
270811af0c4SBarry Smith     Options Database Key:
271cc843e7aSLisandro Dalcin     -viewer_binary_mpiio : Flag for using MPI-IO
272cc843e7aSLisandro Dalcin 
273cc843e7aSLisandro Dalcin     Level: advanced
274cc843e7aSLisandro Dalcin 
275811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`,
276db781477SPatrick Sanan           `PetscViewerBinaryGetUseMPIIO()`
277cc843e7aSLisandro Dalcin @*/
278*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetUseMPIIO(PetscViewer viewer, PetscBool use)
279*d71ae5a4SJacob Faibussowitsch {
280a8aae444SDave May   PetscFunctionBegin;
281cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
282cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, use, 2);
283cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetUseMPIIO_C", (PetscViewer, PetscBool), (viewer, use));
284cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
285cc843e7aSLisandro Dalcin }
286cc843e7aSLisandro Dalcin 
287cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
288*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetUseMPIIO_Binary(PetscViewer viewer, PetscBool use)
289*d71ae5a4SJacob Faibussowitsch {
290cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
291cc843e7aSLisandro Dalcin   PetscFunctionBegin;
292cc73adaaSBarry Smith   PetscCheck(!viewer->setupcalled || vbinary->usempiio == use, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Cannot change MPIIO to %s after setup", PetscBools[use]);
293cc843e7aSLisandro Dalcin   vbinary->usempiio = use;
294a8aae444SDave May   PetscFunctionReturn(0);
295a8aae444SDave May }
296a8aae444SDave May #endif
297a8aae444SDave May 
298cc843e7aSLisandro Dalcin /*@
299bc196f7cSDave May     PetscViewerBinaryGetUseMPIIO - Returns PETSC_TRUE if the binary viewer uses MPI-IO.
3005c6c1daeSBarry Smith 
3015c6c1daeSBarry Smith     Not Collective
3025c6c1daeSBarry Smith 
3035c6c1daeSBarry Smith     Input Parameter:
304811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`; must be a `PETSCVIEWERBINARY`
3055c6c1daeSBarry Smith 
3065c6c1daeSBarry Smith     Output Parameter:
307811af0c4SBarry Smith .   use - `PETSC_TRUE` if MPI-IO is being used
3085c6c1daeSBarry Smith 
3095c6c1daeSBarry Smith     Level: advanced
3105c6c1daeSBarry Smith 
311bc196f7cSDave May     Note:
312bc196f7cSDave May     If MPI-IO is not available, this function will always return PETSC_FALSE
313bc196f7cSDave May 
314811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
3155c6c1daeSBarry Smith @*/
316*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetUseMPIIO(PetscViewer viewer, PetscBool *use)
317*d71ae5a4SJacob Faibussowitsch {
3185c6c1daeSBarry Smith   PetscFunctionBegin;
319cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
320cc843e7aSLisandro Dalcin   PetscValidBoolPointer(use, 2);
321cc843e7aSLisandro Dalcin   *use = PETSC_FALSE;
322cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinaryGetUseMPIIO_C", (PetscViewer, PetscBool *), (viewer, use));
3235c6c1daeSBarry Smith   PetscFunctionReturn(0);
3245c6c1daeSBarry Smith }
3255c6c1daeSBarry Smith 
326cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
327*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetUseMPIIO_Binary(PetscViewer viewer, PetscBool *use)
328*d71ae5a4SJacob Faibussowitsch {
3295c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
3305c6c1daeSBarry Smith 
3315c6c1daeSBarry Smith   PetscFunctionBegin;
332cc843e7aSLisandro Dalcin   *use = vbinary->usempiio;
333cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
334cc843e7aSLisandro Dalcin }
335cc843e7aSLisandro Dalcin #endif
336cc843e7aSLisandro Dalcin 
337cc843e7aSLisandro Dalcin /*@
338cc843e7aSLisandro Dalcin     PetscViewerBinarySetFlowControl - Sets how many messages are allowed to outstanding at the same time during parallel IO reads/writes
339cc843e7aSLisandro Dalcin 
340cc843e7aSLisandro Dalcin     Not Collective
341cc843e7aSLisandro Dalcin 
342d8d19677SJose E. Roman     Input Parameters:
343811af0c4SBarry Smith +   viewer - PetscViewer context, obtained from `PetscViewerBinaryOpen()`
344cc843e7aSLisandro Dalcin -   fc - the number of messages, defaults to 256 if this function was not called
345cc843e7aSLisandro Dalcin 
346cc843e7aSLisandro Dalcin     Level: advanced
347cc843e7aSLisandro Dalcin 
348811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetFlowControl()`
349cc843e7aSLisandro Dalcin @*/
350*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetFlowControl(PetscViewer viewer, PetscInt fc)
351*d71ae5a4SJacob Faibussowitsch {
352cc843e7aSLisandro Dalcin   PetscFunctionBegin;
353cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
354cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, fc, 2);
355cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetFlowControl_C", (PetscViewer, PetscInt), (viewer, fc));
3565c6c1daeSBarry Smith   PetscFunctionReturn(0);
3575c6c1daeSBarry Smith }
3585c6c1daeSBarry Smith 
359*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetFlowControl_Binary(PetscViewer viewer, PetscInt fc)
360*d71ae5a4SJacob Faibussowitsch {
361cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
362cc843e7aSLisandro Dalcin 
363cc843e7aSLisandro Dalcin   PetscFunctionBegin;
36408401ef6SPierre Jolivet   PetscCheck(fc > 1, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_OUTOFRANGE, "Flow control count must be greater than 1, %" PetscInt_FMT " was set", fc);
365cc843e7aSLisandro Dalcin   vbinary->flowcontrol = fc;
366cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
367cc843e7aSLisandro Dalcin }
368cc843e7aSLisandro Dalcin 
369cc843e7aSLisandro Dalcin /*@
3705c6c1daeSBarry Smith     PetscViewerBinaryGetFlowControl - Returns how many messages are allowed to outstanding at the same time during parallel IO reads/writes
3715c6c1daeSBarry Smith 
3725c6c1daeSBarry Smith     Not Collective
3735c6c1daeSBarry Smith 
3745c6c1daeSBarry Smith     Input Parameter:
375811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
3765c6c1daeSBarry Smith 
3775c6c1daeSBarry Smith     Output Parameter:
3785c6c1daeSBarry Smith .   fc - the number of messages
3795c6c1daeSBarry Smith 
3805c6c1daeSBarry Smith     Level: advanced
3815c6c1daeSBarry Smith 
382811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinarySetFlowControl()`
3835c6c1daeSBarry Smith @*/
384*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetFlowControl(PetscViewer viewer, PetscInt *fc)
385*d71ae5a4SJacob Faibussowitsch {
3865c6c1daeSBarry Smith   PetscFunctionBegin;
387cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
388cc843e7aSLisandro Dalcin   PetscValidIntPointer(fc, 2);
389cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetFlowControl_C", (PetscViewer, PetscInt *), (viewer, fc));
3905c6c1daeSBarry Smith   PetscFunctionReturn(0);
3915c6c1daeSBarry Smith }
3925c6c1daeSBarry Smith 
393*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetFlowControl_Binary(PetscViewer viewer, PetscInt *fc)
394*d71ae5a4SJacob Faibussowitsch {
3955c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
3965c6c1daeSBarry Smith 
3975c6c1daeSBarry Smith   PetscFunctionBegin;
398cc843e7aSLisandro Dalcin   *fc = vbinary->flowcontrol;
3995c6c1daeSBarry Smith   PetscFunctionReturn(0);
4005c6c1daeSBarry Smith }
4015c6c1daeSBarry Smith 
4025c6c1daeSBarry Smith /*@C
403811af0c4SBarry Smith     PetscViewerBinaryGetDescriptor - Extracts the file descriptor from a `PetscViewer` of `PetscViewerType` `PETSCVIEWERBINARY`.
4045c6c1daeSBarry Smith 
405811af0c4SBarry Smith     Collective on viewer because it may trigger a `PetscViewerSetUp()` call
4065c6c1daeSBarry Smith 
4075c6c1daeSBarry Smith     Input Parameter:
408811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
4095c6c1daeSBarry Smith 
4105c6c1daeSBarry Smith     Output Parameter:
4115c6c1daeSBarry Smith .   fdes - file descriptor
4125c6c1daeSBarry Smith 
4135c6c1daeSBarry Smith     Level: advanced
4145c6c1daeSBarry Smith 
415811af0c4SBarry Smith     Note:
416811af0c4SBarry Smith       For writable binary `PetscViewer`s, the descriptor will only be valid for the
417811af0c4SBarry Smith     first processor in the communicator that shares the `PetscViewer`. For readable
4185c6c1daeSBarry Smith     files it will only be valid on nodes that have the file. If node 0 does not
4195c6c1daeSBarry Smith     have the file it generates an error even if another node does have the file.
4205c6c1daeSBarry Smith 
4215c6c1daeSBarry Smith     Fortran Note:
4225c6c1daeSBarry Smith     This routine is not supported in Fortran.
4235c6c1daeSBarry Smith 
424811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`
4255c6c1daeSBarry Smith @*/
426*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetDescriptor(PetscViewer viewer, int *fdes)
427*d71ae5a4SJacob Faibussowitsch {
42822a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
4295c6c1daeSBarry Smith 
4305c6c1daeSBarry Smith   PetscFunctionBegin;
43122a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
43222a8f86cSLisandro Dalcin   PetscValidPointer(fdes, 2);
4339566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
43422a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
4355c6c1daeSBarry Smith   *fdes   = vbinary->fdes;
4365c6c1daeSBarry Smith   PetscFunctionReturn(0);
4375c6c1daeSBarry Smith }
4385c6c1daeSBarry Smith 
4395c6c1daeSBarry Smith /*@
4405c6c1daeSBarry Smith     PetscViewerBinarySkipInfo - Binary file will not have .info file created with it
4415c6c1daeSBarry Smith 
4425c6c1daeSBarry Smith     Not Collective
4435c6c1daeSBarry Smith 
444fd292e60Sprj-     Input Parameter:
445811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerCreate()`
4465c6c1daeSBarry Smith 
4475c6c1daeSBarry Smith     Options Database Key:
44810699b91SBarry Smith .   -viewer_binary_skip_info - true indicates do not generate .info file
4495c6c1daeSBarry Smith 
4505c6c1daeSBarry Smith     Level: advanced
4515c6c1daeSBarry Smith 
45295452b02SPatrick Sanan     Notes:
453811af0c4SBarry Smith     This must be called after `PetscViewerSetType()`. If you use `PetscViewerBinaryOpen()` then
4545c6c1daeSBarry Smith     you can only skip the info file with the -viewer_binary_skip_info flag. To use the function you must open the
455811af0c4SBarry Smith     viewer with `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinarySkipInfo()`.
4565c6c1daeSBarry Smith 
4575c6c1daeSBarry Smith     The .info contains meta information about the data in the binary file, for example the block size if it was
4585c6c1daeSBarry Smith     set for a vector or matrix.
4595c6c1daeSBarry Smith 
460811af0c4SBarry Smith     This routine is deprecated, use `PetscViewerBinarySetSkipInfo()`
461811af0c4SBarry Smith 
462811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySetSkipOptions()`,
463db781477SPatrick Sanan           `PetscViewerBinaryGetSkipOptions()`, `PetscViewerBinaryGetSkipInfo()`
4645c6c1daeSBarry Smith @*/
465*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySkipInfo(PetscViewer viewer)
466*d71ae5a4SJacob Faibussowitsch {
4675c6c1daeSBarry Smith   PetscFunctionBegin;
4689566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinarySetSkipInfo(viewer, PETSC_TRUE));
469807ea322SDave May   PetscFunctionReturn(0);
470807ea322SDave May }
471807ea322SDave May 
472807ea322SDave May /*@
473807ea322SDave May     PetscViewerBinarySetSkipInfo - Binary file will not have .info file created with it
474807ea322SDave May 
475807ea322SDave May     Not Collective
476807ea322SDave May 
477d8d19677SJose E. Roman     Input Parameters:
478cc843e7aSLisandro Dalcin +   viewer - PetscViewer context, obtained from PetscViewerCreate()
479cc843e7aSLisandro Dalcin -   skip - PETSC_TRUE implies the .info file will not be generated
480807ea322SDave May 
481807ea322SDave May     Options Database Key:
48210699b91SBarry Smith .   -viewer_binary_skip_info - true indicates do not generate .info file
483807ea322SDave May 
484807ea322SDave May     Level: advanced
485807ea322SDave May 
486811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySetSkipOptions()`,
487db781477SPatrick Sanan           `PetscViewerBinaryGetSkipOptions()`, `PetscViewerBinaryGetSkipInfo()`
488807ea322SDave May @*/
489*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetSkipInfo(PetscViewer viewer, PetscBool skip)
490*d71ae5a4SJacob Faibussowitsch {
491807ea322SDave May   PetscFunctionBegin;
492cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
493cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, skip, 2);
494cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetSkipInfo_C", (PetscViewer, PetscBool), (viewer, skip));
495807ea322SDave May   PetscFunctionReturn(0);
496807ea322SDave May }
497807ea322SDave May 
498*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetSkipInfo_Binary(PetscViewer viewer, PetscBool skip)
499*d71ae5a4SJacob Faibussowitsch {
500807ea322SDave May   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
501807ea322SDave May 
502807ea322SDave May   PetscFunctionBegin;
503cc843e7aSLisandro Dalcin   vbinary->skipinfo = skip;
504807ea322SDave May   PetscFunctionReturn(0);
505807ea322SDave May }
506807ea322SDave May 
507807ea322SDave May /*@
508807ea322SDave May     PetscViewerBinaryGetSkipInfo - check if viewer wrote a .info file
509807ea322SDave May 
510807ea322SDave May     Not Collective
511807ea322SDave May 
512807ea322SDave May     Input Parameter:
513811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
514807ea322SDave May 
515807ea322SDave May     Output Parameter:
516811af0c4SBarry Smith .   skip - `PETSC_TRUE` implies the .info file was not generated
517807ea322SDave May 
518807ea322SDave May     Level: advanced
519807ea322SDave May 
520811af0c4SBarry Smith     Note:
521811af0c4SBarry Smith     This must be called after `PetscViewerSetType()`
522807ea322SDave May 
523811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`,
524db781477SPatrick Sanan           `PetscViewerBinarySetSkipOptions()`, `PetscViewerBinarySetSkipInfo()`
525807ea322SDave May @*/
526*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetSkipInfo(PetscViewer viewer, PetscBool *skip)
527*d71ae5a4SJacob Faibussowitsch {
528807ea322SDave May   PetscFunctionBegin;
529cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
530cc843e7aSLisandro Dalcin   PetscValidBoolPointer(skip, 2);
531cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetSkipInfo_C", (PetscViewer, PetscBool *), (viewer, skip));
532807ea322SDave May   PetscFunctionReturn(0);
533807ea322SDave May }
534807ea322SDave May 
535*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetSkipInfo_Binary(PetscViewer viewer, PetscBool *skip)
536*d71ae5a4SJacob Faibussowitsch {
537807ea322SDave May   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
538807ea322SDave May 
539807ea322SDave May   PetscFunctionBegin;
540cc843e7aSLisandro Dalcin   *skip = vbinary->skipinfo;
541807ea322SDave May   PetscFunctionReturn(0);
542807ea322SDave May }
543807ea322SDave May 
5445c6c1daeSBarry Smith /*@
5455c6c1daeSBarry Smith     PetscViewerBinarySetSkipOptions - do not use the PETSc options database when loading objects
5465c6c1daeSBarry Smith 
5475c6c1daeSBarry Smith     Not Collective
5485c6c1daeSBarry Smith 
5495c6c1daeSBarry Smith     Input Parameters:
550811af0c4SBarry Smith +   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
551811af0c4SBarry Smith -   skip - `PETSC_TRUE` means do not use the options from the options database
5525c6c1daeSBarry Smith 
5535c6c1daeSBarry Smith     Options Database Key:
554811af0c4SBarry Smith .   -viewer_binary_skip_options <true or false> - true means do not use the options from the options database
5555c6c1daeSBarry Smith 
5565c6c1daeSBarry Smith     Level: advanced
5575c6c1daeSBarry Smith 
558811af0c4SBarry Smith     Note:
559811af0c4SBarry Smith     This must be called after `PetscViewerSetType()`
5605c6c1daeSBarry Smith 
561811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
562db781477SPatrick Sanan           `PetscViewerBinaryGetSkipOptions()`
5635c6c1daeSBarry Smith @*/
564*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetSkipOptions(PetscViewer viewer, PetscBool skip)
565*d71ae5a4SJacob Faibussowitsch {
5665c6c1daeSBarry Smith   PetscFunctionBegin;
567cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
568cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, skip, 2);
569cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetSkipOptions_C", (PetscViewer, PetscBool), (viewer, skip));
570807ea322SDave May   PetscFunctionReturn(0);
571807ea322SDave May }
572807ea322SDave May 
573*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetSkipOptions_Binary(PetscViewer viewer, PetscBool skip)
574*d71ae5a4SJacob Faibussowitsch {
575cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
576807ea322SDave May 
577807ea322SDave May   PetscFunctionBegin;
578cc843e7aSLisandro Dalcin   vbinary->skipoptions = skip;
5795c6c1daeSBarry Smith   PetscFunctionReturn(0);
5805c6c1daeSBarry Smith }
5815c6c1daeSBarry Smith 
5825c6c1daeSBarry Smith /*@
5835c6c1daeSBarry Smith     PetscViewerBinaryGetSkipOptions - checks if viewer uses the PETSc options database when loading objects
5845c6c1daeSBarry Smith 
5855c6c1daeSBarry Smith     Not Collective
5865c6c1daeSBarry Smith 
5875c6c1daeSBarry Smith     Input Parameter:
588811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
5895c6c1daeSBarry Smith 
5905c6c1daeSBarry Smith     Output Parameter:
591811af0c4SBarry Smith .   skip - `PETSC_TRUE` means do not use
5925c6c1daeSBarry Smith 
5935c6c1daeSBarry Smith     Level: advanced
5945c6c1daeSBarry Smith 
595811af0c4SBarry Smith     Note:
596811af0c4SBarry Smith     This must be called after `PetscViewerSetType()`
5975c6c1daeSBarry Smith 
598811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
599db781477SPatrick Sanan           `PetscViewerBinarySetSkipOptions()`
6005c6c1daeSBarry Smith @*/
601*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetSkipOptions(PetscViewer viewer, PetscBool *skip)
602*d71ae5a4SJacob Faibussowitsch {
6035c6c1daeSBarry Smith   PetscFunctionBegin;
604cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
605cc843e7aSLisandro Dalcin   PetscValidBoolPointer(skip, 2);
606cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetSkipOptions_C", (PetscViewer, PetscBool *), (viewer, skip));
6075c6c1daeSBarry Smith   PetscFunctionReturn(0);
6085c6c1daeSBarry Smith }
6095c6c1daeSBarry Smith 
610*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetSkipOptions_Binary(PetscViewer viewer, PetscBool *skip)
611*d71ae5a4SJacob Faibussowitsch {
612d21b9a37SPierre Jolivet   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
6135c6c1daeSBarry Smith 
6145c6c1daeSBarry Smith   PetscFunctionBegin;
615cc843e7aSLisandro Dalcin   *skip = vbinary->skipoptions;
6165c6c1daeSBarry Smith   PetscFunctionReturn(0);
6175c6c1daeSBarry Smith }
6185c6c1daeSBarry Smith 
6195c6c1daeSBarry Smith /*@
6205c6c1daeSBarry Smith     PetscViewerBinarySetSkipHeader - do not write a header with size information on output, just raw data
6215c6c1daeSBarry Smith 
6225c6c1daeSBarry Smith     Not Collective
6235c6c1daeSBarry Smith 
6245c6c1daeSBarry Smith     Input Parameters:
625811af0c4SBarry Smith +   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
626811af0c4SBarry Smith -   skip - `PETSC_TRUE `means do not write header
6275c6c1daeSBarry Smith 
6285c6c1daeSBarry Smith     Options Database Key:
629811af0c4SBarry Smith .   -viewer_binary_skip_header <true or false> - true means do not write header
6305c6c1daeSBarry Smith 
6315c6c1daeSBarry Smith     Level: advanced
6325c6c1daeSBarry Smith 
63395452b02SPatrick Sanan     Notes:
634811af0c4SBarry Smith       This must be called after `PetscViewerSetType()`
6355c6c1daeSBarry Smith 
63610699b91SBarry Smith       Is ignored on anything but a binary viewer
6375c6c1daeSBarry Smith 
638811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
639db781477SPatrick Sanan           `PetscViewerBinaryGetSkipHeader()`
6405c6c1daeSBarry Smith @*/
641*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetSkipHeader(PetscViewer viewer, PetscBool skip)
642*d71ae5a4SJacob Faibussowitsch {
6435c6c1daeSBarry Smith   PetscFunctionBegin;
644cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
645cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, skip, 2);
646cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetSkipHeader_C", (PetscViewer, PetscBool), (viewer, skip));
6475c6c1daeSBarry Smith   PetscFunctionReturn(0);
6485c6c1daeSBarry Smith }
6495c6c1daeSBarry Smith 
650*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetSkipHeader_Binary(PetscViewer viewer, PetscBool skip)
651*d71ae5a4SJacob Faibussowitsch {
6525c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
6535c6c1daeSBarry Smith 
6545c6c1daeSBarry Smith   PetscFunctionBegin;
655cc843e7aSLisandro Dalcin   vbinary->skipheader = skip;
6565c6c1daeSBarry Smith   PetscFunctionReturn(0);
6575c6c1daeSBarry Smith }
6585c6c1daeSBarry Smith 
6595c6c1daeSBarry Smith /*@
6605c6c1daeSBarry Smith     PetscViewerBinaryGetSkipHeader - checks whether to write a header with size information on output, or just raw data
6615c6c1daeSBarry Smith 
6625c6c1daeSBarry Smith     Not Collective
6635c6c1daeSBarry Smith 
6645c6c1daeSBarry Smith     Input Parameter:
665811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
6665c6c1daeSBarry Smith 
6675c6c1daeSBarry Smith     Output Parameter:
668811af0c4SBarry Smith .   skip - `PETSC_TRUE` means do not write header
6695c6c1daeSBarry Smith 
6705c6c1daeSBarry Smith     Level: advanced
6715c6c1daeSBarry Smith 
67295452b02SPatrick Sanan     Notes:
67395452b02SPatrick Sanan     This must be called after PetscViewerSetType()
6745c6c1daeSBarry Smith 
675811af0c4SBarry Smith     Returns `PETSC_FALSE` for `PETSCSOCKETVIEWER`, you cannot skip the header for it.
6765c6c1daeSBarry Smith 
677811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
678db781477SPatrick Sanan           `PetscViewerBinarySetSkipHeader()`
6795c6c1daeSBarry Smith @*/
680*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetSkipHeader(PetscViewer viewer, PetscBool *skip)
681*d71ae5a4SJacob Faibussowitsch {
6825c6c1daeSBarry Smith   PetscFunctionBegin;
683cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
684cc843e7aSLisandro Dalcin   PetscValidBoolPointer(skip, 2);
685cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetSkipHeader_C", (PetscViewer, PetscBool *), (viewer, skip));
6865c6c1daeSBarry Smith   PetscFunctionReturn(0);
6875c6c1daeSBarry Smith }
6885c6c1daeSBarry Smith 
689*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetSkipHeader_Binary(PetscViewer viewer, PetscBool *skip)
690*d71ae5a4SJacob Faibussowitsch {
691f3b211e4SSatish Balay   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
6925c6c1daeSBarry Smith 
6935c6c1daeSBarry Smith   PetscFunctionBegin;
694cc843e7aSLisandro Dalcin   *skip = vbinary->skipheader;
6955c6c1daeSBarry Smith   PetscFunctionReturn(0);
6965c6c1daeSBarry Smith }
6975c6c1daeSBarry Smith 
6985c6c1daeSBarry Smith /*@C
6995c6c1daeSBarry Smith     PetscViewerBinaryGetInfoPointer - Extracts the file pointer for the ASCII
7005c6c1daeSBarry Smith           info file associated with a binary file.
7015c6c1daeSBarry Smith 
7025c6c1daeSBarry Smith     Not Collective
7035c6c1daeSBarry Smith 
7045c6c1daeSBarry Smith     Input Parameter:
705811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
7065c6c1daeSBarry Smith 
7075c6c1daeSBarry Smith     Output Parameter:
7080298fd71SBarry Smith .   file - file pointer  Always returns NULL if not a binary viewer
7095c6c1daeSBarry Smith 
7105c6c1daeSBarry Smith     Level: advanced
7115c6c1daeSBarry Smith 
712811af0c4SBarry Smith     Note:
713811af0c4SBarry Smith       For writable binary `PetscViewer`s, the file pointer will only be valid for the
714811af0c4SBarry Smith     first processor in the communicator that shares the `PetscViewer`.
7155c6c1daeSBarry Smith 
7165c6c1daeSBarry Smith     Fortran Note:
7175c6c1daeSBarry Smith     This routine is not supported in Fortran.
7185c6c1daeSBarry Smith 
719811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`
7205c6c1daeSBarry Smith @*/
721*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetInfoPointer(PetscViewer viewer, FILE **file)
722*d71ae5a4SJacob Faibussowitsch {
7235c6c1daeSBarry Smith   PetscFunctionBegin;
724cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
725cc843e7aSLisandro Dalcin   PetscValidPointer(file, 2);
7260298fd71SBarry Smith   *file = NULL;
727cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinaryGetInfoPointer_C", (PetscViewer, FILE **), (viewer, file));
7285c6c1daeSBarry Smith   PetscFunctionReturn(0);
7295c6c1daeSBarry Smith }
7305c6c1daeSBarry Smith 
731*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetInfoPointer_Binary(PetscViewer viewer, FILE **file)
732*d71ae5a4SJacob Faibussowitsch {
733cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
7345c6c1daeSBarry Smith 
7355c6c1daeSBarry Smith   PetscFunctionBegin;
7369566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
737cc843e7aSLisandro Dalcin   *file = vbinary->fdes_info;
738cc843e7aSLisandro Dalcin   if (viewer->format == PETSC_VIEWER_BINARY_MATLAB && !vbinary->matlabheaderwritten) {
7395c6c1daeSBarry Smith     if (vbinary->fdes_info) {
740cc843e7aSLisandro Dalcin       FILE *info = vbinary->fdes_info;
7419566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- begin code written by PetscViewerBinary for MATLAB format ---#\n"));
7429566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#$$ Set.filename = '%s';\n", vbinary->filename));
7439566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#$$ fd = PetscOpenFile(Set.filename);\n"));
7449566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n"));
745cc843e7aSLisandro Dalcin     }
746cc843e7aSLisandro Dalcin     vbinary->matlabheaderwritten = PETSC_TRUE;
7475c6c1daeSBarry Smith   }
7485c6c1daeSBarry Smith   PetscFunctionReturn(0);
7495c6c1daeSBarry Smith }
7505c6c1daeSBarry Smith 
751e0385b85SDave May #if defined(PETSC_HAVE_MPIIO)
752*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_BinaryMPIIO(PetscViewer v)
753*d71ae5a4SJacob Faibussowitsch {
754e0385b85SDave May   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
755e0385b85SDave May 
756e0385b85SDave May   PetscFunctionBegin;
75748a46eb9SPierre Jolivet   if (vbinary->mfdes != MPI_FILE_NULL) PetscCallMPI(MPI_File_close(&vbinary->mfdes));
75848a46eb9SPierre Jolivet   if (vbinary->mfsub != MPI_FILE_NULL) PetscCallMPI(MPI_File_close(&vbinary->mfsub));
759cc843e7aSLisandro Dalcin   vbinary->moff = 0;
760e0385b85SDave May   PetscFunctionReturn(0);
761e0385b85SDave May }
762e0385b85SDave May #endif
763e0385b85SDave May 
764*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_BinarySTDIO(PetscViewer v)
765*d71ae5a4SJacob Faibussowitsch {
766cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
767cc843e7aSLisandro Dalcin 
768cc843e7aSLisandro Dalcin   PetscFunctionBegin;
769cc843e7aSLisandro Dalcin   if (vbinary->fdes != -1) {
7709566063dSJacob Faibussowitsch     PetscCall(PetscBinaryClose(vbinary->fdes));
771cc843e7aSLisandro Dalcin     vbinary->fdes = -1;
772cc843e7aSLisandro Dalcin     if (vbinary->storecompressed) {
773cc843e7aSLisandro Dalcin       char        cmd[8 + PETSC_MAX_PATH_LEN], out[64 + PETSC_MAX_PATH_LEN] = "";
774cc843e7aSLisandro Dalcin       const char *gzfilename = vbinary->ogzfilename ? vbinary->ogzfilename : vbinary->filename;
775cc843e7aSLisandro Dalcin       /* compress the file */
7769566063dSJacob Faibussowitsch       PetscCall(PetscStrncpy(cmd, "gzip -f ", sizeof(cmd)));
7779566063dSJacob Faibussowitsch       PetscCall(PetscStrlcat(cmd, gzfilename, sizeof(cmd)));
778cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_POPEN)
779cc843e7aSLisandro Dalcin       {
780cc843e7aSLisandro Dalcin         FILE *fp;
7819566063dSJacob Faibussowitsch         PetscCall(PetscPOpen(PETSC_COMM_SELF, NULL, cmd, "r", &fp));
782cc73adaaSBarry Smith         PetscCheck(!fgets(out, (int)(sizeof(out) - 1), fp), PETSC_COMM_SELF, PETSC_ERR_LIB, "Error from command %s\n%s", cmd, out);
7839566063dSJacob Faibussowitsch         PetscCall(PetscPClose(PETSC_COMM_SELF, fp));
784cc843e7aSLisandro Dalcin       }
785cc843e7aSLisandro Dalcin #endif
786cc843e7aSLisandro Dalcin     }
787cc843e7aSLisandro Dalcin   }
7889566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary->ogzfilename));
789cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
790cc843e7aSLisandro Dalcin }
791cc843e7aSLisandro Dalcin 
792*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_BinaryInfo(PetscViewer v)
793*d71ae5a4SJacob Faibussowitsch {
794cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
795cc843e7aSLisandro Dalcin 
796cc843e7aSLisandro Dalcin   PetscFunctionBegin;
797cc843e7aSLisandro Dalcin   if (v->format == PETSC_VIEWER_BINARY_MATLAB && vbinary->matlabheaderwritten) {
798cc843e7aSLisandro Dalcin     if (vbinary->fdes_info) {
799cc843e7aSLisandro Dalcin       FILE *info = vbinary->fdes_info;
8009566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- begin code written by PetscViewerBinary for MATLAB format ---#\n"));
8019566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#$$ close(fd);\n"));
8029566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n"));
803cc843e7aSLisandro Dalcin     }
804cc843e7aSLisandro Dalcin   }
805cc843e7aSLisandro Dalcin   if (vbinary->fdes_info) {
806cc843e7aSLisandro Dalcin     FILE *info         = vbinary->fdes_info;
807cc843e7aSLisandro Dalcin     vbinary->fdes_info = NULL;
808cc73adaaSBarry Smith     PetscCheck(!fclose(info), PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file");
809cc843e7aSLisandro Dalcin   }
810cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
811cc843e7aSLisandro Dalcin }
812cc843e7aSLisandro Dalcin 
813*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_Binary(PetscViewer v)
814*d71ae5a4SJacob Faibussowitsch {
815cc843e7aSLisandro Dalcin   PetscFunctionBegin;
816cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
8179566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_BinaryMPIIO(v));
818cc843e7aSLisandro Dalcin #endif
8199566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_BinarySTDIO(v));
8209566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_BinaryInfo(v));
821cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
822cc843e7aSLisandro Dalcin }
823cc843e7aSLisandro Dalcin 
824*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerDestroy_Binary(PetscViewer v)
825*d71ae5a4SJacob Faibussowitsch {
8265c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
8275c6c1daeSBarry Smith 
8285c6c1daeSBarry Smith   PetscFunctionBegin;
8299566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_Binary(v));
8309566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary->filename));
8319566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary));
8322e956fe4SStefano Zampini   PetscCall(PetscViewerBinaryClearFunctionList(v));
833e0385b85SDave May   PetscFunctionReturn(0);
834e0385b85SDave May }
8355c6c1daeSBarry Smith 
8365c6c1daeSBarry Smith /*@C
8375c6c1daeSBarry Smith    PetscViewerBinaryOpen - Opens a file for binary input/output.
8385c6c1daeSBarry Smith 
839d083f849SBarry Smith    Collective
8405c6c1daeSBarry Smith 
8415c6c1daeSBarry Smith    Input Parameters:
8425c6c1daeSBarry Smith +  comm - MPI communicator
8435c6c1daeSBarry Smith .  name - name of file
844cc843e7aSLisandro Dalcin -  mode - open mode of file
8455c6c1daeSBarry Smith $    FILE_MODE_WRITE - create new file for binary output
8465c6c1daeSBarry Smith $    FILE_MODE_READ - open existing file for binary input
8475c6c1daeSBarry Smith $    FILE_MODE_APPEND - open existing file for binary output
8485c6c1daeSBarry Smith 
8495c6c1daeSBarry Smith    Output Parameter:
850cc843e7aSLisandro Dalcin .  viewer - PetscViewer for binary input/output to use with the specified file
8515c6c1daeSBarry Smith 
8525c6c1daeSBarry Smith     Options Database Keys:
853811af0c4SBarry Smith +    -viewer_binary_filename <name> - name of file to use
854811af0c4SBarry Smith .    -viewer_binary_skip_info - true to skip opening an info file
855811af0c4SBarry Smith .    -viewer_binary_skip_options - true to not use options database while creating viewer
856811af0c4SBarry Smith .    -viewer_binary_skip_header - true to skip output object headers to the file
857811af0c4SBarry Smith -    -viewer_binary_mpiio - true to use MPI-IO for input and output to the file (more scalable for large problems)
8585c6c1daeSBarry Smith 
8595c6c1daeSBarry Smith    Level: beginner
8605c6c1daeSBarry Smith 
8615c6c1daeSBarry Smith    Note:
862811af0c4SBarry Smith    This `PetscViewer` should be destroyed with `PetscViewerDestroy()`.
8635c6c1daeSBarry Smith 
8645c6c1daeSBarry Smith     For reading files, the filename may begin with ftp:// or http:// and/or
8655c6c1daeSBarry Smith     end with .gz; in this case file is brought over and uncompressed.
8665c6c1daeSBarry Smith 
8675c6c1daeSBarry Smith     For creating files, if the file name ends with .gz it is automatically
8685c6c1daeSBarry Smith     compressed when closed.
8695c6c1daeSBarry Smith 
870811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
871db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
872db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`, `PetscViewerBinarySetUseMPIIO()`,
873db781477SPatrick Sanan           `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
8745c6c1daeSBarry Smith @*/
875*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryOpen(MPI_Comm comm, const char name[], PetscFileMode mode, PetscViewer *viewer)
876*d71ae5a4SJacob Faibussowitsch {
8775c6c1daeSBarry Smith   PetscFunctionBegin;
8789566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(comm, viewer));
8799566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(*viewer, PETSCVIEWERBINARY));
8809566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(*viewer, mode));
8819566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetName(*viewer, name));
8829566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions(*viewer));
8835c6c1daeSBarry Smith   PetscFunctionReturn(0);
8845c6c1daeSBarry Smith }
8855c6c1daeSBarry Smith 
8865c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
887*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryWriteReadMPIIO(PetscViewer viewer, void *data, PetscInt num, PetscInt *count, PetscDataType dtype, PetscBool write)
888*d71ae5a4SJacob Faibussowitsch {
88922a8f86cSLisandro Dalcin   MPI_Comm            comm    = PetscObjectComm((PetscObject)viewer);
8905c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
89122a8f86cSLisandro Dalcin   MPI_File            mfdes   = vbinary->mfdes;
8925c6c1daeSBarry Smith   MPI_Datatype        mdtype;
89369e10bbaSLisandro Dalcin   PetscMPIInt         rank, cnt;
8945c6c1daeSBarry Smith   MPI_Status          status;
8955c6c1daeSBarry Smith   MPI_Aint            ul, dsize;
8965c6c1daeSBarry Smith 
8975c6c1daeSBarry Smith   PetscFunctionBegin;
8989566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
8999566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(num, &cnt));
9009566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToMPIDataType(dtype, &mdtype));
9015c6c1daeSBarry Smith   if (write) {
90248a46eb9SPierre Jolivet     if (rank == 0) PetscCall(MPIU_File_write_at(mfdes, vbinary->moff, data, cnt, mdtype, &status));
9035c6c1daeSBarry Smith   } else {
904dd400576SPatrick Sanan     if (rank == 0) {
9059566063dSJacob Faibussowitsch       PetscCall(MPIU_File_read_at(mfdes, vbinary->moff, data, cnt, mdtype, &status));
9069566063dSJacob Faibussowitsch       if (cnt > 0) PetscCallMPI(MPI_Get_count(&status, mdtype, &cnt));
9075c6c1daeSBarry Smith     }
9089566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&cnt, 1, MPI_INT, 0, comm));
9099566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(data, cnt, mdtype, 0, comm));
91069e10bbaSLisandro Dalcin   }
9119566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_get_extent(mdtype, &ul, &dsize));
9125c6c1daeSBarry Smith   vbinary->moff += dsize * cnt;
9139860990eSLisandro Dalcin   if (count) *count = cnt;
9145c6c1daeSBarry Smith   PetscFunctionReturn(0);
9155c6c1daeSBarry Smith }
9165c6c1daeSBarry Smith #endif
9175c6c1daeSBarry Smith 
9185c6c1daeSBarry Smith /*@C
9195c6c1daeSBarry Smith    PetscViewerBinaryRead - Reads from a binary file, all processors get the same result
9205c6c1daeSBarry Smith 
921d083f849SBarry Smith    Collective
9225c6c1daeSBarry Smith 
9235c6c1daeSBarry Smith    Input Parameters:
9245c6c1daeSBarry Smith +  viewer - the binary viewer
925907376e6SBarry Smith .  data - location of the data to be written
926060da220SMatthew G. Knepley .  num - number of items of data to read
927907376e6SBarry Smith -  dtype - type of data to read
9285c6c1daeSBarry Smith 
929f8e4bde8SMatthew G. Knepley    Output Parameters:
9305972f5f3SLisandro Dalcin .  count - number of items of data actually read, or NULL.
931f8e4bde8SMatthew G. Knepley 
9325c6c1daeSBarry Smith    Level: beginner
9335c6c1daeSBarry Smith 
934811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
935db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
936db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
9375c6c1daeSBarry Smith @*/
938*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryRead(PetscViewer viewer, void *data, PetscInt num, PetscInt *count, PetscDataType dtype)
939*d71ae5a4SJacob Faibussowitsch {
94022a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
9415c6c1daeSBarry Smith 
94222a8f86cSLisandro Dalcin   PetscFunctionBegin;
94322a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
94422a8f86cSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, num, 3);
9459566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
94622a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
9475c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
948bc196f7cSDave May   if (vbinary->usempiio) {
9499566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWriteReadMPIIO(viewer, data, num, count, dtype, PETSC_FALSE));
9505c6c1daeSBarry Smith   } else {
9515c6c1daeSBarry Smith #endif
9529566063dSJacob Faibussowitsch     PetscCall(PetscBinarySynchronizedRead(PetscObjectComm((PetscObject)viewer), vbinary->fdes, data, num, count, dtype));
9535c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
9545c6c1daeSBarry Smith   }
9555c6c1daeSBarry Smith #endif
9565c6c1daeSBarry Smith   PetscFunctionReturn(0);
9575c6c1daeSBarry Smith }
9585c6c1daeSBarry Smith 
9595c6c1daeSBarry Smith /*@C
960811af0c4SBarry Smith    PetscViewerBinaryWrite - writes to a binary file, only from the first MPI rank
9615c6c1daeSBarry Smith 
962d083f849SBarry Smith    Collective
9635c6c1daeSBarry Smith 
9645c6c1daeSBarry Smith    Input Parameters:
9655c6c1daeSBarry Smith +  viewer - the binary viewer
9665c6c1daeSBarry Smith .  data - location of data
9675824c9d0SPatrick Sanan .  count - number of items of data to write
968f253e43cSLisandro Dalcin -  dtype - type of data to write
9695c6c1daeSBarry Smith 
9705c6c1daeSBarry Smith    Level: beginner
9715c6c1daeSBarry Smith 
972811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
973db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`, `PetscDataType`
974db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
9755c6c1daeSBarry Smith @*/
976*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryWrite(PetscViewer viewer, const void *data, PetscInt count, PetscDataType dtype)
977*d71ae5a4SJacob Faibussowitsch {
97822a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
9795c6c1daeSBarry Smith 
9805c6c1daeSBarry Smith   PetscFunctionBegin;
98122a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
98222a8f86cSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, count, 3);
9839566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
98422a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
9855c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
986bc196f7cSDave May   if (vbinary->usempiio) {
9879566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWriteReadMPIIO(viewer, (void *)data, count, NULL, dtype, PETSC_TRUE));
9885c6c1daeSBarry Smith   } else {
9895c6c1daeSBarry Smith #endif
9909566063dSJacob Faibussowitsch     PetscCall(PetscBinarySynchronizedWrite(PetscObjectComm((PetscObject)viewer), vbinary->fdes, data, count, dtype));
9915c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
9925c6c1daeSBarry Smith   }
9935c6c1daeSBarry Smith #endif
9945c6c1daeSBarry Smith   PetscFunctionReturn(0);
9955c6c1daeSBarry Smith }
9965c6c1daeSBarry Smith 
997*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryWriteReadAll(PetscViewer viewer, PetscBool write, void *data, PetscInt count, PetscInt start, PetscInt total, PetscDataType dtype)
998*d71ae5a4SJacob Faibussowitsch {
9995972f5f3SLisandro Dalcin   MPI_Comm              comm = PetscObjectComm((PetscObject)viewer);
10005972f5f3SLisandro Dalcin   PetscMPIInt           size, rank;
10015972f5f3SLisandro Dalcin   MPI_Datatype          mdtype;
1002ec4bef21SJose E. Roman   PETSC_UNUSED MPI_Aint lb;
1003ec4bef21SJose E. Roman   MPI_Aint              dsize;
10045972f5f3SLisandro Dalcin   PetscBool             useMPIIO;
10055972f5f3SLisandro Dalcin 
10065972f5f3SLisandro Dalcin   PetscFunctionBegin;
10075972f5f3SLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
1008064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveBool(viewer, ((start >= 0) || (start == PETSC_DETERMINE)), 5);
1009064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveBool(viewer, ((total >= 0) || (total == PETSC_DETERMINE)), 6);
1010064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveInt(viewer, total, 6);
10119566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
10125972f5f3SLisandro Dalcin 
10139566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToMPIDataType(dtype, &mdtype));
10149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_get_extent(mdtype, &lb, &dsize));
10159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
10169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
10175972f5f3SLisandro Dalcin 
10189566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetUseMPIIO(viewer, &useMPIIO));
10195972f5f3SLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
10205972f5f3SLisandro Dalcin   if (useMPIIO) {
10215972f5f3SLisandro Dalcin     MPI_File    mfdes;
10225972f5f3SLisandro Dalcin     MPI_Offset  off;
10235972f5f3SLisandro Dalcin     PetscMPIInt cnt;
10245972f5f3SLisandro Dalcin 
10255972f5f3SLisandro Dalcin     if (start == PETSC_DETERMINE) {
10269566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Scan(&count, &start, 1, MPIU_INT, MPI_SUM, comm));
10275972f5f3SLisandro Dalcin       start -= count;
10285972f5f3SLisandro Dalcin     }
10295972f5f3SLisandro Dalcin     if (total == PETSC_DETERMINE) {
10305972f5f3SLisandro Dalcin       total = start + count;
10319566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Bcast(&total, 1, MPIU_INT, size - 1, comm));
10325972f5f3SLisandro Dalcin     }
10339566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(count, &cnt));
10349566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryGetMPIIODescriptor(viewer, &mfdes));
10359566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryGetMPIIOOffset(viewer, &off));
10365972f5f3SLisandro Dalcin     off += (MPI_Offset)(start * dsize);
10375972f5f3SLisandro Dalcin     if (write) {
10389566063dSJacob Faibussowitsch       PetscCall(MPIU_File_write_at_all(mfdes, off, data, cnt, mdtype, MPI_STATUS_IGNORE));
10395972f5f3SLisandro Dalcin     } else {
10409566063dSJacob Faibussowitsch       PetscCall(MPIU_File_read_at_all(mfdes, off, data, cnt, mdtype, MPI_STATUS_IGNORE));
10415972f5f3SLisandro Dalcin     }
10425972f5f3SLisandro Dalcin     off = (MPI_Offset)(total * dsize);
10439566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryAddMPIIOOffset(viewer, off));
10445972f5f3SLisandro Dalcin     PetscFunctionReturn(0);
10455972f5f3SLisandro Dalcin   }
10465972f5f3SLisandro Dalcin #endif
10475972f5f3SLisandro Dalcin   {
10485972f5f3SLisandro Dalcin     int         fdes;
10495972f5f3SLisandro Dalcin     char       *workbuf = NULL;
1050dd400576SPatrick Sanan     PetscInt    tcount = rank == 0 ? 0 : count, maxcount = 0, message_count, flowcontrolcount;
10515972f5f3SLisandro Dalcin     PetscMPIInt tag, cnt, maxcnt, scnt = 0, rcnt = 0, j;
10525972f5f3SLisandro Dalcin     MPI_Status  status;
10535972f5f3SLisandro Dalcin 
10549566063dSJacob Faibussowitsch     PetscCall(PetscCommGetNewTag(comm, &tag));
10559566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Reduce(&tcount, &maxcount, 1, MPIU_INT, MPI_MAX, 0, comm));
10569566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(maxcount, &maxcnt));
10579566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(count, &cnt));
10585972f5f3SLisandro Dalcin 
10599566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryGetDescriptor(viewer, &fdes));
10609566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlowControlStart(viewer, &message_count, &flowcontrolcount));
1061dd400576SPatrick Sanan     if (rank == 0) {
10629566063dSJacob Faibussowitsch       PetscCall(PetscMalloc(maxcnt * dsize, &workbuf));
10635972f5f3SLisandro Dalcin       if (write) {
10649566063dSJacob Faibussowitsch         PetscCall(PetscBinaryWrite(fdes, data, cnt, dtype));
10655972f5f3SLisandro Dalcin       } else {
10669566063dSJacob Faibussowitsch         PetscCall(PetscBinaryRead(fdes, data, cnt, NULL, dtype));
10675972f5f3SLisandro Dalcin       }
10685972f5f3SLisandro Dalcin       for (j = 1; j < size; j++) {
10699566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlowControlStepMain(viewer, j, &message_count, flowcontrolcount));
10705972f5f3SLisandro Dalcin         if (write) {
10719566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Recv(workbuf, maxcnt, mdtype, j, tag, comm, &status));
10729566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Get_count(&status, mdtype, &rcnt));
10739566063dSJacob Faibussowitsch           PetscCall(PetscBinaryWrite(fdes, workbuf, rcnt, dtype));
10745972f5f3SLisandro Dalcin         } else {
10759566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Recv(&scnt, 1, MPI_INT, j, tag, comm, MPI_STATUS_IGNORE));
10769566063dSJacob Faibussowitsch           PetscCall(PetscBinaryRead(fdes, workbuf, scnt, NULL, dtype));
10779566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Send(workbuf, scnt, mdtype, j, tag, comm));
10785972f5f3SLisandro Dalcin         }
10795972f5f3SLisandro Dalcin       }
10809566063dSJacob Faibussowitsch       PetscCall(PetscFree(workbuf));
10819566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlowControlEndMain(viewer, &message_count));
10825972f5f3SLisandro Dalcin     } else {
10839566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlowControlStepWorker(viewer, rank, &message_count));
10845972f5f3SLisandro Dalcin       if (write) {
10859566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(data, cnt, mdtype, 0, tag, comm));
10865972f5f3SLisandro Dalcin       } else {
10879566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(&cnt, 1, MPI_INT, 0, tag, comm));
10889566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Recv(data, cnt, mdtype, 0, tag, comm, MPI_STATUS_IGNORE));
10895972f5f3SLisandro Dalcin       }
10909566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlowControlEndWorker(viewer, &message_count));
10915972f5f3SLisandro Dalcin     }
10925972f5f3SLisandro Dalcin   }
10935972f5f3SLisandro Dalcin   PetscFunctionReturn(0);
10945972f5f3SLisandro Dalcin }
10955972f5f3SLisandro Dalcin 
10965972f5f3SLisandro Dalcin /*@C
1097811af0c4SBarry Smith    PetscViewerBinaryReadAll - reads from a binary file from all MPI ranks, each rank receives its own portion of the data
10985972f5f3SLisandro Dalcin 
10995972f5f3SLisandro Dalcin    Collective
11005972f5f3SLisandro Dalcin 
11015972f5f3SLisandro Dalcin    Input Parameters:
11025972f5f3SLisandro Dalcin +  viewer - the binary viewer
11035972f5f3SLisandro Dalcin .  data - location of data
11045972f5f3SLisandro Dalcin .  count - local number of items of data to read
1105811af0c4SBarry Smith .  start - local start, can be `PETSC_DETERMINE`
1106811af0c4SBarry Smith .  total - global number of items of data to read, can be `PETSC_DETERMINE`
11075972f5f3SLisandro Dalcin -  dtype - type of data to read
11085972f5f3SLisandro Dalcin 
11095972f5f3SLisandro Dalcin    Level: advanced
11105972f5f3SLisandro Dalcin 
1111811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryRead()`, `PetscViewerBinaryWriteAll()`
11125972f5f3SLisandro Dalcin @*/
1113*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryReadAll(PetscViewer viewer, void *data, PetscInt count, PetscInt start, PetscInt total, PetscDataType dtype)
1114*d71ae5a4SJacob Faibussowitsch {
11155972f5f3SLisandro Dalcin   PetscFunctionBegin;
11169566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWriteReadAll(viewer, PETSC_FALSE, data, count, start, total, dtype));
11175972f5f3SLisandro Dalcin   PetscFunctionReturn(0);
11185972f5f3SLisandro Dalcin }
11195972f5f3SLisandro Dalcin 
11205972f5f3SLisandro Dalcin /*@C
1121811af0c4SBarry Smith    PetscViewerBinaryWriteAll - writes to a binary file from all MPI ranks, each rank writes its own portion of the data
11225972f5f3SLisandro Dalcin 
11235972f5f3SLisandro Dalcin    Collective
11245972f5f3SLisandro Dalcin 
11255972f5f3SLisandro Dalcin    Input Parameters:
11265972f5f3SLisandro Dalcin +  viewer - the binary viewer
11275972f5f3SLisandro Dalcin .  data - location of data
11285972f5f3SLisandro Dalcin .  count - local number of items of data to write
1129811af0c4SBarry Smith .  start - local start, can be `PETSC_DETERMINE`
1130811af0c4SBarry Smith .  total - global number of items of data to write, can be `PETSC_DETERMINE`
11315972f5f3SLisandro Dalcin -  dtype - type of data to write
11325972f5f3SLisandro Dalcin 
11335972f5f3SLisandro Dalcin    Level: advanced
11345972f5f3SLisandro Dalcin 
1135811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryWriteAll()`, `PetscViewerBinaryReadAll()`
11365972f5f3SLisandro Dalcin @*/
1137*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryWriteAll(PetscViewer viewer, const void *data, PetscInt count, PetscInt start, PetscInt total, PetscDataType dtype)
1138*d71ae5a4SJacob Faibussowitsch {
11395972f5f3SLisandro Dalcin   PetscFunctionBegin;
11409566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWriteReadAll(viewer, PETSC_TRUE, (void *)data, count, start, total, dtype));
11415972f5f3SLisandro Dalcin   PetscFunctionReturn(0);
11425972f5f3SLisandro Dalcin }
11435972f5f3SLisandro Dalcin 
11445c6c1daeSBarry Smith /*@C
1145811af0c4SBarry Smith    PetscViewerBinaryWriteStringArray - writes to a binary file, only from the first MPI rank, an array of strings
11465c6c1daeSBarry Smith 
1147d083f849SBarry Smith    Collective
11485c6c1daeSBarry Smith 
11495c6c1daeSBarry Smith    Input Parameters:
11505c6c1daeSBarry Smith +  viewer - the binary viewer
11515c6c1daeSBarry Smith -  data - location of the array of strings
11525c6c1daeSBarry Smith 
11535c6c1daeSBarry Smith    Level: intermediate
11545c6c1daeSBarry Smith 
1155811af0c4SBarry Smith     Note:
1156811af0c4SBarry Smith     The array of strings must be null terminated
11575c6c1daeSBarry Smith 
1158811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
1159db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
1160db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
11615c6c1daeSBarry Smith @*/
1162*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryWriteStringArray(PetscViewer viewer, const char *const *data)
1163*d71ae5a4SJacob Faibussowitsch {
11645c6c1daeSBarry Smith   PetscInt i, n = 0, *sizes;
1165cc843e7aSLisandro Dalcin   size_t   len;
11665c6c1daeSBarry Smith 
116722a8f86cSLisandro Dalcin   PetscFunctionBegin;
11689566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
11695c6c1daeSBarry Smith   /* count number of strings */
11709371c9d4SSatish Balay   while (data[n++])
11719371c9d4SSatish Balay     ;
11725c6c1daeSBarry Smith   n--;
11739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n + 1, &sizes));
11745c6c1daeSBarry Smith   sizes[0] = n;
11755c6c1daeSBarry Smith   for (i = 0; i < n; i++) {
11769566063dSJacob Faibussowitsch     PetscCall(PetscStrlen(data[i], &len));
1177cc843e7aSLisandro Dalcin     sizes[i + 1] = (PetscInt)len + 1; /* size includes space for the null terminator */
11785c6c1daeSBarry Smith   }
11799566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWrite(viewer, sizes, n + 1, PETSC_INT));
118048a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall(PetscViewerBinaryWrite(viewer, (void *)data[i], sizes[i + 1], PETSC_CHAR));
11819566063dSJacob Faibussowitsch   PetscCall(PetscFree(sizes));
11825c6c1daeSBarry Smith   PetscFunctionReturn(0);
11835c6c1daeSBarry Smith }
11845c6c1daeSBarry Smith 
11855c6c1daeSBarry Smith /*@C
1186811af0c4SBarry Smith    PetscViewerBinaryReadStringArray - reads a binary file an array of strings to all MPI ranks
11875c6c1daeSBarry Smith 
1188d083f849SBarry Smith    Collective
11895c6c1daeSBarry Smith 
11905c6c1daeSBarry Smith    Input Parameter:
11915c6c1daeSBarry Smith .  viewer - the binary viewer
11925c6c1daeSBarry Smith 
11935c6c1daeSBarry Smith    Output Parameter:
11945c6c1daeSBarry Smith .  data - location of the array of strings
11955c6c1daeSBarry Smith 
11965c6c1daeSBarry Smith    Level: intermediate
11975c6c1daeSBarry Smith 
1198811af0c4SBarry Smith     Note:
1199811af0c4SBarry Smith     The array of strings must null terminated
12005c6c1daeSBarry Smith 
1201811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
1202db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
1203db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
12045c6c1daeSBarry Smith @*/
1205*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryReadStringArray(PetscViewer viewer, char ***data)
1206*d71ae5a4SJacob Faibussowitsch {
1207060da220SMatthew G. Knepley   PetscInt i, n, *sizes, N = 0;
12085c6c1daeSBarry Smith 
120922a8f86cSLisandro Dalcin   PetscFunctionBegin;
12109566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
12115c6c1daeSBarry Smith   /* count number of strings */
12129566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &n, 1, NULL, PETSC_INT));
12139566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &sizes));
12149566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, sizes, n, NULL, PETSC_INT));
1215a297a907SKarl Rupp   for (i = 0; i < n; i++) N += sizes[i];
12169566063dSJacob Faibussowitsch   PetscCall(PetscMalloc((n + 1) * sizeof(char *) + N * sizeof(char), data));
12175c6c1daeSBarry Smith   (*data)[0] = (char *)((*data) + n + 1);
1218a297a907SKarl Rupp   for (i = 1; i < n; i++) (*data)[i] = (*data)[i - 1] + sizes[i - 1];
12199566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, (*data)[0], N, NULL, PETSC_CHAR));
122002c9f0b5SLisandro Dalcin   (*data)[n] = NULL;
12219566063dSJacob Faibussowitsch   PetscCall(PetscFree(sizes));
12225c6c1daeSBarry Smith   PetscFunctionReturn(0);
12235c6c1daeSBarry Smith }
12245c6c1daeSBarry Smith 
1225cc843e7aSLisandro Dalcin /*@C
1226cc843e7aSLisandro Dalcin      PetscViewerFileSetMode - Sets the open mode of file
1227cc843e7aSLisandro Dalcin 
1228811af0c4SBarry Smith     Logically Collective on viewer
1229cc843e7aSLisandro Dalcin 
1230cc843e7aSLisandro Dalcin   Input Parameters:
1231811af0c4SBarry Smith +  viewer - the `PetscViewer`; must be a a `PETSCVIEWERBINARY`, `PETSCVIEWERMATLAB`, `PETSCVIEWERHDF5`, or `PETSCVIEWERASCII`  `PetscViewer`
1232cc843e7aSLisandro Dalcin -  mode - open mode of file
1233cc843e7aSLisandro Dalcin $    FILE_MODE_WRITE - create new file for output
1234cc843e7aSLisandro Dalcin $    FILE_MODE_READ - open existing file for input
1235cc843e7aSLisandro Dalcin $    FILE_MODE_APPEND - open existing file for output
1236cc843e7aSLisandro Dalcin 
1237cc843e7aSLisandro Dalcin   Level: advanced
1238cc843e7aSLisandro Dalcin 
1239db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
1240cc843e7aSLisandro Dalcin @*/
1241*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerFileSetMode(PetscViewer viewer, PetscFileMode mode)
1242*d71ae5a4SJacob Faibussowitsch {
1243cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1244cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1245cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveEnum(viewer, mode, 2);
124608401ef6SPierre Jolivet   PetscCheck(mode != FILE_MODE_UNDEFINED, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Cannot set FILE_MODE_UNDEFINED");
1247f7d195e4SLawrence Mitchell   PetscCheck(mode >= FILE_MODE_UNDEFINED && mode <= FILE_MODE_APPEND_UPDATE, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_OUTOFRANGE, "Invalid file mode %d", (int)mode);
1248cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerFileSetMode_C", (PetscViewer, PetscFileMode), (viewer, mode));
1249cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1250cc843e7aSLisandro Dalcin }
1251cc843e7aSLisandro Dalcin 
1252*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetMode_Binary(PetscViewer viewer, PetscFileMode mode)
1253*d71ae5a4SJacob Faibussowitsch {
1254cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1255cc843e7aSLisandro Dalcin 
1256cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1257cc73adaaSBarry Smith   PetscCheck(!viewer->setupcalled || vbinary->filemode == mode, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Cannot change mode to %s after setup", PetscFileModes[mode]);
1258cc843e7aSLisandro Dalcin   vbinary->filemode = mode;
1259cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1260cc843e7aSLisandro Dalcin }
1261cc843e7aSLisandro Dalcin 
1262cc843e7aSLisandro Dalcin /*@C
1263cc843e7aSLisandro Dalcin      PetscViewerFileGetMode - Gets the open mode of file
1264cc843e7aSLisandro Dalcin 
1265cc843e7aSLisandro Dalcin     Not Collective
1266cc843e7aSLisandro Dalcin 
1267cc843e7aSLisandro Dalcin   Input Parameter:
1268811af0c4SBarry Smith .  viewer - the `PetscViewer`; must be a `PETSCVIEWERBINARY`, `PETSCVIEWERMATLAB`, `PETSCVIEWERHDF5`, or `PETSCVIEWERASCII`  `PetscViewer`
1269cc843e7aSLisandro Dalcin 
1270cc843e7aSLisandro Dalcin   Output Parameter:
1271cc843e7aSLisandro Dalcin .  mode - open mode of file
1272cc843e7aSLisandro Dalcin $    FILE_MODE_WRITE - create new file for binary output
1273cc843e7aSLisandro Dalcin $    FILE_MODE_READ - open existing file for binary input
1274cc843e7aSLisandro Dalcin $    FILE_MODE_APPEND - open existing file for binary output
1275cc843e7aSLisandro Dalcin 
1276cc843e7aSLisandro Dalcin   Level: advanced
1277cc843e7aSLisandro Dalcin 
1278db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
1279cc843e7aSLisandro Dalcin @*/
1280*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerFileGetMode(PetscViewer viewer, PetscFileMode *mode)
1281*d71ae5a4SJacob Faibussowitsch {
1282cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1283cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1284cc843e7aSLisandro Dalcin   PetscValidPointer(mode, 2);
1285cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerFileGetMode_C", (PetscViewer, PetscFileMode *), (viewer, mode));
1286cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1287cc843e7aSLisandro Dalcin }
1288cc843e7aSLisandro Dalcin 
1289*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetMode_Binary(PetscViewer viewer, PetscFileMode *mode)
1290*d71ae5a4SJacob Faibussowitsch {
1291cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1292cc843e7aSLisandro Dalcin 
1293cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1294cc843e7aSLisandro Dalcin   *mode = vbinary->filemode;
1295cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1296cc843e7aSLisandro Dalcin }
1297cc843e7aSLisandro Dalcin 
1298*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetName_Binary(PetscViewer viewer, const char name[])
1299*d71ae5a4SJacob Faibussowitsch {
1300cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1301cc843e7aSLisandro Dalcin 
1302cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1303cc843e7aSLisandro Dalcin   if (viewer->setupcalled && vbinary->filename) {
1304cc843e7aSLisandro Dalcin     /* gzip can be run after the file with the previous filename has been closed */
13059566063dSJacob Faibussowitsch     PetscCall(PetscFree(vbinary->ogzfilename));
13069566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(vbinary->filename, &vbinary->ogzfilename));
1307cc843e7aSLisandro Dalcin   }
13089566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary->filename));
13099566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &vbinary->filename));
1310cc843e7aSLisandro Dalcin   viewer->setupcalled = PETSC_FALSE;
1311cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1312cc843e7aSLisandro Dalcin }
1313cc843e7aSLisandro Dalcin 
1314*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetName_Binary(PetscViewer viewer, const char **name)
1315*d71ae5a4SJacob Faibussowitsch {
13165c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
13175c6c1daeSBarry Smith 
13185c6c1daeSBarry Smith   PetscFunctionBegin;
13195c6c1daeSBarry Smith   *name = vbinary->filename;
13205c6c1daeSBarry Smith   PetscFunctionReturn(0);
13215c6c1daeSBarry Smith }
13225c6c1daeSBarry Smith 
13235c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
1324*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetUp_BinaryMPIIO(PetscViewer viewer)
1325*d71ae5a4SJacob Faibussowitsch {
13265c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1327e8a65b7dSLisandro Dalcin   int                 amode;
13285c6c1daeSBarry Smith 
13295c6c1daeSBarry Smith   PetscFunctionBegin;
1330a297a907SKarl Rupp   vbinary->storecompressed = PETSC_FALSE;
13315c6c1daeSBarry Smith 
1332cc843e7aSLisandro Dalcin   vbinary->moff = 0;
1333cc843e7aSLisandro Dalcin   switch (vbinary->filemode) {
1334*d71ae5a4SJacob Faibussowitsch   case FILE_MODE_READ:
1335*d71ae5a4SJacob Faibussowitsch     amode = MPI_MODE_RDONLY;
1336*d71ae5a4SJacob Faibussowitsch     break;
1337*d71ae5a4SJacob Faibussowitsch   case FILE_MODE_WRITE:
1338*d71ae5a4SJacob Faibussowitsch     amode = MPI_MODE_WRONLY | MPI_MODE_CREATE;
1339*d71ae5a4SJacob Faibussowitsch     break;
1340*d71ae5a4SJacob Faibussowitsch   case FILE_MODE_APPEND:
1341*d71ae5a4SJacob Faibussowitsch     amode = MPI_MODE_WRONLY | MPI_MODE_CREATE | MPI_MODE_APPEND;
1342*d71ae5a4SJacob Faibussowitsch     break;
1343*d71ae5a4SJacob Faibussowitsch   case FILE_MODE_UNDEFINED:
1344*d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerSetUp()");
1345*d71ae5a4SJacob Faibussowitsch   default:
1346*d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[vbinary->filemode]);
13475c6c1daeSBarry Smith   }
13489566063dSJacob Faibussowitsch   PetscCallMPI(MPI_File_open(PetscObjectComm((PetscObject)viewer), vbinary->filename, amode, MPI_INFO_NULL, &vbinary->mfdes));
134922a8f86cSLisandro Dalcin   /*
135022a8f86cSLisandro Dalcin       The MPI standard does not have MPI_MODE_TRUNCATE. We emulate this behavior by setting the file size to zero.
135122a8f86cSLisandro Dalcin   */
13529566063dSJacob Faibussowitsch   if (vbinary->filemode == FILE_MODE_WRITE) PetscCallMPI(MPI_File_set_size(vbinary->mfdes, 0));
135322a8f86cSLisandro Dalcin   /*
135422a8f86cSLisandro Dalcin       Initially, all processes view the file as a linear byte stream. Therefore, for files opened with MPI_MODE_APPEND,
135522a8f86cSLisandro Dalcin       MPI_File_get_position[_shared](fh, &offset) returns the absolute byte position at the end of file.
135622a8f86cSLisandro Dalcin       Otherwise, we would need to call MPI_File_get_byte_offset(fh, offset, &byte_offset) to convert
135722a8f86cSLisandro Dalcin       the offset in etype units to an absolute byte position.
135822a8f86cSLisandro Dalcin    */
13599566063dSJacob Faibussowitsch   if (vbinary->filemode == FILE_MODE_APPEND) PetscCallMPI(MPI_File_get_position(vbinary->mfdes, &vbinary->moff));
1360cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1361cc843e7aSLisandro Dalcin }
1362cc843e7aSLisandro Dalcin #endif
13635c6c1daeSBarry Smith 
1364*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetUp_BinarySTDIO(PetscViewer viewer)
1365*d71ae5a4SJacob Faibussowitsch {
1366cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1367cc843e7aSLisandro Dalcin   const char         *fname;
1368cc843e7aSLisandro Dalcin   char                bname[PETSC_MAX_PATH_LEN], *gz;
1369cc843e7aSLisandro Dalcin   PetscBool           found;
1370cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
13715c6c1daeSBarry Smith 
1372cc843e7aSLisandro Dalcin   PetscFunctionBegin;
13739566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
13745c6c1daeSBarry Smith 
1375cc843e7aSLisandro Dalcin   /* if file name ends in .gz strip that off and note user wants file compressed */
1376cc843e7aSLisandro Dalcin   vbinary->storecompressed = PETSC_FALSE;
1377cc843e7aSLisandro Dalcin   if (vbinary->filemode == FILE_MODE_WRITE) {
13789566063dSJacob Faibussowitsch     PetscCall(PetscStrstr(vbinary->filename, ".gz", &gz));
13799371c9d4SSatish Balay     if (gz && gz[3] == 0) {
13809371c9d4SSatish Balay       *gz                      = 0;
13819371c9d4SSatish Balay       vbinary->storecompressed = PETSC_TRUE;
13829371c9d4SSatish Balay     }
1383cc843e7aSLisandro Dalcin   }
1384cc843e7aSLisandro Dalcin #if !defined(PETSC_HAVE_POPEN)
138528b400f6SJacob Faibussowitsch   PetscCheck(!vbinary->storecompressed, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP_SYS, "Cannot run gzip on this machine");
1386cc843e7aSLisandro Dalcin #endif
1387cc843e7aSLisandro Dalcin 
1388cc843e7aSLisandro Dalcin   fname = vbinary->filename;
1389cc843e7aSLisandro Dalcin   if (vbinary->filemode == FILE_MODE_READ) { /* possibly get the file from remote site or compressed file */
13909566063dSJacob Faibussowitsch     PetscCall(PetscFileRetrieve(PetscObjectComm((PetscObject)viewer), fname, bname, PETSC_MAX_PATH_LEN, &found));
139128b400f6SJacob Faibussowitsch     PetscCheck(found, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_OPEN, "Cannot locate file: %s", fname);
1392cc843e7aSLisandro Dalcin     fname = bname;
13935c6c1daeSBarry Smith   }
13945c6c1daeSBarry Smith 
1395cc843e7aSLisandro Dalcin   vbinary->fdes = -1;
1396dd400576SPatrick Sanan   if (rank == 0) { /* only first processor opens file*/
1397cc843e7aSLisandro Dalcin     PetscFileMode mode = vbinary->filemode;
1398cc843e7aSLisandro Dalcin     if (mode == FILE_MODE_APPEND) {
1399cc843e7aSLisandro Dalcin       /* check if asked to append to a non-existing file */
14009566063dSJacob Faibussowitsch       PetscCall(PetscTestFile(fname, '\0', &found));
1401cc843e7aSLisandro Dalcin       if (!found) mode = FILE_MODE_WRITE;
1402cc843e7aSLisandro Dalcin     }
14039566063dSJacob Faibussowitsch     PetscCall(PetscBinaryOpen(fname, mode, &vbinary->fdes));
1404cc843e7aSLisandro Dalcin   }
1405cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1406cc843e7aSLisandro Dalcin }
1407cc843e7aSLisandro Dalcin 
1408*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetUp_BinaryInfo(PetscViewer viewer)
1409*d71ae5a4SJacob Faibussowitsch {
1410cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1411cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
1412cc843e7aSLisandro Dalcin   PetscBool           found;
1413cc843e7aSLisandro Dalcin 
1414cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1415cc843e7aSLisandro Dalcin   vbinary->fdes_info = NULL;
14169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
1417dd400576SPatrick Sanan   if (!vbinary->skipinfo && (vbinary->filemode == FILE_MODE_READ || rank == 0)) {
1418cc843e7aSLisandro Dalcin     char infoname[PETSC_MAX_PATH_LEN], iname[PETSC_MAX_PATH_LEN], *gz;
1419cc843e7aSLisandro Dalcin 
14209566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(infoname, vbinary->filename, sizeof(infoname)));
1421cc843e7aSLisandro Dalcin     /* remove .gz if it ends file name */
14229566063dSJacob Faibussowitsch     PetscCall(PetscStrstr(infoname, ".gz", &gz));
1423cc843e7aSLisandro Dalcin     if (gz && gz[3] == 0) *gz = 0;
1424cc843e7aSLisandro Dalcin 
14259566063dSJacob Faibussowitsch     PetscCall(PetscStrlcat(infoname, ".info", sizeof(infoname)));
1426cc843e7aSLisandro Dalcin     if (vbinary->filemode == FILE_MODE_READ) {
14279566063dSJacob Faibussowitsch       PetscCall(PetscFixFilename(infoname, iname));
14289566063dSJacob Faibussowitsch       PetscCall(PetscFileRetrieve(PetscObjectComm((PetscObject)viewer), iname, infoname, PETSC_MAX_PATH_LEN, &found));
14299566063dSJacob Faibussowitsch       if (found) PetscCall(PetscOptionsInsertFile(PetscObjectComm((PetscObject)viewer), ((PetscObject)viewer)->options, infoname, PETSC_FALSE));
1430dd400576SPatrick Sanan     } else if (rank == 0) { /* write or append */
1431cc843e7aSLisandro Dalcin       const char *omode  = (vbinary->filemode == FILE_MODE_APPEND) ? "a" : "w";
1432cc843e7aSLisandro Dalcin       vbinary->fdes_info = fopen(infoname, omode);
143328b400f6SJacob Faibussowitsch       PetscCheck(vbinary->fdes_info, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open .info file %s for writing", infoname);
14345c6c1daeSBarry Smith     }
14355c6c1daeSBarry Smith   }
14365c6c1daeSBarry Smith   PetscFunctionReturn(0);
14375c6c1daeSBarry Smith }
14385c6c1daeSBarry Smith 
1439*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerSetUp_Binary(PetscViewer viewer)
1440*d71ae5a4SJacob Faibussowitsch {
14415c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1442cc843e7aSLisandro Dalcin   PetscBool           usempiio;
1443cc843e7aSLisandro Dalcin 
14445c6c1daeSBarry Smith   PetscFunctionBegin;
14459566063dSJacob Faibussowitsch   if (!vbinary->setfromoptionscalled) PetscCall(PetscViewerSetFromOptions(viewer));
144628b400f6SJacob Faibussowitsch   PetscCheck(vbinary->filename, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetName()");
144708401ef6SPierre Jolivet   PetscCheck(vbinary->filemode != (PetscFileMode)-1, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode()");
14489566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_Binary(viewer));
1449cc843e7aSLisandro Dalcin 
14509566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetUseMPIIO(viewer, &usempiio));
1451cc843e7aSLisandro Dalcin   if (usempiio) {
1452cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
14539566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetUp_BinaryMPIIO(viewer));
1454cc843e7aSLisandro Dalcin #endif
1455cc843e7aSLisandro Dalcin   } else {
14569566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetUp_BinarySTDIO(viewer));
1457cc843e7aSLisandro Dalcin   }
14589566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetUp_BinaryInfo(viewer));
1459cc843e7aSLisandro Dalcin 
14609566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectState((PetscObject)viewer, "File: %s", vbinary->filename));
14615c6c1daeSBarry Smith   PetscFunctionReturn(0);
14625c6c1daeSBarry Smith }
14635c6c1daeSBarry Smith 
1464*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerView_Binary(PetscViewer v, PetscViewer viewer)
1465*d71ae5a4SJacob Faibussowitsch {
1466cb6ad94fSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
1467cb6ad94fSLisandro Dalcin   const char         *fname   = vbinary->filename ? vbinary->filename : "not yet set";
1468cc843e7aSLisandro Dalcin   const char         *fmode   = vbinary->filemode != (PetscFileMode)-1 ? PetscFileModes[vbinary->filemode] : "not yet set";
1469cb6ad94fSLisandro Dalcin   PetscBool           usempiio;
14702bf49c77SBarry Smith 
14712bf49c77SBarry Smith   PetscFunctionBegin;
14729566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetUseMPIIO(v, &usempiio));
14739566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Filename: %s\n", fname));
14749566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Mode: %s (%s)\n", fmode, usempiio ? "mpiio" : "stdio"));
14752bf49c77SBarry Smith   PetscFunctionReturn(0);
14762bf49c77SBarry Smith }
14772bf49c77SBarry Smith 
1478*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerSetFromOptions_Binary(PetscViewer viewer, PetscOptionItems *PetscOptionsObject)
1479*d71ae5a4SJacob Faibussowitsch {
148022a8f86cSLisandro Dalcin   PetscViewer_Binary *binary = (PetscViewer_Binary *)viewer->data;
1481d22fd6bcSDave May   char                defaultname[PETSC_MAX_PATH_LEN];
148203a1d59fSDave May   PetscBool           flg;
1483e0385b85SDave May 
148403a1d59fSDave May   PetscFunctionBegin;
148522a8f86cSLisandro Dalcin   if (viewer->setupcalled) PetscFunctionReturn(0);
1486d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Binary PetscViewer Options");
14879566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(defaultname, PETSC_MAX_PATH_LEN - 1, "binaryoutput"));
14889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-viewer_binary_filename", "Specify filename", "PetscViewerFileSetName", defaultname, defaultname, sizeof(defaultname), &flg));
14899566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscViewerFileSetName_Binary(viewer, defaultname));
14909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_skip_info", "Skip writing/reading .info file", "PetscViewerBinarySetSkipInfo", binary->skipinfo, &binary->skipinfo, NULL));
14919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_skip_options", "Skip parsing Vec/Mat load options", "PetscViewerBinarySetSkipOptions", binary->skipoptions, &binary->skipoptions, NULL));
14929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_skip_header", "Skip writing/reading header information", "PetscViewerBinarySetSkipHeader", binary->skipheader, &binary->skipheader, NULL));
149303a1d59fSDave May #if defined(PETSC_HAVE_MPIIO)
14949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_mpiio", "Use MPI-IO functionality to write/read binary file", "PetscViewerBinarySetUseMPIIO", binary->usempiio, &binary->usempiio, NULL));
14955972f5f3SLisandro Dalcin #else
14969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_mpiio", "Use MPI-IO functionality to write/read binary file (NOT AVAILABLE)", "PetscViewerBinarySetUseMPIIO", PETSC_FALSE, NULL, NULL));
149703a1d59fSDave May #endif
1498d0609cedSBarry Smith   PetscOptionsHeadEnd();
1499bc196f7cSDave May   binary->setfromoptionscalled = PETSC_TRUE;
150003a1d59fSDave May   PetscFunctionReturn(0);
150103a1d59fSDave May }
150203a1d59fSDave May 
15038556b5ebSBarry Smith /*MC
15048556b5ebSBarry Smith    PETSCVIEWERBINARY - A viewer that saves to binary files
15058556b5ebSBarry Smith 
15061b266c99SBarry Smith   Level: beginner
15071b266c99SBarry Smith 
1508811af0c4SBarry Smith .seealso: `PetscViewerBinaryOpen()`, `PETSC_VIEWER_STDOUT_()`, `PETSC_VIEWER_STDOUT_SELF`, `PETSC_VIEWER_STDOUT_WORLD`, `PetscViewerCreate()`, `PetscViewerASCIIOpen()`,
1509811af0c4SBarry Smith           `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`, `PETSCVIEWERDRAW`, `PETSCVIEWERSOCKET`
1510811af0c4SBarry Smith           `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()`,
1511811af0c4SBarry Smith           `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`
15128556b5ebSBarry Smith M*/
15138556b5ebSBarry Smith 
1514*d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscViewerCreate_Binary(PetscViewer v)
1515*d71ae5a4SJacob Faibussowitsch {
15165c6c1daeSBarry Smith   PetscViewer_Binary *vbinary;
15175c6c1daeSBarry Smith 
15185c6c1daeSBarry Smith   PetscFunctionBegin;
15194dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&vbinary));
15205c6c1daeSBarry Smith   v->data = (void *)vbinary;
1521cc843e7aSLisandro Dalcin 
152203a1d59fSDave May   v->ops->setfromoptions   = PetscViewerSetFromOptions_Binary;
15235c6c1daeSBarry Smith   v->ops->destroy          = PetscViewerDestroy_Binary;
15242bf49c77SBarry Smith   v->ops->view             = PetscViewerView_Binary;
1525c98fd787SBarry Smith   v->ops->setup            = PetscViewerSetUp_Binary;
1526cc843e7aSLisandro Dalcin   v->ops->flush            = NULL; /* Should we support Flush() ? */
1527cc843e7aSLisandro Dalcin   v->ops->getsubviewer     = PetscViewerGetSubViewer_Binary;
1528cc843e7aSLisandro Dalcin   v->ops->restoresubviewer = PetscViewerRestoreSubViewer_Binary;
1529cc843e7aSLisandro Dalcin   v->ops->read             = PetscViewerBinaryRead;
1530cc843e7aSLisandro Dalcin 
1531cc843e7aSLisandro Dalcin   vbinary->fdes = -1;
1532e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
1533cc843e7aSLisandro Dalcin   vbinary->usempiio = PETSC_FALSE;
1534e8a65b7dSLisandro Dalcin   vbinary->mfdes    = MPI_FILE_NULL;
1535e8a65b7dSLisandro Dalcin   vbinary->mfsub    = MPI_FILE_NULL;
1536e8a65b7dSLisandro Dalcin #endif
1537cc843e7aSLisandro Dalcin   vbinary->filename        = NULL;
15387e4fd573SVaclav Hapla   vbinary->filemode        = FILE_MODE_UNDEFINED;
153902c9f0b5SLisandro Dalcin   vbinary->fdes_info       = NULL;
15405c6c1daeSBarry Smith   vbinary->skipinfo        = PETSC_FALSE;
15415c6c1daeSBarry Smith   vbinary->skipoptions     = PETSC_TRUE;
15425c6c1daeSBarry Smith   vbinary->skipheader      = PETSC_FALSE;
15435c6c1daeSBarry Smith   vbinary->storecompressed = PETSC_FALSE;
1544f90597f1SStefano Zampini   vbinary->ogzfilename     = NULL;
15455c6c1daeSBarry Smith   vbinary->flowcontrol     = 256; /* seems a good number for Cray XT-5 */
15465c6c1daeSBarry Smith 
1547cc843e7aSLisandro Dalcin   vbinary->setfromoptionscalled = PETSC_FALSE;
1548cc843e7aSLisandro Dalcin 
15499566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetFlowControl_C", PetscViewerBinaryGetFlowControl_Binary));
15509566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetFlowControl_C", PetscViewerBinarySetFlowControl_Binary));
15519566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipHeader_C", PetscViewerBinaryGetSkipHeader_Binary));
15529566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipHeader_C", PetscViewerBinarySetSkipHeader_Binary));
15539566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipOptions_C", PetscViewerBinaryGetSkipOptions_Binary));
15549566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipOptions_C", PetscViewerBinarySetSkipOptions_Binary));
15559566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipInfo_C", PetscViewerBinaryGetSkipInfo_Binary));
15569566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipInfo_C", PetscViewerBinarySetSkipInfo_Binary));
15579566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetInfoPointer_C", PetscViewerBinaryGetInfoPointer_Binary));
15589566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_Binary));
15599566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_Binary));
15609566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_Binary));
15619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_Binary));
15625c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
15639566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetUseMPIIO_C", PetscViewerBinaryGetUseMPIIO_Binary));
15649566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetUseMPIIO_C", PetscViewerBinarySetUseMPIIO_Binary));
15655c6c1daeSBarry Smith #endif
15665c6c1daeSBarry Smith   PetscFunctionReturn(0);
15675c6c1daeSBarry Smith }
15685c6c1daeSBarry Smith 
15695c6c1daeSBarry Smith /* ---------------------------------------------------------------------*/
15705c6c1daeSBarry Smith /*
15715c6c1daeSBarry Smith     The variable Petsc_Viewer_Binary_keyval is used to indicate an MPI attribute that
15725c6c1daeSBarry Smith   is attached to a communicator, in this case the attribute is a PetscViewer.
15735c6c1daeSBarry Smith */
1574d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_Binary_keyval = MPI_KEYVAL_INVALID;
15755c6c1daeSBarry Smith 
15765c6c1daeSBarry Smith /*@C
1577811af0c4SBarry Smith      PETSC_VIEWER_BINARY_ - Creates a `PETSCVIEWERBINARY` `PetscViewer` shared by all processors
15785c6c1daeSBarry Smith                      in a communicator.
15795c6c1daeSBarry Smith 
1580d083f849SBarry Smith      Collective
15815c6c1daeSBarry Smith 
15825c6c1daeSBarry Smith      Input Parameter:
1583811af0c4SBarry Smith .    comm - the MPI communicator to share the `PETSCVIEWERBINARY`
15845c6c1daeSBarry Smith 
15855c6c1daeSBarry Smith      Level: intermediate
15865c6c1daeSBarry Smith 
15875c6c1daeSBarry Smith    Options Database Keys:
158810699b91SBarry Smith +    -viewer_binary_filename <name> - filename in which to store the binary data, defaults to binaryoutput
158910699b91SBarry Smith .    -viewer_binary_skip_info - true means do not create .info file for this viewer
159010699b91SBarry Smith .    -viewer_binary_skip_options - true means do not use the options database for this viewer
159110699b91SBarry Smith .    -viewer_binary_skip_header - true means do not store the usual header information in the binary file
159210699b91SBarry Smith -    -viewer_binary_mpiio - true means use the file via MPI-IO, maybe faster for large files and many MPI ranks
15935c6c1daeSBarry Smith 
1594811af0c4SBarry Smith    Environmental variable:
159510699b91SBarry Smith -   PETSC_VIEWER_BINARY_FILENAME - filename in which to store the binary data, defaults to binaryoutput
15965c6c1daeSBarry Smith 
1597811af0c4SBarry Smith      Note:
1598811af0c4SBarry Smith      Unlike almost all other PETSc routines, `PETSC_VIEWER_BINARY_` does not return
15995c6c1daeSBarry Smith      an error code.  The binary PetscViewer is usually used in the form
16005c6c1daeSBarry Smith $       XXXView(XXX object,PETSC_VIEWER_BINARY_(comm));
16015c6c1daeSBarry Smith 
1602811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PETSC_VIEWER_BINARY_WORLD`, `PETSC_VIEWER_BINARY_SELF`, `PetscViewerBinaryOpen()`, `PetscViewerCreate()`,
1603db781477SPatrick Sanan           `PetscViewerDestroy()`
16045c6c1daeSBarry Smith @*/
1605*d71ae5a4SJacob Faibussowitsch PetscViewer PETSC_VIEWER_BINARY_(MPI_Comm comm)
1606*d71ae5a4SJacob Faibussowitsch {
16075c6c1daeSBarry Smith   PetscErrorCode ierr;
16085c6c1daeSBarry Smith   PetscBool      flg;
16095c6c1daeSBarry Smith   PetscViewer    viewer;
16105c6c1daeSBarry Smith   char           fname[PETSC_MAX_PATH_LEN];
16115c6c1daeSBarry Smith   MPI_Comm       ncomm;
16125c6c1daeSBarry Smith 
16135c6c1daeSBarry Smith   PetscFunctionBegin;
16149371c9d4SSatish Balay   ierr = PetscCommDuplicate(comm, &ncomm, NULL);
16159371c9d4SSatish Balay   if (ierr) {
16169371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16179371c9d4SSatish Balay     PetscFunctionReturn(NULL);
16189371c9d4SSatish Balay   }
16195c6c1daeSBarry Smith   if (Petsc_Viewer_Binary_keyval == MPI_KEYVAL_INVALID) {
162002c9f0b5SLisandro Dalcin     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Viewer_Binary_keyval, NULL);
16219371c9d4SSatish Balay     if (ierr) {
16229371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16239371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16249371c9d4SSatish Balay     }
16255c6c1daeSBarry Smith   }
162647435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm, Petsc_Viewer_Binary_keyval, (void **)&viewer, (int *)&flg);
16279371c9d4SSatish Balay   if (ierr) {
16289371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16299371c9d4SSatish Balay     PetscFunctionReturn(NULL);
16309371c9d4SSatish Balay   }
16315c6c1daeSBarry Smith   if (!flg) { /* PetscViewer not yet created */
16325c6c1daeSBarry Smith     ierr = PetscOptionsGetenv(ncomm, "PETSC_VIEWER_BINARY_FILENAME", fname, PETSC_MAX_PATH_LEN, &flg);
16339371c9d4SSatish Balay     if (ierr) {
16349371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16359371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16369371c9d4SSatish Balay     }
16375c6c1daeSBarry Smith     if (!flg) {
16385c6c1daeSBarry Smith       ierr = PetscStrcpy(fname, "binaryoutput");
16399371c9d4SSatish Balay       if (ierr) {
16409371c9d4SSatish Balay         PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16419371c9d4SSatish Balay         PetscFunctionReturn(NULL);
16429371c9d4SSatish Balay       }
16435c6c1daeSBarry Smith     }
16445c6c1daeSBarry Smith     ierr = PetscViewerBinaryOpen(ncomm, fname, FILE_MODE_WRITE, &viewer);
16459371c9d4SSatish Balay     if (ierr) {
16469371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16479371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16489371c9d4SSatish Balay     }
16495c6c1daeSBarry Smith     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
16509371c9d4SSatish Balay     if (ierr) {
16519371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16529371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16539371c9d4SSatish Balay     }
165447435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm, Petsc_Viewer_Binary_keyval, (void *)viewer);
16559371c9d4SSatish Balay     if (ierr) {
16569371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16579371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16589371c9d4SSatish Balay     }
16595c6c1daeSBarry Smith   }
16605c6c1daeSBarry Smith   ierr = PetscCommDestroy(&ncomm);
16619371c9d4SSatish Balay   if (ierr) {
16629371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16639371c9d4SSatish Balay     PetscFunctionReturn(NULL);
16649371c9d4SSatish Balay   }
16655c6c1daeSBarry Smith   PetscFunctionReturn(viewer);
16665c6c1daeSBarry Smith }
1667