xref: /petsc/src/sys/classes/viewer/impls/binary/binv.c (revision 7666791831c4618c39250d1ff427002f918abde9)
1af0996ceSBarry Smith #include <petsc/private/viewerimpl.h> /*I   "petscviewer.h"   I*/
25c6c1daeSBarry Smith 
3*76667918SBarry Smith /*
4*76667918SBarry Smith    This needs to start the same as PetscViewer_Socket.
5*76667918SBarry Smith */
65c6c1daeSBarry Smith typedef struct {
75c6c1daeSBarry Smith   int       fdes;        /* file descriptor, ignored if using MPI IO */
8*76667918SBarry Smith   PetscInt  flowcontrol; /* allow only <flowcontrol> messages outstanding at a time while doing IO */
9*76667918SBarry Smith   PetscBool skipheader;  /* don't write header, only raw data */
105c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
11bc196f7cSDave May   PetscBool  usempiio;
125c6c1daeSBarry Smith   MPI_File   mfdes; /* ignored unless using MPI IO */
13e8a65b7dSLisandro Dalcin   MPI_File   mfsub; /* subviewer support */
145c6c1daeSBarry Smith   MPI_Offset moff;
155c6c1daeSBarry Smith #endif
16cc843e7aSLisandro Dalcin   char         *filename;            /* file name */
17cc843e7aSLisandro Dalcin   PetscFileMode filemode;            /* read/write/append mode */
185c6c1daeSBarry Smith   FILE         *fdes_info;           /* optional file containing info on binary file*/
195c6c1daeSBarry Smith   PetscBool     storecompressed;     /* gzip the write binary file when closing it*/
20f90597f1SStefano Zampini   char         *ogzfilename;         /* gzip can be run after the filename has been updated */
215c6c1daeSBarry Smith   PetscBool     skipinfo;            /* Don't create info file for writing; don't use for reading */
225c6c1daeSBarry Smith   PetscBool     skipoptions;         /* don't use PETSc options database when loading */
23a261c58fSBarry Smith   PetscBool     matlabheaderwritten; /* if format is PETSC_VIEWER_BINARY_MATLAB has the MATLAB .info header been written yet */
24c98fd787SBarry Smith   PetscBool     setfromoptionscalled;
255c6c1daeSBarry Smith } PetscViewer_Binary;
265c6c1daeSBarry Smith 
27d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryClearFunctionList(PetscViewer v)
28d71ae5a4SJacob Faibussowitsch {
292e956fe4SStefano Zampini   PetscFunctionBegin;
302e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetFlowControl_C", NULL));
312e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetFlowControl_C", NULL));
322e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipHeader_C", NULL));
332e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipHeader_C", NULL));
342e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipOptions_C", NULL));
352e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipOptions_C", NULL));
362e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipInfo_C", NULL));
372e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipInfo_C", NULL));
382e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetInfoPointer_C", NULL));
392e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", NULL));
402e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", NULL));
412e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", NULL));
422e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", NULL));
432e956fe4SStefano Zampini #if defined(PETSC_HAVE_MPIIO)
442e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetUseMPIIO_C", NULL));
452e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetUseMPIIO_C", NULL));
462e956fe4SStefano Zampini #endif
472e956fe4SStefano Zampini   PetscFunctionReturn(0);
482e956fe4SStefano Zampini }
492e956fe4SStefano Zampini 
50cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
51d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySyncMPIIO(PetscViewer viewer)
52d71ae5a4SJacob Faibussowitsch {
53cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
54cc843e7aSLisandro Dalcin 
55cc843e7aSLisandro Dalcin   PetscFunctionBegin;
56cc843e7aSLisandro Dalcin   if (vbinary->filemode == FILE_MODE_READ) PetscFunctionReturn(0);
5748a46eb9SPierre Jolivet   if (vbinary->mfsub != MPI_FILE_NULL) PetscCallMPI(MPI_File_sync(vbinary->mfsub));
58cc843e7aSLisandro Dalcin   if (vbinary->mfdes != MPI_FILE_NULL) {
599566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer)));
609566063dSJacob Faibussowitsch     PetscCallMPI(MPI_File_sync(vbinary->mfdes));
61cc843e7aSLisandro Dalcin   }
62cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
63cc843e7aSLisandro Dalcin }
64cc843e7aSLisandro Dalcin #endif
65cc843e7aSLisandro Dalcin 
66d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerGetSubViewer_Binary(PetscViewer viewer, MPI_Comm comm, PetscViewer *outviewer)
67d71ae5a4SJacob Faibussowitsch {
68e8a65b7dSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
69cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
705c6c1daeSBarry Smith 
715c6c1daeSBarry Smith   PetscFunctionBegin;
729566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
73e8a65b7dSLisandro Dalcin 
74e8a65b7dSLisandro Dalcin   /* Return subviewer in process zero */
759566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
76dd400576SPatrick Sanan   if (rank == 0) {
77e5afcf28SBarry Smith     PetscMPIInt flg;
78e5afcf28SBarry Smith 
799566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_compare(PETSC_COMM_SELF, comm, &flg));
80cc73adaaSBarry Smith     PetscCheck(flg == MPI_IDENT || flg == MPI_CONGRUENT, PETSC_COMM_SELF, PETSC_ERR_SUP, "PetscViewerGetSubViewer() for PETSCVIEWERBINARY requires a singleton MPI_Comm");
819566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(comm, outviewer));
829566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(*outviewer, PETSCVIEWERBINARY));
839566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy((*outviewer)->data, vbinary, sizeof(PetscViewer_Binary)));
8412f4c3a9SLisandro Dalcin     (*outviewer)->setupcalled = PETSC_TRUE;
8596509d17SLisandro Dalcin   } else {
8696509d17SLisandro Dalcin     *outviewer = NULL;
8796509d17SLisandro Dalcin   }
88e8a65b7dSLisandro Dalcin 
89e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
90e8a65b7dSLisandro Dalcin   if (vbinary->usempiio && *outviewer) {
91e8a65b7dSLisandro Dalcin     PetscViewer_Binary *obinary = (PetscViewer_Binary *)(*outviewer)->data;
92e8a65b7dSLisandro Dalcin     /* Parent viewer opens a new MPI file handle on PETSC_COMM_SELF and keeps track of it for future reuse */
93e8a65b7dSLisandro Dalcin     if (vbinary->mfsub == MPI_FILE_NULL) {
94e8a65b7dSLisandro Dalcin       int amode;
95cc843e7aSLisandro Dalcin       switch (vbinary->filemode) {
96d71ae5a4SJacob Faibussowitsch       case FILE_MODE_READ:
97d71ae5a4SJacob Faibussowitsch         amode = MPI_MODE_RDONLY;
98d71ae5a4SJacob Faibussowitsch         break;
99d71ae5a4SJacob Faibussowitsch       case FILE_MODE_WRITE:
100d71ae5a4SJacob Faibussowitsch         amode = MPI_MODE_WRONLY;
101d71ae5a4SJacob Faibussowitsch         break;
102d71ae5a4SJacob Faibussowitsch       case FILE_MODE_APPEND:
103d71ae5a4SJacob Faibussowitsch         amode = MPI_MODE_WRONLY;
104d71ae5a4SJacob Faibussowitsch         break;
105d71ae5a4SJacob Faibussowitsch       default:
106d71ae5a4SJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[vbinary->filemode]);
107e8a65b7dSLisandro Dalcin       }
1089566063dSJacob Faibussowitsch       PetscCallMPI(MPI_File_open(PETSC_COMM_SELF, vbinary->filename, amode, MPI_INFO_NULL, &vbinary->mfsub));
109e8a65b7dSLisandro Dalcin     }
110e8a65b7dSLisandro Dalcin     /* Subviewer gets the MPI file handle on PETSC_COMM_SELF */
111e8a65b7dSLisandro Dalcin     obinary->mfdes = vbinary->mfsub;
112e8a65b7dSLisandro Dalcin     obinary->mfsub = MPI_FILE_NULL;
113e8a65b7dSLisandro Dalcin     obinary->moff  = vbinary->moff;
114e8a65b7dSLisandro Dalcin   }
115e8a65b7dSLisandro Dalcin #endif
116cc843e7aSLisandro Dalcin 
117cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
1189566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinarySyncMPIIO(viewer));
119cc843e7aSLisandro Dalcin #endif
1205c6c1daeSBarry Smith   PetscFunctionReturn(0);
1215c6c1daeSBarry Smith }
1225c6c1daeSBarry Smith 
123d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerRestoreSubViewer_Binary(PetscViewer viewer, MPI_Comm comm, PetscViewer *outviewer)
124d71ae5a4SJacob Faibussowitsch {
125e8a65b7dSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
126cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
127e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
128e8a65b7dSLisandro Dalcin   MPI_Offset moff = 0;
129e8a65b7dSLisandro Dalcin #endif
1305c6c1daeSBarry Smith 
1315c6c1daeSBarry Smith   PetscFunctionBegin;
1329566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
133c5853193SPierre Jolivet   PetscCheck(rank == 0 || !*outviewer, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Subviewer not obtained from viewer");
134e8a65b7dSLisandro Dalcin 
135e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
136e8a65b7dSLisandro Dalcin   if (vbinary->usempiio && *outviewer) {
137e8a65b7dSLisandro Dalcin     PetscViewer_Binary *obinary = (PetscViewer_Binary *)(*outviewer)->data;
13808401ef6SPierre Jolivet     PetscCheck(obinary->mfdes == vbinary->mfsub, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Subviewer not obtained from viewer");
1399566063dSJacob Faibussowitsch     if (obinary->mfsub != MPI_FILE_NULL) PetscCallMPI(MPI_File_close(&obinary->mfsub));
140e8a65b7dSLisandro Dalcin     moff = obinary->moff;
141e8a65b7dSLisandro Dalcin   }
142e8a65b7dSLisandro Dalcin #endif
143e8a65b7dSLisandro Dalcin 
144e8a65b7dSLisandro Dalcin   if (*outviewer) {
145e8a65b7dSLisandro Dalcin     PetscViewer_Binary *obinary = (PetscViewer_Binary *)(*outviewer)->data;
14608401ef6SPierre Jolivet     PetscCheck(obinary->fdes == vbinary->fdes, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Subviewer not obtained from viewer");
1479566063dSJacob Faibussowitsch     PetscCall(PetscFree((*outviewer)->data));
1482e956fe4SStefano Zampini     PetscCall(PetscViewerBinaryClearFunctionList(*outviewer));
1499566063dSJacob Faibussowitsch     PetscCall(PetscHeaderDestroy(outviewer));
1505c6c1daeSBarry Smith   }
151e8a65b7dSLisandro Dalcin 
152e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
153e8a65b7dSLisandro Dalcin   if (vbinary->usempiio) {
154e8a65b7dSLisandro Dalcin     PetscInt64 ioff = (PetscInt64)moff; /* We could use MPI_OFFSET datatype (requires MPI 2.2) */
1559566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&ioff, 1, MPIU_INT64, 0, PetscObjectComm((PetscObject)viewer)));
156e8a65b7dSLisandro Dalcin     vbinary->moff = (MPI_Offset)ioff;
157e8a65b7dSLisandro Dalcin   }
158e8a65b7dSLisandro Dalcin #endif
159cc843e7aSLisandro Dalcin 
160cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
1619566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinarySyncMPIIO(viewer));
162cc843e7aSLisandro Dalcin #endif
1635c6c1daeSBarry Smith   PetscFunctionReturn(0);
1645c6c1daeSBarry Smith }
1655c6c1daeSBarry Smith 
1665c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
1675c6c1daeSBarry Smith /*@C
16885b8072bSPatrick Sanan     PetscViewerBinaryGetMPIIOOffset - Gets the current global offset that should be passed to `MPI_File_set_view()` or `MPI_File_{write|read}_at[_all]()`
1695c6c1daeSBarry Smith 
1705c6c1daeSBarry Smith     Not Collective
1715c6c1daeSBarry Smith 
1725c6c1daeSBarry Smith     Input Parameter:
17385b8072bSPatrick Sanan .   viewer - PetscViewer context, obtained from `PetscViewerBinaryOpen()`
1745c6c1daeSBarry Smith 
1755c6c1daeSBarry Smith     Output Parameter:
17622a8f86cSLisandro Dalcin .   off - the current global offset
1775c6c1daeSBarry Smith 
1785c6c1daeSBarry Smith     Level: advanced
1795c6c1daeSBarry Smith 
180811af0c4SBarry Smith     Note:
181811af0c4SBarry Smith     Use `PetscViewerBinaryAddMPIIOOffset()` to increase this value after you have written a view.
182811af0c4SBarry Smith 
1835c6c1daeSBarry Smith     Fortran Note:
1845c6c1daeSBarry Smith     This routine is not supported in Fortran.
1855c6c1daeSBarry Smith 
186811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryAddMPIIOOffset()`
1875c6c1daeSBarry Smith @*/
188d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetMPIIOOffset(PetscViewer viewer, MPI_Offset *off)
189d71ae5a4SJacob Faibussowitsch {
19022a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
1915c6c1daeSBarry Smith 
1925c6c1daeSBarry Smith   PetscFunctionBegin;
19322a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
19422a8f86cSLisandro Dalcin   PetscValidPointer(off, 2);
19522a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
1965c6c1daeSBarry Smith   *off    = vbinary->moff;
1975c6c1daeSBarry Smith   PetscFunctionReturn(0);
1985c6c1daeSBarry Smith }
1995c6c1daeSBarry Smith 
2005c6c1daeSBarry Smith /*@C
20122a8f86cSLisandro Dalcin     PetscViewerBinaryAddMPIIOOffset - Adds to the current global offset
2025c6c1daeSBarry Smith 
20322a8f86cSLisandro Dalcin     Logically Collective
2045c6c1daeSBarry Smith 
2055c6c1daeSBarry Smith     Input Parameters:
206811af0c4SBarry Smith +   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
20722a8f86cSLisandro Dalcin -   off - the addition to the global offset
2085c6c1daeSBarry Smith 
2095c6c1daeSBarry Smith     Level: advanced
2105c6c1daeSBarry Smith 
211811af0c4SBarry Smith     Note:
212811af0c4SBarry Smith     Use `PetscViewerBinaryGetMPIIOOffset()` to get the value that you should pass to `MPI_File_set_view()` or `MPI_File_{write|read}_at[_all]()`
213811af0c4SBarry Smith 
2145c6c1daeSBarry Smith     Fortran Note:
2155c6c1daeSBarry Smith     This routine is not supported in Fortran.
2165c6c1daeSBarry Smith 
217811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
2185c6c1daeSBarry Smith @*/
219d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryAddMPIIOOffset(PetscViewer viewer, MPI_Offset off)
220d71ae5a4SJacob Faibussowitsch {
22122a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
2225c6c1daeSBarry Smith 
2235c6c1daeSBarry Smith   PetscFunctionBegin;
22422a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
22522a8f86cSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, (PetscInt)off, 2);
22622a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
2275c6c1daeSBarry Smith   vbinary->moff += off;
2285c6c1daeSBarry Smith   PetscFunctionReturn(0);
2295c6c1daeSBarry Smith }
2305c6c1daeSBarry Smith 
2315c6c1daeSBarry Smith /*@C
232811af0c4SBarry Smith     PetscViewerBinaryGetMPIIODescriptor - Extracts the MPI IO file descriptor from a `PetscViewer`.
2335c6c1daeSBarry Smith 
2345c6c1daeSBarry Smith     Not Collective
2355c6c1daeSBarry Smith 
2365c6c1daeSBarry Smith     Input Parameter:
237811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
2385c6c1daeSBarry Smith 
2395c6c1daeSBarry Smith     Output Parameter:
2405c6c1daeSBarry Smith .   fdes - file descriptor
2415c6c1daeSBarry Smith 
2425c6c1daeSBarry Smith     Level: advanced
2435c6c1daeSBarry Smith 
2445c6c1daeSBarry Smith     Fortran Note:
2455c6c1daeSBarry Smith     This routine is not supported in Fortran.
2465c6c1daeSBarry Smith 
247811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
2485c6c1daeSBarry Smith @*/
249d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetMPIIODescriptor(PetscViewer viewer, MPI_File *fdes)
250d71ae5a4SJacob Faibussowitsch {
25122a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
2525c6c1daeSBarry Smith 
2535c6c1daeSBarry Smith   PetscFunctionBegin;
25422a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
25522a8f86cSLisandro Dalcin   PetscValidPointer(fdes, 2);
2569566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
25722a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
2585c6c1daeSBarry Smith   *fdes   = vbinary->mfdes;
2595c6c1daeSBarry Smith   PetscFunctionReturn(0);
2605c6c1daeSBarry Smith }
261cc843e7aSLisandro Dalcin #endif
2625c6c1daeSBarry Smith 
263cc843e7aSLisandro Dalcin /*@
264cc843e7aSLisandro Dalcin     PetscViewerBinarySetUseMPIIO - Sets a binary viewer to use MPI-IO for reading/writing. Must be called
265811af0c4SBarry Smith         before `PetscViewerFileSetName()`
266cc843e7aSLisandro Dalcin 
267811af0c4SBarry Smith     Logically Collective on viewer
268cc843e7aSLisandro Dalcin 
269cc843e7aSLisandro Dalcin     Input Parameters:
270811af0c4SBarry Smith +   viewer - the `PetscViewer`; must be a `PETSCVIEWERBINARY`
271811af0c4SBarry Smith -   use - `PETSC_TRUE` means MPI-IO will be used
272cc843e7aSLisandro Dalcin 
273811af0c4SBarry Smith     Options Database Key:
274cc843e7aSLisandro Dalcin     -viewer_binary_mpiio : Flag for using MPI-IO
275cc843e7aSLisandro Dalcin 
276cc843e7aSLisandro Dalcin     Level: advanced
277cc843e7aSLisandro Dalcin 
278811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`,
279db781477SPatrick Sanan           `PetscViewerBinaryGetUseMPIIO()`
280cc843e7aSLisandro Dalcin @*/
281d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetUseMPIIO(PetscViewer viewer, PetscBool use)
282d71ae5a4SJacob Faibussowitsch {
283a8aae444SDave May   PetscFunctionBegin;
284cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
285cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, use, 2);
286cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetUseMPIIO_C", (PetscViewer, PetscBool), (viewer, use));
287cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
288cc843e7aSLisandro Dalcin }
289cc843e7aSLisandro Dalcin 
290cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
291d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetUseMPIIO_Binary(PetscViewer viewer, PetscBool use)
292d71ae5a4SJacob Faibussowitsch {
293cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
294cc843e7aSLisandro Dalcin   PetscFunctionBegin;
295cc73adaaSBarry Smith   PetscCheck(!viewer->setupcalled || vbinary->usempiio == use, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Cannot change MPIIO to %s after setup", PetscBools[use]);
296cc843e7aSLisandro Dalcin   vbinary->usempiio = use;
297a8aae444SDave May   PetscFunctionReturn(0);
298a8aae444SDave May }
299a8aae444SDave May #endif
300a8aae444SDave May 
301cc843e7aSLisandro Dalcin /*@
302bc196f7cSDave May     PetscViewerBinaryGetUseMPIIO - Returns PETSC_TRUE if the binary viewer uses MPI-IO.
3035c6c1daeSBarry Smith 
3045c6c1daeSBarry Smith     Not Collective
3055c6c1daeSBarry Smith 
3065c6c1daeSBarry Smith     Input Parameter:
307811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`; must be a `PETSCVIEWERBINARY`
3085c6c1daeSBarry Smith 
3095c6c1daeSBarry Smith     Output Parameter:
310811af0c4SBarry Smith .   use - `PETSC_TRUE` if MPI-IO is being used
3115c6c1daeSBarry Smith 
3125c6c1daeSBarry Smith     Level: advanced
3135c6c1daeSBarry Smith 
314bc196f7cSDave May     Note:
315bc196f7cSDave May     If MPI-IO is not available, this function will always return PETSC_FALSE
316bc196f7cSDave May 
317811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
3185c6c1daeSBarry Smith @*/
319d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetUseMPIIO(PetscViewer viewer, PetscBool *use)
320d71ae5a4SJacob Faibussowitsch {
3215c6c1daeSBarry Smith   PetscFunctionBegin;
322cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
323cc843e7aSLisandro Dalcin   PetscValidBoolPointer(use, 2);
324cc843e7aSLisandro Dalcin   *use = PETSC_FALSE;
325cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinaryGetUseMPIIO_C", (PetscViewer, PetscBool *), (viewer, use));
3265c6c1daeSBarry Smith   PetscFunctionReturn(0);
3275c6c1daeSBarry Smith }
3285c6c1daeSBarry Smith 
329cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
330d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetUseMPIIO_Binary(PetscViewer viewer, PetscBool *use)
331d71ae5a4SJacob Faibussowitsch {
3325c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
3335c6c1daeSBarry Smith 
3345c6c1daeSBarry Smith   PetscFunctionBegin;
335cc843e7aSLisandro Dalcin   *use = vbinary->usempiio;
336cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
337cc843e7aSLisandro Dalcin }
338cc843e7aSLisandro Dalcin #endif
339cc843e7aSLisandro Dalcin 
340cc843e7aSLisandro Dalcin /*@
341cc843e7aSLisandro Dalcin     PetscViewerBinarySetFlowControl - Sets how many messages are allowed to outstanding at the same time during parallel IO reads/writes
342cc843e7aSLisandro Dalcin 
343cc843e7aSLisandro Dalcin     Not Collective
344cc843e7aSLisandro Dalcin 
345d8d19677SJose E. Roman     Input Parameters:
346811af0c4SBarry Smith +   viewer - PetscViewer context, obtained from `PetscViewerBinaryOpen()`
347cc843e7aSLisandro Dalcin -   fc - the number of messages, defaults to 256 if this function was not called
348cc843e7aSLisandro Dalcin 
349cc843e7aSLisandro Dalcin     Level: advanced
350cc843e7aSLisandro Dalcin 
351811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetFlowControl()`
352cc843e7aSLisandro Dalcin @*/
353d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetFlowControl(PetscViewer viewer, PetscInt fc)
354d71ae5a4SJacob Faibussowitsch {
355cc843e7aSLisandro Dalcin   PetscFunctionBegin;
356cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
357cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, fc, 2);
358cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetFlowControl_C", (PetscViewer, PetscInt), (viewer, fc));
3595c6c1daeSBarry Smith   PetscFunctionReturn(0);
3605c6c1daeSBarry Smith }
3615c6c1daeSBarry Smith 
362d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetFlowControl_Binary(PetscViewer viewer, PetscInt fc)
363d71ae5a4SJacob Faibussowitsch {
364cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
365cc843e7aSLisandro Dalcin 
366cc843e7aSLisandro Dalcin   PetscFunctionBegin;
36708401ef6SPierre Jolivet   PetscCheck(fc > 1, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_OUTOFRANGE, "Flow control count must be greater than 1, %" PetscInt_FMT " was set", fc);
368cc843e7aSLisandro Dalcin   vbinary->flowcontrol = fc;
369cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
370cc843e7aSLisandro Dalcin }
371cc843e7aSLisandro Dalcin 
372cc843e7aSLisandro Dalcin /*@
3735c6c1daeSBarry Smith     PetscViewerBinaryGetFlowControl - Returns how many messages are allowed to outstanding at the same time during parallel IO reads/writes
3745c6c1daeSBarry Smith 
3755c6c1daeSBarry Smith     Not Collective
3765c6c1daeSBarry Smith 
3775c6c1daeSBarry Smith     Input Parameter:
378811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
3795c6c1daeSBarry Smith 
3805c6c1daeSBarry Smith     Output Parameter:
3815c6c1daeSBarry Smith .   fc - the number of messages
3825c6c1daeSBarry Smith 
3835c6c1daeSBarry Smith     Level: advanced
3845c6c1daeSBarry Smith 
385811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinarySetFlowControl()`
3865c6c1daeSBarry Smith @*/
387d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetFlowControl(PetscViewer viewer, PetscInt *fc)
388d71ae5a4SJacob Faibussowitsch {
3895c6c1daeSBarry Smith   PetscFunctionBegin;
390cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
391cc843e7aSLisandro Dalcin   PetscValidIntPointer(fc, 2);
392cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetFlowControl_C", (PetscViewer, PetscInt *), (viewer, fc));
3935c6c1daeSBarry Smith   PetscFunctionReturn(0);
3945c6c1daeSBarry Smith }
3955c6c1daeSBarry Smith 
396*76667918SBarry Smith PETSC_INTERN PetscErrorCode PetscViewerBinaryGetFlowControl_Binary(PetscViewer viewer, PetscInt *fc)
397d71ae5a4SJacob Faibussowitsch {
3985c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
3995c6c1daeSBarry Smith 
4005c6c1daeSBarry Smith   PetscFunctionBegin;
401cc843e7aSLisandro Dalcin   *fc = vbinary->flowcontrol;
4025c6c1daeSBarry Smith   PetscFunctionReturn(0);
4035c6c1daeSBarry Smith }
4045c6c1daeSBarry Smith 
4055c6c1daeSBarry Smith /*@C
406811af0c4SBarry Smith     PetscViewerBinaryGetDescriptor - Extracts the file descriptor from a `PetscViewer` of `PetscViewerType` `PETSCVIEWERBINARY`.
4075c6c1daeSBarry Smith 
408811af0c4SBarry Smith     Collective on viewer because it may trigger a `PetscViewerSetUp()` call
4095c6c1daeSBarry Smith 
4105c6c1daeSBarry Smith     Input Parameter:
411811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
4125c6c1daeSBarry Smith 
4135c6c1daeSBarry Smith     Output Parameter:
4145c6c1daeSBarry Smith .   fdes - file descriptor
4155c6c1daeSBarry Smith 
4165c6c1daeSBarry Smith     Level: advanced
4175c6c1daeSBarry Smith 
418811af0c4SBarry Smith     Note:
419811af0c4SBarry Smith       For writable binary `PetscViewer`s, the descriptor will only be valid for the
420811af0c4SBarry Smith     first processor in the communicator that shares the `PetscViewer`. For readable
4215c6c1daeSBarry Smith     files it will only be valid on nodes that have the file. If node 0 does not
4225c6c1daeSBarry Smith     have the file it generates an error even if another node does have the file.
4235c6c1daeSBarry Smith 
4245c6c1daeSBarry Smith     Fortran Note:
4255c6c1daeSBarry Smith     This routine is not supported in Fortran.
4265c6c1daeSBarry Smith 
427811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`
4285c6c1daeSBarry Smith @*/
429d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetDescriptor(PetscViewer viewer, int *fdes)
430d71ae5a4SJacob Faibussowitsch {
43122a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
4325c6c1daeSBarry Smith 
4335c6c1daeSBarry Smith   PetscFunctionBegin;
43422a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
43522a8f86cSLisandro Dalcin   PetscValidPointer(fdes, 2);
4369566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
43722a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
4385c6c1daeSBarry Smith   *fdes   = vbinary->fdes;
4395c6c1daeSBarry Smith   PetscFunctionReturn(0);
4405c6c1daeSBarry Smith }
4415c6c1daeSBarry Smith 
4425c6c1daeSBarry Smith /*@
4435c6c1daeSBarry Smith     PetscViewerBinarySkipInfo - Binary file will not have .info file created with it
4445c6c1daeSBarry Smith 
4455c6c1daeSBarry Smith     Not Collective
4465c6c1daeSBarry Smith 
447fd292e60Sprj-     Input Parameter:
448811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerCreate()`
4495c6c1daeSBarry Smith 
4505c6c1daeSBarry Smith     Options Database Key:
45110699b91SBarry Smith .   -viewer_binary_skip_info - true indicates do not generate .info file
4525c6c1daeSBarry Smith 
4535c6c1daeSBarry Smith     Level: advanced
4545c6c1daeSBarry Smith 
45595452b02SPatrick Sanan     Notes:
456811af0c4SBarry Smith     This must be called after `PetscViewerSetType()`. If you use `PetscViewerBinaryOpen()` then
4575c6c1daeSBarry Smith     you can only skip the info file with the -viewer_binary_skip_info flag. To use the function you must open the
458811af0c4SBarry Smith     viewer with `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinarySkipInfo()`.
4595c6c1daeSBarry Smith 
4605c6c1daeSBarry Smith     The .info contains meta information about the data in the binary file, for example the block size if it was
4615c6c1daeSBarry Smith     set for a vector or matrix.
4625c6c1daeSBarry Smith 
463811af0c4SBarry Smith     This routine is deprecated, use `PetscViewerBinarySetSkipInfo()`
464811af0c4SBarry Smith 
465811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySetSkipOptions()`,
466db781477SPatrick Sanan           `PetscViewerBinaryGetSkipOptions()`, `PetscViewerBinaryGetSkipInfo()`
4675c6c1daeSBarry Smith @*/
468d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySkipInfo(PetscViewer viewer)
469d71ae5a4SJacob Faibussowitsch {
4705c6c1daeSBarry Smith   PetscFunctionBegin;
4719566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinarySetSkipInfo(viewer, PETSC_TRUE));
472807ea322SDave May   PetscFunctionReturn(0);
473807ea322SDave May }
474807ea322SDave May 
475807ea322SDave May /*@
476807ea322SDave May     PetscViewerBinarySetSkipInfo - Binary file will not have .info file created with it
477807ea322SDave May 
478807ea322SDave May     Not Collective
479807ea322SDave May 
480d8d19677SJose E. Roman     Input Parameters:
481cc843e7aSLisandro Dalcin +   viewer - PetscViewer context, obtained from PetscViewerCreate()
482cc843e7aSLisandro Dalcin -   skip - PETSC_TRUE implies the .info file will not be generated
483807ea322SDave May 
484807ea322SDave May     Options Database Key:
48510699b91SBarry Smith .   -viewer_binary_skip_info - true indicates do not generate .info file
486807ea322SDave May 
487807ea322SDave May     Level: advanced
488807ea322SDave May 
489811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySetSkipOptions()`,
490db781477SPatrick Sanan           `PetscViewerBinaryGetSkipOptions()`, `PetscViewerBinaryGetSkipInfo()`
491807ea322SDave May @*/
492d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetSkipInfo(PetscViewer viewer, PetscBool skip)
493d71ae5a4SJacob Faibussowitsch {
494807ea322SDave May   PetscFunctionBegin;
495cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
496cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, skip, 2);
497cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetSkipInfo_C", (PetscViewer, PetscBool), (viewer, skip));
498807ea322SDave May   PetscFunctionReturn(0);
499807ea322SDave May }
500807ea322SDave May 
501d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetSkipInfo_Binary(PetscViewer viewer, PetscBool skip)
502d71ae5a4SJacob Faibussowitsch {
503807ea322SDave May   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
504807ea322SDave May 
505807ea322SDave May   PetscFunctionBegin;
506cc843e7aSLisandro Dalcin   vbinary->skipinfo = skip;
507807ea322SDave May   PetscFunctionReturn(0);
508807ea322SDave May }
509807ea322SDave May 
510807ea322SDave May /*@
511807ea322SDave May     PetscViewerBinaryGetSkipInfo - check if viewer wrote a .info file
512807ea322SDave May 
513807ea322SDave May     Not Collective
514807ea322SDave May 
515807ea322SDave May     Input Parameter:
516811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
517807ea322SDave May 
518807ea322SDave May     Output Parameter:
519811af0c4SBarry Smith .   skip - `PETSC_TRUE` implies the .info file was not generated
520807ea322SDave May 
521807ea322SDave May     Level: advanced
522807ea322SDave May 
523811af0c4SBarry Smith     Note:
524811af0c4SBarry Smith     This must be called after `PetscViewerSetType()`
525807ea322SDave May 
526811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`,
527db781477SPatrick Sanan           `PetscViewerBinarySetSkipOptions()`, `PetscViewerBinarySetSkipInfo()`
528807ea322SDave May @*/
529d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetSkipInfo(PetscViewer viewer, PetscBool *skip)
530d71ae5a4SJacob Faibussowitsch {
531807ea322SDave May   PetscFunctionBegin;
532cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
533cc843e7aSLisandro Dalcin   PetscValidBoolPointer(skip, 2);
534cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetSkipInfo_C", (PetscViewer, PetscBool *), (viewer, skip));
535807ea322SDave May   PetscFunctionReturn(0);
536807ea322SDave May }
537807ea322SDave May 
538d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetSkipInfo_Binary(PetscViewer viewer, PetscBool *skip)
539d71ae5a4SJacob Faibussowitsch {
540807ea322SDave May   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
541807ea322SDave May 
542807ea322SDave May   PetscFunctionBegin;
543cc843e7aSLisandro Dalcin   *skip = vbinary->skipinfo;
544807ea322SDave May   PetscFunctionReturn(0);
545807ea322SDave May }
546807ea322SDave May 
5475c6c1daeSBarry Smith /*@
5485c6c1daeSBarry Smith     PetscViewerBinarySetSkipOptions - do not use the PETSc options database when loading objects
5495c6c1daeSBarry Smith 
5505c6c1daeSBarry Smith     Not Collective
5515c6c1daeSBarry Smith 
5525c6c1daeSBarry Smith     Input Parameters:
553811af0c4SBarry Smith +   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
554811af0c4SBarry Smith -   skip - `PETSC_TRUE` means do not use the options from the options database
5555c6c1daeSBarry Smith 
5565c6c1daeSBarry Smith     Options Database Key:
557811af0c4SBarry Smith .   -viewer_binary_skip_options <true or false> - true means do not use the options from the options database
5585c6c1daeSBarry Smith 
5595c6c1daeSBarry Smith     Level: advanced
5605c6c1daeSBarry Smith 
561811af0c4SBarry Smith     Note:
562811af0c4SBarry Smith     This must be called after `PetscViewerSetType()`
5635c6c1daeSBarry Smith 
564811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
565db781477SPatrick Sanan           `PetscViewerBinaryGetSkipOptions()`
5665c6c1daeSBarry Smith @*/
567d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetSkipOptions(PetscViewer viewer, PetscBool skip)
568d71ae5a4SJacob Faibussowitsch {
5695c6c1daeSBarry Smith   PetscFunctionBegin;
570cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
571cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, skip, 2);
572cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetSkipOptions_C", (PetscViewer, PetscBool), (viewer, skip));
573807ea322SDave May   PetscFunctionReturn(0);
574807ea322SDave May }
575807ea322SDave May 
576d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetSkipOptions_Binary(PetscViewer viewer, PetscBool skip)
577d71ae5a4SJacob Faibussowitsch {
578cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
579807ea322SDave May 
580807ea322SDave May   PetscFunctionBegin;
581cc843e7aSLisandro Dalcin   vbinary->skipoptions = skip;
5825c6c1daeSBarry Smith   PetscFunctionReturn(0);
5835c6c1daeSBarry Smith }
5845c6c1daeSBarry Smith 
5855c6c1daeSBarry Smith /*@
5865c6c1daeSBarry Smith     PetscViewerBinaryGetSkipOptions - checks if viewer uses the PETSc options database when loading objects
5875c6c1daeSBarry Smith 
5885c6c1daeSBarry Smith     Not Collective
5895c6c1daeSBarry Smith 
5905c6c1daeSBarry Smith     Input Parameter:
591811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
5925c6c1daeSBarry Smith 
5935c6c1daeSBarry Smith     Output Parameter:
594811af0c4SBarry Smith .   skip - `PETSC_TRUE` means do not use
5955c6c1daeSBarry Smith 
5965c6c1daeSBarry Smith     Level: advanced
5975c6c1daeSBarry Smith 
598811af0c4SBarry Smith     Note:
599811af0c4SBarry Smith     This must be called after `PetscViewerSetType()`
6005c6c1daeSBarry Smith 
601811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
602db781477SPatrick Sanan           `PetscViewerBinarySetSkipOptions()`
6035c6c1daeSBarry Smith @*/
604d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetSkipOptions(PetscViewer viewer, PetscBool *skip)
605d71ae5a4SJacob Faibussowitsch {
6065c6c1daeSBarry Smith   PetscFunctionBegin;
607cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
608cc843e7aSLisandro Dalcin   PetscValidBoolPointer(skip, 2);
609cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetSkipOptions_C", (PetscViewer, PetscBool *), (viewer, skip));
6105c6c1daeSBarry Smith   PetscFunctionReturn(0);
6115c6c1daeSBarry Smith }
6125c6c1daeSBarry Smith 
613d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetSkipOptions_Binary(PetscViewer viewer, PetscBool *skip)
614d71ae5a4SJacob Faibussowitsch {
615d21b9a37SPierre Jolivet   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
6165c6c1daeSBarry Smith 
6175c6c1daeSBarry Smith   PetscFunctionBegin;
618cc843e7aSLisandro Dalcin   *skip = vbinary->skipoptions;
6195c6c1daeSBarry Smith   PetscFunctionReturn(0);
6205c6c1daeSBarry Smith }
6215c6c1daeSBarry Smith 
6225c6c1daeSBarry Smith /*@
6235c6c1daeSBarry Smith     PetscViewerBinarySetSkipHeader - do not write a header with size information on output, just raw data
6245c6c1daeSBarry Smith 
6255c6c1daeSBarry Smith     Not Collective
6265c6c1daeSBarry Smith 
6275c6c1daeSBarry Smith     Input Parameters:
628811af0c4SBarry Smith +   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
629811af0c4SBarry Smith -   skip - `PETSC_TRUE `means do not write header
6305c6c1daeSBarry Smith 
6315c6c1daeSBarry Smith     Options Database Key:
632811af0c4SBarry Smith .   -viewer_binary_skip_header <true or false> - true means do not write header
6335c6c1daeSBarry Smith 
6345c6c1daeSBarry Smith     Level: advanced
6355c6c1daeSBarry Smith 
63695452b02SPatrick Sanan     Notes:
637811af0c4SBarry Smith       This must be called after `PetscViewerSetType()`
6385c6c1daeSBarry Smith 
63910699b91SBarry Smith       Is ignored on anything but a binary viewer
6405c6c1daeSBarry Smith 
641811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
642db781477SPatrick Sanan           `PetscViewerBinaryGetSkipHeader()`
6435c6c1daeSBarry Smith @*/
644d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetSkipHeader(PetscViewer viewer, PetscBool skip)
645d71ae5a4SJacob Faibussowitsch {
6465c6c1daeSBarry Smith   PetscFunctionBegin;
647cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
648cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, skip, 2);
649cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetSkipHeader_C", (PetscViewer, PetscBool), (viewer, skip));
6505c6c1daeSBarry Smith   PetscFunctionReturn(0);
6515c6c1daeSBarry Smith }
6525c6c1daeSBarry Smith 
653d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetSkipHeader_Binary(PetscViewer viewer, PetscBool skip)
654d71ae5a4SJacob Faibussowitsch {
6555c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
6565c6c1daeSBarry Smith 
6575c6c1daeSBarry Smith   PetscFunctionBegin;
658cc843e7aSLisandro Dalcin   vbinary->skipheader = skip;
6595c6c1daeSBarry Smith   PetscFunctionReturn(0);
6605c6c1daeSBarry Smith }
6615c6c1daeSBarry Smith 
6625c6c1daeSBarry Smith /*@
6635c6c1daeSBarry Smith     PetscViewerBinaryGetSkipHeader - checks whether to write a header with size information on output, or just raw data
6645c6c1daeSBarry Smith 
6655c6c1daeSBarry Smith     Not Collective
6665c6c1daeSBarry Smith 
6675c6c1daeSBarry Smith     Input Parameter:
668811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
6695c6c1daeSBarry Smith 
6705c6c1daeSBarry Smith     Output Parameter:
671811af0c4SBarry Smith .   skip - `PETSC_TRUE` means do not write header
6725c6c1daeSBarry Smith 
6735c6c1daeSBarry Smith     Level: advanced
6745c6c1daeSBarry Smith 
67595452b02SPatrick Sanan     Notes:
67695452b02SPatrick Sanan     This must be called after PetscViewerSetType()
6775c6c1daeSBarry Smith 
678811af0c4SBarry Smith     Returns `PETSC_FALSE` for `PETSCSOCKETVIEWER`, you cannot skip the header for it.
6795c6c1daeSBarry Smith 
680811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
681db781477SPatrick Sanan           `PetscViewerBinarySetSkipHeader()`
6825c6c1daeSBarry Smith @*/
683d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetSkipHeader(PetscViewer viewer, PetscBool *skip)
684d71ae5a4SJacob Faibussowitsch {
6855c6c1daeSBarry Smith   PetscFunctionBegin;
686cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
687cc843e7aSLisandro Dalcin   PetscValidBoolPointer(skip, 2);
688cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetSkipHeader_C", (PetscViewer, PetscBool *), (viewer, skip));
6895c6c1daeSBarry Smith   PetscFunctionReturn(0);
6905c6c1daeSBarry Smith }
6915c6c1daeSBarry Smith 
692d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetSkipHeader_Binary(PetscViewer viewer, PetscBool *skip)
693d71ae5a4SJacob Faibussowitsch {
694f3b211e4SSatish Balay   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
6955c6c1daeSBarry Smith 
6965c6c1daeSBarry Smith   PetscFunctionBegin;
697cc843e7aSLisandro Dalcin   *skip = vbinary->skipheader;
6985c6c1daeSBarry Smith   PetscFunctionReturn(0);
6995c6c1daeSBarry Smith }
7005c6c1daeSBarry Smith 
7015c6c1daeSBarry Smith /*@C
7025c6c1daeSBarry Smith     PetscViewerBinaryGetInfoPointer - Extracts the file pointer for the ASCII
7035c6c1daeSBarry Smith           info file associated with a binary file.
7045c6c1daeSBarry Smith 
7055c6c1daeSBarry Smith     Not Collective
7065c6c1daeSBarry Smith 
7075c6c1daeSBarry Smith     Input Parameter:
708811af0c4SBarry Smith .   viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
7095c6c1daeSBarry Smith 
7105c6c1daeSBarry Smith     Output Parameter:
7110298fd71SBarry Smith .   file - file pointer  Always returns NULL if not a binary viewer
7125c6c1daeSBarry Smith 
7135c6c1daeSBarry Smith     Level: advanced
7145c6c1daeSBarry Smith 
715811af0c4SBarry Smith     Note:
716811af0c4SBarry Smith       For writable binary `PetscViewer`s, the file pointer will only be valid for the
717811af0c4SBarry Smith     first processor in the communicator that shares the `PetscViewer`.
7185c6c1daeSBarry Smith 
7195c6c1daeSBarry Smith     Fortran Note:
7205c6c1daeSBarry Smith     This routine is not supported in Fortran.
7215c6c1daeSBarry Smith 
722811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`
7235c6c1daeSBarry Smith @*/
724d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetInfoPointer(PetscViewer viewer, FILE **file)
725d71ae5a4SJacob Faibussowitsch {
7265c6c1daeSBarry Smith   PetscFunctionBegin;
727cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
728cc843e7aSLisandro Dalcin   PetscValidPointer(file, 2);
7290298fd71SBarry Smith   *file = NULL;
730cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinaryGetInfoPointer_C", (PetscViewer, FILE **), (viewer, file));
7315c6c1daeSBarry Smith   PetscFunctionReturn(0);
7325c6c1daeSBarry Smith }
7335c6c1daeSBarry Smith 
734d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetInfoPointer_Binary(PetscViewer viewer, FILE **file)
735d71ae5a4SJacob Faibussowitsch {
736cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
7375c6c1daeSBarry Smith 
7385c6c1daeSBarry Smith   PetscFunctionBegin;
7399566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
740cc843e7aSLisandro Dalcin   *file = vbinary->fdes_info;
741cc843e7aSLisandro Dalcin   if (viewer->format == PETSC_VIEWER_BINARY_MATLAB && !vbinary->matlabheaderwritten) {
7425c6c1daeSBarry Smith     if (vbinary->fdes_info) {
743cc843e7aSLisandro Dalcin       FILE *info = vbinary->fdes_info;
7449566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- begin code written by PetscViewerBinary for MATLAB format ---#\n"));
7459566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#$$ Set.filename = '%s';\n", vbinary->filename));
7469566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#$$ fd = PetscOpenFile(Set.filename);\n"));
7479566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n"));
748cc843e7aSLisandro Dalcin     }
749cc843e7aSLisandro Dalcin     vbinary->matlabheaderwritten = PETSC_TRUE;
7505c6c1daeSBarry Smith   }
7515c6c1daeSBarry Smith   PetscFunctionReturn(0);
7525c6c1daeSBarry Smith }
7535c6c1daeSBarry Smith 
754e0385b85SDave May #if defined(PETSC_HAVE_MPIIO)
755d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_BinaryMPIIO(PetscViewer v)
756d71ae5a4SJacob Faibussowitsch {
757e0385b85SDave May   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
758e0385b85SDave May 
759e0385b85SDave May   PetscFunctionBegin;
76048a46eb9SPierre Jolivet   if (vbinary->mfdes != MPI_FILE_NULL) PetscCallMPI(MPI_File_close(&vbinary->mfdes));
76148a46eb9SPierre Jolivet   if (vbinary->mfsub != MPI_FILE_NULL) PetscCallMPI(MPI_File_close(&vbinary->mfsub));
762cc843e7aSLisandro Dalcin   vbinary->moff = 0;
763e0385b85SDave May   PetscFunctionReturn(0);
764e0385b85SDave May }
765e0385b85SDave May #endif
766e0385b85SDave May 
767d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_BinarySTDIO(PetscViewer v)
768d71ae5a4SJacob Faibussowitsch {
769cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
770cc843e7aSLisandro Dalcin 
771cc843e7aSLisandro Dalcin   PetscFunctionBegin;
772cc843e7aSLisandro Dalcin   if (vbinary->fdes != -1) {
7739566063dSJacob Faibussowitsch     PetscCall(PetscBinaryClose(vbinary->fdes));
774cc843e7aSLisandro Dalcin     vbinary->fdes = -1;
775cc843e7aSLisandro Dalcin     if (vbinary->storecompressed) {
776cc843e7aSLisandro Dalcin       char        cmd[8 + PETSC_MAX_PATH_LEN], out[64 + PETSC_MAX_PATH_LEN] = "";
777cc843e7aSLisandro Dalcin       const char *gzfilename = vbinary->ogzfilename ? vbinary->ogzfilename : vbinary->filename;
778cc843e7aSLisandro Dalcin       /* compress the file */
7799566063dSJacob Faibussowitsch       PetscCall(PetscStrncpy(cmd, "gzip -f ", sizeof(cmd)));
7809566063dSJacob Faibussowitsch       PetscCall(PetscStrlcat(cmd, gzfilename, sizeof(cmd)));
781cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_POPEN)
782cc843e7aSLisandro Dalcin       {
783cc843e7aSLisandro Dalcin         FILE *fp;
7849566063dSJacob Faibussowitsch         PetscCall(PetscPOpen(PETSC_COMM_SELF, NULL, cmd, "r", &fp));
785cc73adaaSBarry Smith         PetscCheck(!fgets(out, (int)(sizeof(out) - 1), fp), PETSC_COMM_SELF, PETSC_ERR_LIB, "Error from command %s\n%s", cmd, out);
7869566063dSJacob Faibussowitsch         PetscCall(PetscPClose(PETSC_COMM_SELF, fp));
787cc843e7aSLisandro Dalcin       }
788cc843e7aSLisandro Dalcin #endif
789cc843e7aSLisandro Dalcin     }
790cc843e7aSLisandro Dalcin   }
7919566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary->ogzfilename));
792cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
793cc843e7aSLisandro Dalcin }
794cc843e7aSLisandro Dalcin 
795d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_BinaryInfo(PetscViewer v)
796d71ae5a4SJacob Faibussowitsch {
797cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
798cc843e7aSLisandro Dalcin 
799cc843e7aSLisandro Dalcin   PetscFunctionBegin;
800cc843e7aSLisandro Dalcin   if (v->format == PETSC_VIEWER_BINARY_MATLAB && vbinary->matlabheaderwritten) {
801cc843e7aSLisandro Dalcin     if (vbinary->fdes_info) {
802cc843e7aSLisandro Dalcin       FILE *info = vbinary->fdes_info;
8039566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- begin code written by PetscViewerBinary for MATLAB format ---#\n"));
8049566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#$$ close(fd);\n"));
8059566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n"));
806cc843e7aSLisandro Dalcin     }
807cc843e7aSLisandro Dalcin   }
808cc843e7aSLisandro Dalcin   if (vbinary->fdes_info) {
809cc843e7aSLisandro Dalcin     FILE *info         = vbinary->fdes_info;
810cc843e7aSLisandro Dalcin     vbinary->fdes_info = NULL;
811cc73adaaSBarry Smith     PetscCheck(!fclose(info), PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file");
812cc843e7aSLisandro Dalcin   }
813cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
814cc843e7aSLisandro Dalcin }
815cc843e7aSLisandro Dalcin 
816d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_Binary(PetscViewer v)
817d71ae5a4SJacob Faibussowitsch {
818cc843e7aSLisandro Dalcin   PetscFunctionBegin;
819cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
8209566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_BinaryMPIIO(v));
821cc843e7aSLisandro Dalcin #endif
8229566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_BinarySTDIO(v));
8239566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_BinaryInfo(v));
824cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
825cc843e7aSLisandro Dalcin }
826cc843e7aSLisandro Dalcin 
827d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerDestroy_Binary(PetscViewer v)
828d71ae5a4SJacob Faibussowitsch {
8295c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
8305c6c1daeSBarry Smith 
8315c6c1daeSBarry Smith   PetscFunctionBegin;
8329566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_Binary(v));
8339566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary->filename));
8349566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary));
8352e956fe4SStefano Zampini   PetscCall(PetscViewerBinaryClearFunctionList(v));
836e0385b85SDave May   PetscFunctionReturn(0);
837e0385b85SDave May }
8385c6c1daeSBarry Smith 
8395c6c1daeSBarry Smith /*@C
8405c6c1daeSBarry Smith    PetscViewerBinaryOpen - Opens a file for binary input/output.
8415c6c1daeSBarry Smith 
842d083f849SBarry Smith    Collective
8435c6c1daeSBarry Smith 
8445c6c1daeSBarry Smith    Input Parameters:
8455c6c1daeSBarry Smith +  comm - MPI communicator
8465c6c1daeSBarry Smith .  name - name of file
847cc843e7aSLisandro Dalcin -  mode - open mode of file
8485c6c1daeSBarry Smith $    FILE_MODE_WRITE - create new file for binary output
8495c6c1daeSBarry Smith $    FILE_MODE_READ - open existing file for binary input
8505c6c1daeSBarry Smith $    FILE_MODE_APPEND - open existing file for binary output
8515c6c1daeSBarry Smith 
8525c6c1daeSBarry Smith    Output Parameter:
853cc843e7aSLisandro Dalcin .  viewer - PetscViewer for binary input/output to use with the specified file
8545c6c1daeSBarry Smith 
8555c6c1daeSBarry Smith     Options Database Keys:
856811af0c4SBarry Smith +    -viewer_binary_filename <name> - name of file to use
857811af0c4SBarry Smith .    -viewer_binary_skip_info - true to skip opening an info file
858811af0c4SBarry Smith .    -viewer_binary_skip_options - true to not use options database while creating viewer
859811af0c4SBarry Smith .    -viewer_binary_skip_header - true to skip output object headers to the file
860811af0c4SBarry Smith -    -viewer_binary_mpiio - true to use MPI-IO for input and output to the file (more scalable for large problems)
8615c6c1daeSBarry Smith 
8625c6c1daeSBarry Smith    Level: beginner
8635c6c1daeSBarry Smith 
8645c6c1daeSBarry Smith    Note:
865811af0c4SBarry Smith    This `PetscViewer` should be destroyed with `PetscViewerDestroy()`.
8665c6c1daeSBarry Smith 
8675c6c1daeSBarry Smith     For reading files, the filename may begin with ftp:// or http:// and/or
8685c6c1daeSBarry Smith     end with .gz; in this case file is brought over and uncompressed.
8695c6c1daeSBarry Smith 
8705c6c1daeSBarry Smith     For creating files, if the file name ends with .gz it is automatically
8715c6c1daeSBarry Smith     compressed when closed.
8725c6c1daeSBarry Smith 
873811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
874db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
875db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`, `PetscViewerBinarySetUseMPIIO()`,
876db781477SPatrick Sanan           `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
8775c6c1daeSBarry Smith @*/
878d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryOpen(MPI_Comm comm, const char name[], PetscFileMode mode, PetscViewer *viewer)
879d71ae5a4SJacob Faibussowitsch {
8805c6c1daeSBarry Smith   PetscFunctionBegin;
8819566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(comm, viewer));
8829566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(*viewer, PETSCVIEWERBINARY));
8839566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(*viewer, mode));
8849566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetName(*viewer, name));
8859566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions(*viewer));
8865c6c1daeSBarry Smith   PetscFunctionReturn(0);
8875c6c1daeSBarry Smith }
8885c6c1daeSBarry Smith 
8895c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
890d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryWriteReadMPIIO(PetscViewer viewer, void *data, PetscInt num, PetscInt *count, PetscDataType dtype, PetscBool write)
891d71ae5a4SJacob Faibussowitsch {
89222a8f86cSLisandro Dalcin   MPI_Comm            comm    = PetscObjectComm((PetscObject)viewer);
8935c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
89422a8f86cSLisandro Dalcin   MPI_File            mfdes   = vbinary->mfdes;
8955c6c1daeSBarry Smith   MPI_Datatype        mdtype;
89669e10bbaSLisandro Dalcin   PetscMPIInt         rank, cnt;
8975c6c1daeSBarry Smith   MPI_Status          status;
8985c6c1daeSBarry Smith   MPI_Aint            ul, dsize;
8995c6c1daeSBarry Smith 
9005c6c1daeSBarry Smith   PetscFunctionBegin;
9019566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
9029566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(num, &cnt));
9039566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToMPIDataType(dtype, &mdtype));
9045c6c1daeSBarry Smith   if (write) {
90548a46eb9SPierre Jolivet     if (rank == 0) PetscCall(MPIU_File_write_at(mfdes, vbinary->moff, data, cnt, mdtype, &status));
9065c6c1daeSBarry Smith   } else {
907dd400576SPatrick Sanan     if (rank == 0) {
9089566063dSJacob Faibussowitsch       PetscCall(MPIU_File_read_at(mfdes, vbinary->moff, data, cnt, mdtype, &status));
9099566063dSJacob Faibussowitsch       if (cnt > 0) PetscCallMPI(MPI_Get_count(&status, mdtype, &cnt));
9105c6c1daeSBarry Smith     }
9119566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&cnt, 1, MPI_INT, 0, comm));
9129566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(data, cnt, mdtype, 0, comm));
91369e10bbaSLisandro Dalcin   }
9149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_get_extent(mdtype, &ul, &dsize));
9155c6c1daeSBarry Smith   vbinary->moff += dsize * cnt;
9169860990eSLisandro Dalcin   if (count) *count = cnt;
9175c6c1daeSBarry Smith   PetscFunctionReturn(0);
9185c6c1daeSBarry Smith }
9195c6c1daeSBarry Smith #endif
9205c6c1daeSBarry Smith 
9215c6c1daeSBarry Smith /*@C
9225c6c1daeSBarry Smith    PetscViewerBinaryRead - Reads from a binary file, all processors get the same result
9235c6c1daeSBarry Smith 
924d083f849SBarry Smith    Collective
9255c6c1daeSBarry Smith 
9265c6c1daeSBarry Smith    Input Parameters:
9275c6c1daeSBarry Smith +  viewer - the binary viewer
928907376e6SBarry Smith .  data - location of the data to be written
929060da220SMatthew G. Knepley .  num - number of items of data to read
930907376e6SBarry Smith -  dtype - type of data to read
9315c6c1daeSBarry Smith 
932f8e4bde8SMatthew G. Knepley    Output Parameters:
9335972f5f3SLisandro Dalcin .  count - number of items of data actually read, or NULL.
934f8e4bde8SMatthew G. Knepley 
9355c6c1daeSBarry Smith    Level: beginner
9365c6c1daeSBarry Smith 
937811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
938db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
939db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
9405c6c1daeSBarry Smith @*/
941d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryRead(PetscViewer viewer, void *data, PetscInt num, PetscInt *count, PetscDataType dtype)
942d71ae5a4SJacob Faibussowitsch {
94322a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
9445c6c1daeSBarry Smith 
94522a8f86cSLisandro Dalcin   PetscFunctionBegin;
94622a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
94722a8f86cSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, num, 3);
9489566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
94922a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
9505c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
951bc196f7cSDave May   if (vbinary->usempiio) {
9529566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWriteReadMPIIO(viewer, data, num, count, dtype, PETSC_FALSE));
9535c6c1daeSBarry Smith   } else {
9545c6c1daeSBarry Smith #endif
9559566063dSJacob Faibussowitsch     PetscCall(PetscBinarySynchronizedRead(PetscObjectComm((PetscObject)viewer), vbinary->fdes, data, num, count, dtype));
9565c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
9575c6c1daeSBarry Smith   }
9585c6c1daeSBarry Smith #endif
9595c6c1daeSBarry Smith   PetscFunctionReturn(0);
9605c6c1daeSBarry Smith }
9615c6c1daeSBarry Smith 
9625c6c1daeSBarry Smith /*@C
963811af0c4SBarry Smith    PetscViewerBinaryWrite - writes to a binary file, only from the first MPI rank
9645c6c1daeSBarry Smith 
965d083f849SBarry Smith    Collective
9665c6c1daeSBarry Smith 
9675c6c1daeSBarry Smith    Input Parameters:
9685c6c1daeSBarry Smith +  viewer - the binary viewer
9695c6c1daeSBarry Smith .  data - location of data
9705824c9d0SPatrick Sanan .  count - number of items of data to write
971f253e43cSLisandro Dalcin -  dtype - type of data to write
9725c6c1daeSBarry Smith 
9735c6c1daeSBarry Smith    Level: beginner
9745c6c1daeSBarry Smith 
975811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
976db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`, `PetscDataType`
977db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
9785c6c1daeSBarry Smith @*/
979d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryWrite(PetscViewer viewer, const void *data, PetscInt count, PetscDataType dtype)
980d71ae5a4SJacob Faibussowitsch {
98122a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
9825c6c1daeSBarry Smith 
9835c6c1daeSBarry Smith   PetscFunctionBegin;
98422a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
98522a8f86cSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, count, 3);
9869566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
98722a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
9885c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
989bc196f7cSDave May   if (vbinary->usempiio) {
9909566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWriteReadMPIIO(viewer, (void *)data, count, NULL, dtype, PETSC_TRUE));
9915c6c1daeSBarry Smith   } else {
9925c6c1daeSBarry Smith #endif
9939566063dSJacob Faibussowitsch     PetscCall(PetscBinarySynchronizedWrite(PetscObjectComm((PetscObject)viewer), vbinary->fdes, data, count, dtype));
9945c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
9955c6c1daeSBarry Smith   }
9965c6c1daeSBarry Smith #endif
9975c6c1daeSBarry Smith   PetscFunctionReturn(0);
9985c6c1daeSBarry Smith }
9995c6c1daeSBarry Smith 
1000d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryWriteReadAll(PetscViewer viewer, PetscBool write, void *data, PetscInt count, PetscInt start, PetscInt total, PetscDataType dtype)
1001d71ae5a4SJacob Faibussowitsch {
10025972f5f3SLisandro Dalcin   MPI_Comm              comm = PetscObjectComm((PetscObject)viewer);
10035972f5f3SLisandro Dalcin   PetscMPIInt           size, rank;
10045972f5f3SLisandro Dalcin   MPI_Datatype          mdtype;
1005ec4bef21SJose E. Roman   PETSC_UNUSED MPI_Aint lb;
1006ec4bef21SJose E. Roman   MPI_Aint              dsize;
10075972f5f3SLisandro Dalcin   PetscBool             useMPIIO;
10085972f5f3SLisandro Dalcin 
10095972f5f3SLisandro Dalcin   PetscFunctionBegin;
10105972f5f3SLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
1011064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveBool(viewer, ((start >= 0) || (start == PETSC_DETERMINE)), 5);
1012064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveBool(viewer, ((total >= 0) || (total == PETSC_DETERMINE)), 6);
1013064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveInt(viewer, total, 6);
10149566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
10155972f5f3SLisandro Dalcin 
10169566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToMPIDataType(dtype, &mdtype));
10179566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_get_extent(mdtype, &lb, &dsize));
10189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
10199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
10205972f5f3SLisandro Dalcin 
10219566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetUseMPIIO(viewer, &useMPIIO));
10225972f5f3SLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
10235972f5f3SLisandro Dalcin   if (useMPIIO) {
10245972f5f3SLisandro Dalcin     MPI_File    mfdes;
10255972f5f3SLisandro Dalcin     MPI_Offset  off;
10265972f5f3SLisandro Dalcin     PetscMPIInt cnt;
10275972f5f3SLisandro Dalcin 
10285972f5f3SLisandro Dalcin     if (start == PETSC_DETERMINE) {
10299566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Scan(&count, &start, 1, MPIU_INT, MPI_SUM, comm));
10305972f5f3SLisandro Dalcin       start -= count;
10315972f5f3SLisandro Dalcin     }
10325972f5f3SLisandro Dalcin     if (total == PETSC_DETERMINE) {
10335972f5f3SLisandro Dalcin       total = start + count;
10349566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Bcast(&total, 1, MPIU_INT, size - 1, comm));
10355972f5f3SLisandro Dalcin     }
10369566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(count, &cnt));
10379566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryGetMPIIODescriptor(viewer, &mfdes));
10389566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryGetMPIIOOffset(viewer, &off));
10395972f5f3SLisandro Dalcin     off += (MPI_Offset)(start * dsize);
10405972f5f3SLisandro Dalcin     if (write) {
10419566063dSJacob Faibussowitsch       PetscCall(MPIU_File_write_at_all(mfdes, off, data, cnt, mdtype, MPI_STATUS_IGNORE));
10425972f5f3SLisandro Dalcin     } else {
10439566063dSJacob Faibussowitsch       PetscCall(MPIU_File_read_at_all(mfdes, off, data, cnt, mdtype, MPI_STATUS_IGNORE));
10445972f5f3SLisandro Dalcin     }
10455972f5f3SLisandro Dalcin     off = (MPI_Offset)(total * dsize);
10469566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryAddMPIIOOffset(viewer, off));
10475972f5f3SLisandro Dalcin     PetscFunctionReturn(0);
10485972f5f3SLisandro Dalcin   }
10495972f5f3SLisandro Dalcin #endif
10505972f5f3SLisandro Dalcin   {
10515972f5f3SLisandro Dalcin     int         fdes;
10525972f5f3SLisandro Dalcin     char       *workbuf = NULL;
1053dd400576SPatrick Sanan     PetscInt    tcount = rank == 0 ? 0 : count, maxcount = 0, message_count, flowcontrolcount;
10545972f5f3SLisandro Dalcin     PetscMPIInt tag, cnt, maxcnt, scnt = 0, rcnt = 0, j;
10555972f5f3SLisandro Dalcin     MPI_Status  status;
10565972f5f3SLisandro Dalcin 
10579566063dSJacob Faibussowitsch     PetscCall(PetscCommGetNewTag(comm, &tag));
10589566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Reduce(&tcount, &maxcount, 1, MPIU_INT, MPI_MAX, 0, comm));
10599566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(maxcount, &maxcnt));
10609566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(count, &cnt));
10615972f5f3SLisandro Dalcin 
10629566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryGetDescriptor(viewer, &fdes));
10639566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlowControlStart(viewer, &message_count, &flowcontrolcount));
1064dd400576SPatrick Sanan     if (rank == 0) {
10659566063dSJacob Faibussowitsch       PetscCall(PetscMalloc(maxcnt * dsize, &workbuf));
10665972f5f3SLisandro Dalcin       if (write) {
10679566063dSJacob Faibussowitsch         PetscCall(PetscBinaryWrite(fdes, data, cnt, dtype));
10685972f5f3SLisandro Dalcin       } else {
10699566063dSJacob Faibussowitsch         PetscCall(PetscBinaryRead(fdes, data, cnt, NULL, dtype));
10705972f5f3SLisandro Dalcin       }
10715972f5f3SLisandro Dalcin       for (j = 1; j < size; j++) {
10729566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlowControlStepMain(viewer, j, &message_count, flowcontrolcount));
10735972f5f3SLisandro Dalcin         if (write) {
10749566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Recv(workbuf, maxcnt, mdtype, j, tag, comm, &status));
10759566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Get_count(&status, mdtype, &rcnt));
10769566063dSJacob Faibussowitsch           PetscCall(PetscBinaryWrite(fdes, workbuf, rcnt, dtype));
10775972f5f3SLisandro Dalcin         } else {
10789566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Recv(&scnt, 1, MPI_INT, j, tag, comm, MPI_STATUS_IGNORE));
10799566063dSJacob Faibussowitsch           PetscCall(PetscBinaryRead(fdes, workbuf, scnt, NULL, dtype));
10809566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Send(workbuf, scnt, mdtype, j, tag, comm));
10815972f5f3SLisandro Dalcin         }
10825972f5f3SLisandro Dalcin       }
10839566063dSJacob Faibussowitsch       PetscCall(PetscFree(workbuf));
10849566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlowControlEndMain(viewer, &message_count));
10855972f5f3SLisandro Dalcin     } else {
10869566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlowControlStepWorker(viewer, rank, &message_count));
10875972f5f3SLisandro Dalcin       if (write) {
10889566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(data, cnt, mdtype, 0, tag, comm));
10895972f5f3SLisandro Dalcin       } else {
10909566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(&cnt, 1, MPI_INT, 0, tag, comm));
10919566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Recv(data, cnt, mdtype, 0, tag, comm, MPI_STATUS_IGNORE));
10925972f5f3SLisandro Dalcin       }
10939566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlowControlEndWorker(viewer, &message_count));
10945972f5f3SLisandro Dalcin     }
10955972f5f3SLisandro Dalcin   }
10965972f5f3SLisandro Dalcin   PetscFunctionReturn(0);
10975972f5f3SLisandro Dalcin }
10985972f5f3SLisandro Dalcin 
10995972f5f3SLisandro Dalcin /*@C
1100811af0c4SBarry Smith    PetscViewerBinaryReadAll - reads from a binary file from all MPI ranks, each rank receives its own portion of the data
11015972f5f3SLisandro Dalcin 
11025972f5f3SLisandro Dalcin    Collective
11035972f5f3SLisandro Dalcin 
11045972f5f3SLisandro Dalcin    Input Parameters:
11055972f5f3SLisandro Dalcin +  viewer - the binary viewer
11065972f5f3SLisandro Dalcin .  data - location of data
11075972f5f3SLisandro Dalcin .  count - local number of items of data to read
1108811af0c4SBarry Smith .  start - local start, can be `PETSC_DETERMINE`
1109811af0c4SBarry Smith .  total - global number of items of data to read, can be `PETSC_DETERMINE`
11105972f5f3SLisandro Dalcin -  dtype - type of data to read
11115972f5f3SLisandro Dalcin 
11125972f5f3SLisandro Dalcin    Level: advanced
11135972f5f3SLisandro Dalcin 
1114811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryRead()`, `PetscViewerBinaryWriteAll()`
11155972f5f3SLisandro Dalcin @*/
1116d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryReadAll(PetscViewer viewer, void *data, PetscInt count, PetscInt start, PetscInt total, PetscDataType dtype)
1117d71ae5a4SJacob Faibussowitsch {
11185972f5f3SLisandro Dalcin   PetscFunctionBegin;
11199566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWriteReadAll(viewer, PETSC_FALSE, data, count, start, total, dtype));
11205972f5f3SLisandro Dalcin   PetscFunctionReturn(0);
11215972f5f3SLisandro Dalcin }
11225972f5f3SLisandro Dalcin 
11235972f5f3SLisandro Dalcin /*@C
1124811af0c4SBarry Smith    PetscViewerBinaryWriteAll - writes to a binary file from all MPI ranks, each rank writes its own portion of the data
11255972f5f3SLisandro Dalcin 
11265972f5f3SLisandro Dalcin    Collective
11275972f5f3SLisandro Dalcin 
11285972f5f3SLisandro Dalcin    Input Parameters:
11295972f5f3SLisandro Dalcin +  viewer - the binary viewer
11305972f5f3SLisandro Dalcin .  data - location of data
11315972f5f3SLisandro Dalcin .  count - local number of items of data to write
1132811af0c4SBarry Smith .  start - local start, can be `PETSC_DETERMINE`
1133811af0c4SBarry Smith .  total - global number of items of data to write, can be `PETSC_DETERMINE`
11345972f5f3SLisandro Dalcin -  dtype - type of data to write
11355972f5f3SLisandro Dalcin 
11365972f5f3SLisandro Dalcin    Level: advanced
11375972f5f3SLisandro Dalcin 
1138811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryWriteAll()`, `PetscViewerBinaryReadAll()`
11395972f5f3SLisandro Dalcin @*/
1140d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryWriteAll(PetscViewer viewer, const void *data, PetscInt count, PetscInt start, PetscInt total, PetscDataType dtype)
1141d71ae5a4SJacob Faibussowitsch {
11425972f5f3SLisandro Dalcin   PetscFunctionBegin;
11439566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWriteReadAll(viewer, PETSC_TRUE, (void *)data, count, start, total, dtype));
11445972f5f3SLisandro Dalcin   PetscFunctionReturn(0);
11455972f5f3SLisandro Dalcin }
11465972f5f3SLisandro Dalcin 
11475c6c1daeSBarry Smith /*@C
1148811af0c4SBarry Smith    PetscViewerBinaryWriteStringArray - writes to a binary file, only from the first MPI rank, an array of strings
11495c6c1daeSBarry Smith 
1150d083f849SBarry Smith    Collective
11515c6c1daeSBarry Smith 
11525c6c1daeSBarry Smith    Input Parameters:
11535c6c1daeSBarry Smith +  viewer - the binary viewer
11545c6c1daeSBarry Smith -  data - location of the array of strings
11555c6c1daeSBarry Smith 
11565c6c1daeSBarry Smith    Level: intermediate
11575c6c1daeSBarry Smith 
1158811af0c4SBarry Smith     Note:
1159811af0c4SBarry Smith     The array of strings must be null terminated
11605c6c1daeSBarry Smith 
1161811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
1162db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
1163db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
11645c6c1daeSBarry Smith @*/
1165d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryWriteStringArray(PetscViewer viewer, const char *const *data)
1166d71ae5a4SJacob Faibussowitsch {
11675c6c1daeSBarry Smith   PetscInt i, n = 0, *sizes;
1168cc843e7aSLisandro Dalcin   size_t   len;
11695c6c1daeSBarry Smith 
117022a8f86cSLisandro Dalcin   PetscFunctionBegin;
11719566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
11725c6c1daeSBarry Smith   /* count number of strings */
11739371c9d4SSatish Balay   while (data[n++])
11749371c9d4SSatish Balay     ;
11755c6c1daeSBarry Smith   n--;
11769566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n + 1, &sizes));
11775c6c1daeSBarry Smith   sizes[0] = n;
11785c6c1daeSBarry Smith   for (i = 0; i < n; i++) {
11799566063dSJacob Faibussowitsch     PetscCall(PetscStrlen(data[i], &len));
1180cc843e7aSLisandro Dalcin     sizes[i + 1] = (PetscInt)len + 1; /* size includes space for the null terminator */
11815c6c1daeSBarry Smith   }
11829566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWrite(viewer, sizes, n + 1, PETSC_INT));
118348a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall(PetscViewerBinaryWrite(viewer, (void *)data[i], sizes[i + 1], PETSC_CHAR));
11849566063dSJacob Faibussowitsch   PetscCall(PetscFree(sizes));
11855c6c1daeSBarry Smith   PetscFunctionReturn(0);
11865c6c1daeSBarry Smith }
11875c6c1daeSBarry Smith 
11885c6c1daeSBarry Smith /*@C
1189811af0c4SBarry Smith    PetscViewerBinaryReadStringArray - reads a binary file an array of strings to all MPI ranks
11905c6c1daeSBarry Smith 
1191d083f849SBarry Smith    Collective
11925c6c1daeSBarry Smith 
11935c6c1daeSBarry Smith    Input Parameter:
11945c6c1daeSBarry Smith .  viewer - the binary viewer
11955c6c1daeSBarry Smith 
11965c6c1daeSBarry Smith    Output Parameter:
11975c6c1daeSBarry Smith .  data - location of the array of strings
11985c6c1daeSBarry Smith 
11995c6c1daeSBarry Smith    Level: intermediate
12005c6c1daeSBarry Smith 
1201811af0c4SBarry Smith     Note:
1202811af0c4SBarry Smith     The array of strings must null terminated
12035c6c1daeSBarry Smith 
1204811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
1205db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
1206db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
12075c6c1daeSBarry Smith @*/
1208d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryReadStringArray(PetscViewer viewer, char ***data)
1209d71ae5a4SJacob Faibussowitsch {
1210060da220SMatthew G. Knepley   PetscInt i, n, *sizes, N = 0;
12115c6c1daeSBarry Smith 
121222a8f86cSLisandro Dalcin   PetscFunctionBegin;
12139566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
12145c6c1daeSBarry Smith   /* count number of strings */
12159566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &n, 1, NULL, PETSC_INT));
12169566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &sizes));
12179566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, sizes, n, NULL, PETSC_INT));
1218a297a907SKarl Rupp   for (i = 0; i < n; i++) N += sizes[i];
12199566063dSJacob Faibussowitsch   PetscCall(PetscMalloc((n + 1) * sizeof(char *) + N * sizeof(char), data));
12205c6c1daeSBarry Smith   (*data)[0] = (char *)((*data) + n + 1);
1221a297a907SKarl Rupp   for (i = 1; i < n; i++) (*data)[i] = (*data)[i - 1] + sizes[i - 1];
12229566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, (*data)[0], N, NULL, PETSC_CHAR));
122302c9f0b5SLisandro Dalcin   (*data)[n] = NULL;
12249566063dSJacob Faibussowitsch   PetscCall(PetscFree(sizes));
12255c6c1daeSBarry Smith   PetscFunctionReturn(0);
12265c6c1daeSBarry Smith }
12275c6c1daeSBarry Smith 
1228cc843e7aSLisandro Dalcin /*@C
1229cc843e7aSLisandro Dalcin      PetscViewerFileSetMode - Sets the open mode of file
1230cc843e7aSLisandro Dalcin 
1231811af0c4SBarry Smith     Logically Collective on viewer
1232cc843e7aSLisandro Dalcin 
1233cc843e7aSLisandro Dalcin   Input Parameters:
1234811af0c4SBarry Smith +  viewer - the `PetscViewer`; must be a a `PETSCVIEWERBINARY`, `PETSCVIEWERMATLAB`, `PETSCVIEWERHDF5`, or `PETSCVIEWERASCII`  `PetscViewer`
1235cc843e7aSLisandro Dalcin -  mode - open mode of file
1236cc843e7aSLisandro Dalcin $    FILE_MODE_WRITE - create new file for output
1237cc843e7aSLisandro Dalcin $    FILE_MODE_READ - open existing file for input
1238cc843e7aSLisandro Dalcin $    FILE_MODE_APPEND - open existing file for output
1239cc843e7aSLisandro Dalcin 
1240cc843e7aSLisandro Dalcin   Level: advanced
1241cc843e7aSLisandro Dalcin 
1242db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
1243cc843e7aSLisandro Dalcin @*/
1244d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerFileSetMode(PetscViewer viewer, PetscFileMode mode)
1245d71ae5a4SJacob Faibussowitsch {
1246cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1247cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1248cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveEnum(viewer, mode, 2);
124908401ef6SPierre Jolivet   PetscCheck(mode != FILE_MODE_UNDEFINED, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Cannot set FILE_MODE_UNDEFINED");
1250f7d195e4SLawrence Mitchell   PetscCheck(mode >= FILE_MODE_UNDEFINED && mode <= FILE_MODE_APPEND_UPDATE, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_OUTOFRANGE, "Invalid file mode %d", (int)mode);
1251cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerFileSetMode_C", (PetscViewer, PetscFileMode), (viewer, mode));
1252cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1253cc843e7aSLisandro Dalcin }
1254cc843e7aSLisandro Dalcin 
1255d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetMode_Binary(PetscViewer viewer, PetscFileMode mode)
1256d71ae5a4SJacob Faibussowitsch {
1257cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1258cc843e7aSLisandro Dalcin 
1259cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1260cc73adaaSBarry Smith   PetscCheck(!viewer->setupcalled || vbinary->filemode == mode, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Cannot change mode to %s after setup", PetscFileModes[mode]);
1261cc843e7aSLisandro Dalcin   vbinary->filemode = mode;
1262cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1263cc843e7aSLisandro Dalcin }
1264cc843e7aSLisandro Dalcin 
1265cc843e7aSLisandro Dalcin /*@C
1266cc843e7aSLisandro Dalcin      PetscViewerFileGetMode - Gets the open mode of file
1267cc843e7aSLisandro Dalcin 
1268cc843e7aSLisandro Dalcin     Not Collective
1269cc843e7aSLisandro Dalcin 
1270cc843e7aSLisandro Dalcin   Input Parameter:
1271811af0c4SBarry Smith .  viewer - the `PetscViewer`; must be a `PETSCVIEWERBINARY`, `PETSCVIEWERMATLAB`, `PETSCVIEWERHDF5`, or `PETSCVIEWERASCII`  `PetscViewer`
1272cc843e7aSLisandro Dalcin 
1273cc843e7aSLisandro Dalcin   Output Parameter:
1274cc843e7aSLisandro Dalcin .  mode - open mode of file
1275cc843e7aSLisandro Dalcin $    FILE_MODE_WRITE - create new file for binary output
1276cc843e7aSLisandro Dalcin $    FILE_MODE_READ - open existing file for binary input
1277cc843e7aSLisandro Dalcin $    FILE_MODE_APPEND - open existing file for binary output
1278cc843e7aSLisandro Dalcin 
1279cc843e7aSLisandro Dalcin   Level: advanced
1280cc843e7aSLisandro Dalcin 
1281db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
1282cc843e7aSLisandro Dalcin @*/
1283d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerFileGetMode(PetscViewer viewer, PetscFileMode *mode)
1284d71ae5a4SJacob Faibussowitsch {
1285cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1286cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1287cc843e7aSLisandro Dalcin   PetscValidPointer(mode, 2);
1288cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerFileGetMode_C", (PetscViewer, PetscFileMode *), (viewer, mode));
1289cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1290cc843e7aSLisandro Dalcin }
1291cc843e7aSLisandro Dalcin 
1292d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetMode_Binary(PetscViewer viewer, PetscFileMode *mode)
1293d71ae5a4SJacob Faibussowitsch {
1294cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1295cc843e7aSLisandro Dalcin 
1296cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1297cc843e7aSLisandro Dalcin   *mode = vbinary->filemode;
1298cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1299cc843e7aSLisandro Dalcin }
1300cc843e7aSLisandro Dalcin 
1301d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetName_Binary(PetscViewer viewer, const char name[])
1302d71ae5a4SJacob Faibussowitsch {
1303cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1304cc843e7aSLisandro Dalcin 
1305cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1306cc843e7aSLisandro Dalcin   if (viewer->setupcalled && vbinary->filename) {
1307cc843e7aSLisandro Dalcin     /* gzip can be run after the file with the previous filename has been closed */
13089566063dSJacob Faibussowitsch     PetscCall(PetscFree(vbinary->ogzfilename));
13099566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(vbinary->filename, &vbinary->ogzfilename));
1310cc843e7aSLisandro Dalcin   }
13119566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary->filename));
13129566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &vbinary->filename));
1313cc843e7aSLisandro Dalcin   viewer->setupcalled = PETSC_FALSE;
1314cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1315cc843e7aSLisandro Dalcin }
1316cc843e7aSLisandro Dalcin 
1317d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetName_Binary(PetscViewer viewer, const char **name)
1318d71ae5a4SJacob Faibussowitsch {
13195c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
13205c6c1daeSBarry Smith 
13215c6c1daeSBarry Smith   PetscFunctionBegin;
13225c6c1daeSBarry Smith   *name = vbinary->filename;
13235c6c1daeSBarry Smith   PetscFunctionReturn(0);
13245c6c1daeSBarry Smith }
13255c6c1daeSBarry Smith 
13265c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
1327d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetUp_BinaryMPIIO(PetscViewer viewer)
1328d71ae5a4SJacob Faibussowitsch {
13295c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1330e8a65b7dSLisandro Dalcin   int                 amode;
13315c6c1daeSBarry Smith 
13325c6c1daeSBarry Smith   PetscFunctionBegin;
1333a297a907SKarl Rupp   vbinary->storecompressed = PETSC_FALSE;
13345c6c1daeSBarry Smith 
1335cc843e7aSLisandro Dalcin   vbinary->moff = 0;
1336cc843e7aSLisandro Dalcin   switch (vbinary->filemode) {
1337d71ae5a4SJacob Faibussowitsch   case FILE_MODE_READ:
1338d71ae5a4SJacob Faibussowitsch     amode = MPI_MODE_RDONLY;
1339d71ae5a4SJacob Faibussowitsch     break;
1340d71ae5a4SJacob Faibussowitsch   case FILE_MODE_WRITE:
1341d71ae5a4SJacob Faibussowitsch     amode = MPI_MODE_WRONLY | MPI_MODE_CREATE;
1342d71ae5a4SJacob Faibussowitsch     break;
1343d71ae5a4SJacob Faibussowitsch   case FILE_MODE_APPEND:
1344d71ae5a4SJacob Faibussowitsch     amode = MPI_MODE_WRONLY | MPI_MODE_CREATE | MPI_MODE_APPEND;
1345d71ae5a4SJacob Faibussowitsch     break;
1346d71ae5a4SJacob Faibussowitsch   case FILE_MODE_UNDEFINED:
1347d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerSetUp()");
1348d71ae5a4SJacob Faibussowitsch   default:
1349d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[vbinary->filemode]);
13505c6c1daeSBarry Smith   }
13519566063dSJacob Faibussowitsch   PetscCallMPI(MPI_File_open(PetscObjectComm((PetscObject)viewer), vbinary->filename, amode, MPI_INFO_NULL, &vbinary->mfdes));
135222a8f86cSLisandro Dalcin   /*
135322a8f86cSLisandro Dalcin       The MPI standard does not have MPI_MODE_TRUNCATE. We emulate this behavior by setting the file size to zero.
135422a8f86cSLisandro Dalcin   */
13559566063dSJacob Faibussowitsch   if (vbinary->filemode == FILE_MODE_WRITE) PetscCallMPI(MPI_File_set_size(vbinary->mfdes, 0));
135622a8f86cSLisandro Dalcin   /*
135722a8f86cSLisandro Dalcin       Initially, all processes view the file as a linear byte stream. Therefore, for files opened with MPI_MODE_APPEND,
135822a8f86cSLisandro Dalcin       MPI_File_get_position[_shared](fh, &offset) returns the absolute byte position at the end of file.
135922a8f86cSLisandro Dalcin       Otherwise, we would need to call MPI_File_get_byte_offset(fh, offset, &byte_offset) to convert
136022a8f86cSLisandro Dalcin       the offset in etype units to an absolute byte position.
136122a8f86cSLisandro Dalcin    */
13629566063dSJacob Faibussowitsch   if (vbinary->filemode == FILE_MODE_APPEND) PetscCallMPI(MPI_File_get_position(vbinary->mfdes, &vbinary->moff));
1363cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1364cc843e7aSLisandro Dalcin }
1365cc843e7aSLisandro Dalcin #endif
13665c6c1daeSBarry Smith 
1367d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetUp_BinarySTDIO(PetscViewer viewer)
1368d71ae5a4SJacob Faibussowitsch {
1369cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1370cc843e7aSLisandro Dalcin   const char         *fname;
1371cc843e7aSLisandro Dalcin   char                bname[PETSC_MAX_PATH_LEN], *gz;
1372cc843e7aSLisandro Dalcin   PetscBool           found;
1373cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
13745c6c1daeSBarry Smith 
1375cc843e7aSLisandro Dalcin   PetscFunctionBegin;
13769566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
13775c6c1daeSBarry Smith 
1378cc843e7aSLisandro Dalcin   /* if file name ends in .gz strip that off and note user wants file compressed */
1379cc843e7aSLisandro Dalcin   vbinary->storecompressed = PETSC_FALSE;
1380cc843e7aSLisandro Dalcin   if (vbinary->filemode == FILE_MODE_WRITE) {
13819566063dSJacob Faibussowitsch     PetscCall(PetscStrstr(vbinary->filename, ".gz", &gz));
13829371c9d4SSatish Balay     if (gz && gz[3] == 0) {
13839371c9d4SSatish Balay       *gz                      = 0;
13849371c9d4SSatish Balay       vbinary->storecompressed = PETSC_TRUE;
13859371c9d4SSatish Balay     }
1386cc843e7aSLisandro Dalcin   }
1387cc843e7aSLisandro Dalcin #if !defined(PETSC_HAVE_POPEN)
138828b400f6SJacob Faibussowitsch   PetscCheck(!vbinary->storecompressed, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP_SYS, "Cannot run gzip on this machine");
1389cc843e7aSLisandro Dalcin #endif
1390cc843e7aSLisandro Dalcin 
1391cc843e7aSLisandro Dalcin   fname = vbinary->filename;
1392cc843e7aSLisandro Dalcin   if (vbinary->filemode == FILE_MODE_READ) { /* possibly get the file from remote site or compressed file */
13939566063dSJacob Faibussowitsch     PetscCall(PetscFileRetrieve(PetscObjectComm((PetscObject)viewer), fname, bname, PETSC_MAX_PATH_LEN, &found));
139428b400f6SJacob Faibussowitsch     PetscCheck(found, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_OPEN, "Cannot locate file: %s", fname);
1395cc843e7aSLisandro Dalcin     fname = bname;
13965c6c1daeSBarry Smith   }
13975c6c1daeSBarry Smith 
1398cc843e7aSLisandro Dalcin   vbinary->fdes = -1;
1399dd400576SPatrick Sanan   if (rank == 0) { /* only first processor opens file*/
1400cc843e7aSLisandro Dalcin     PetscFileMode mode = vbinary->filemode;
1401cc843e7aSLisandro Dalcin     if (mode == FILE_MODE_APPEND) {
1402cc843e7aSLisandro Dalcin       /* check if asked to append to a non-existing file */
14039566063dSJacob Faibussowitsch       PetscCall(PetscTestFile(fname, '\0', &found));
1404cc843e7aSLisandro Dalcin       if (!found) mode = FILE_MODE_WRITE;
1405cc843e7aSLisandro Dalcin     }
14069566063dSJacob Faibussowitsch     PetscCall(PetscBinaryOpen(fname, mode, &vbinary->fdes));
1407cc843e7aSLisandro Dalcin   }
1408cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1409cc843e7aSLisandro Dalcin }
1410cc843e7aSLisandro Dalcin 
1411d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetUp_BinaryInfo(PetscViewer viewer)
1412d71ae5a4SJacob Faibussowitsch {
1413cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1414cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
1415cc843e7aSLisandro Dalcin   PetscBool           found;
1416cc843e7aSLisandro Dalcin 
1417cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1418cc843e7aSLisandro Dalcin   vbinary->fdes_info = NULL;
14199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
1420dd400576SPatrick Sanan   if (!vbinary->skipinfo && (vbinary->filemode == FILE_MODE_READ || rank == 0)) {
1421cc843e7aSLisandro Dalcin     char infoname[PETSC_MAX_PATH_LEN], iname[PETSC_MAX_PATH_LEN], *gz;
1422cc843e7aSLisandro Dalcin 
14239566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(infoname, vbinary->filename, sizeof(infoname)));
1424cc843e7aSLisandro Dalcin     /* remove .gz if it ends file name */
14259566063dSJacob Faibussowitsch     PetscCall(PetscStrstr(infoname, ".gz", &gz));
1426cc843e7aSLisandro Dalcin     if (gz && gz[3] == 0) *gz = 0;
1427cc843e7aSLisandro Dalcin 
14289566063dSJacob Faibussowitsch     PetscCall(PetscStrlcat(infoname, ".info", sizeof(infoname)));
1429cc843e7aSLisandro Dalcin     if (vbinary->filemode == FILE_MODE_READ) {
14309566063dSJacob Faibussowitsch       PetscCall(PetscFixFilename(infoname, iname));
14319566063dSJacob Faibussowitsch       PetscCall(PetscFileRetrieve(PetscObjectComm((PetscObject)viewer), iname, infoname, PETSC_MAX_PATH_LEN, &found));
14329566063dSJacob Faibussowitsch       if (found) PetscCall(PetscOptionsInsertFile(PetscObjectComm((PetscObject)viewer), ((PetscObject)viewer)->options, infoname, PETSC_FALSE));
1433dd400576SPatrick Sanan     } else if (rank == 0) { /* write or append */
1434cc843e7aSLisandro Dalcin       const char *omode  = (vbinary->filemode == FILE_MODE_APPEND) ? "a" : "w";
1435cc843e7aSLisandro Dalcin       vbinary->fdes_info = fopen(infoname, omode);
143628b400f6SJacob Faibussowitsch       PetscCheck(vbinary->fdes_info, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open .info file %s for writing", infoname);
14375c6c1daeSBarry Smith     }
14385c6c1daeSBarry Smith   }
14395c6c1daeSBarry Smith   PetscFunctionReturn(0);
14405c6c1daeSBarry Smith }
14415c6c1daeSBarry Smith 
1442d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerSetUp_Binary(PetscViewer viewer)
1443d71ae5a4SJacob Faibussowitsch {
14445c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1445cc843e7aSLisandro Dalcin   PetscBool           usempiio;
1446cc843e7aSLisandro Dalcin 
14475c6c1daeSBarry Smith   PetscFunctionBegin;
14489566063dSJacob Faibussowitsch   if (!vbinary->setfromoptionscalled) PetscCall(PetscViewerSetFromOptions(viewer));
144928b400f6SJacob Faibussowitsch   PetscCheck(vbinary->filename, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetName()");
145008401ef6SPierre Jolivet   PetscCheck(vbinary->filemode != (PetscFileMode)-1, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode()");
14519566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_Binary(viewer));
1452cc843e7aSLisandro Dalcin 
14539566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetUseMPIIO(viewer, &usempiio));
1454cc843e7aSLisandro Dalcin   if (usempiio) {
1455cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
14569566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetUp_BinaryMPIIO(viewer));
1457cc843e7aSLisandro Dalcin #endif
1458cc843e7aSLisandro Dalcin   } else {
14599566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetUp_BinarySTDIO(viewer));
1460cc843e7aSLisandro Dalcin   }
14619566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetUp_BinaryInfo(viewer));
1462cc843e7aSLisandro Dalcin 
14639566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectState((PetscObject)viewer, "File: %s", vbinary->filename));
14645c6c1daeSBarry Smith   PetscFunctionReturn(0);
14655c6c1daeSBarry Smith }
14665c6c1daeSBarry Smith 
1467d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerView_Binary(PetscViewer v, PetscViewer viewer)
1468d71ae5a4SJacob Faibussowitsch {
1469cb6ad94fSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
1470cb6ad94fSLisandro Dalcin   const char         *fname   = vbinary->filename ? vbinary->filename : "not yet set";
1471cc843e7aSLisandro Dalcin   const char         *fmode   = vbinary->filemode != (PetscFileMode)-1 ? PetscFileModes[vbinary->filemode] : "not yet set";
1472cb6ad94fSLisandro Dalcin   PetscBool           usempiio;
14732bf49c77SBarry Smith 
14742bf49c77SBarry Smith   PetscFunctionBegin;
14759566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetUseMPIIO(v, &usempiio));
14769566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Filename: %s\n", fname));
14779566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Mode: %s (%s)\n", fmode, usempiio ? "mpiio" : "stdio"));
14782bf49c77SBarry Smith   PetscFunctionReturn(0);
14792bf49c77SBarry Smith }
14802bf49c77SBarry Smith 
1481d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerSetFromOptions_Binary(PetscViewer viewer, PetscOptionItems *PetscOptionsObject)
1482d71ae5a4SJacob Faibussowitsch {
148322a8f86cSLisandro Dalcin   PetscViewer_Binary *binary = (PetscViewer_Binary *)viewer->data;
1484d22fd6bcSDave May   char                defaultname[PETSC_MAX_PATH_LEN];
148503a1d59fSDave May   PetscBool           flg;
1486e0385b85SDave May 
148703a1d59fSDave May   PetscFunctionBegin;
148822a8f86cSLisandro Dalcin   if (viewer->setupcalled) PetscFunctionReturn(0);
1489d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Binary PetscViewer Options");
14909566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(defaultname, PETSC_MAX_PATH_LEN - 1, "binaryoutput"));
14919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-viewer_binary_filename", "Specify filename", "PetscViewerFileSetName", defaultname, defaultname, sizeof(defaultname), &flg));
14929566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscViewerFileSetName_Binary(viewer, defaultname));
14939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_skip_info", "Skip writing/reading .info file", "PetscViewerBinarySetSkipInfo", binary->skipinfo, &binary->skipinfo, NULL));
14949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_skip_options", "Skip parsing Vec/Mat load options", "PetscViewerBinarySetSkipOptions", binary->skipoptions, &binary->skipoptions, NULL));
14959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_skip_header", "Skip writing/reading header information", "PetscViewerBinarySetSkipHeader", binary->skipheader, &binary->skipheader, NULL));
149603a1d59fSDave May #if defined(PETSC_HAVE_MPIIO)
14979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_mpiio", "Use MPI-IO functionality to write/read binary file", "PetscViewerBinarySetUseMPIIO", binary->usempiio, &binary->usempiio, NULL));
14985972f5f3SLisandro Dalcin #else
14999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_mpiio", "Use MPI-IO functionality to write/read binary file (NOT AVAILABLE)", "PetscViewerBinarySetUseMPIIO", PETSC_FALSE, NULL, NULL));
150003a1d59fSDave May #endif
1501d0609cedSBarry Smith   PetscOptionsHeadEnd();
1502bc196f7cSDave May   binary->setfromoptionscalled = PETSC_TRUE;
150303a1d59fSDave May   PetscFunctionReturn(0);
150403a1d59fSDave May }
150503a1d59fSDave May 
15068556b5ebSBarry Smith /*MC
15078556b5ebSBarry Smith    PETSCVIEWERBINARY - A viewer that saves to binary files
15088556b5ebSBarry Smith 
15091b266c99SBarry Smith   Level: beginner
15101b266c99SBarry Smith 
1511811af0c4SBarry Smith .seealso: `PetscViewerBinaryOpen()`, `PETSC_VIEWER_STDOUT_()`, `PETSC_VIEWER_STDOUT_SELF`, `PETSC_VIEWER_STDOUT_WORLD`, `PetscViewerCreate()`, `PetscViewerASCIIOpen()`,
1512811af0c4SBarry Smith           `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`, `PETSCVIEWERDRAW`, `PETSCVIEWERSOCKET`
1513811af0c4SBarry Smith           `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()`,
1514811af0c4SBarry Smith           `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`
15158556b5ebSBarry Smith M*/
15168556b5ebSBarry Smith 
1517d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscViewerCreate_Binary(PetscViewer v)
1518d71ae5a4SJacob Faibussowitsch {
15195c6c1daeSBarry Smith   PetscViewer_Binary *vbinary;
15205c6c1daeSBarry Smith 
15215c6c1daeSBarry Smith   PetscFunctionBegin;
15224dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&vbinary));
15235c6c1daeSBarry Smith   v->data = (void *)vbinary;
1524cc843e7aSLisandro Dalcin 
152503a1d59fSDave May   v->ops->setfromoptions   = PetscViewerSetFromOptions_Binary;
15265c6c1daeSBarry Smith   v->ops->destroy          = PetscViewerDestroy_Binary;
15272bf49c77SBarry Smith   v->ops->view             = PetscViewerView_Binary;
1528c98fd787SBarry Smith   v->ops->setup            = PetscViewerSetUp_Binary;
1529cc843e7aSLisandro Dalcin   v->ops->flush            = NULL; /* Should we support Flush() ? */
1530cc843e7aSLisandro Dalcin   v->ops->getsubviewer     = PetscViewerGetSubViewer_Binary;
1531cc843e7aSLisandro Dalcin   v->ops->restoresubviewer = PetscViewerRestoreSubViewer_Binary;
1532cc843e7aSLisandro Dalcin   v->ops->read             = PetscViewerBinaryRead;
1533cc843e7aSLisandro Dalcin 
1534cc843e7aSLisandro Dalcin   vbinary->fdes = -1;
1535e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
1536cc843e7aSLisandro Dalcin   vbinary->usempiio = PETSC_FALSE;
1537e8a65b7dSLisandro Dalcin   vbinary->mfdes    = MPI_FILE_NULL;
1538e8a65b7dSLisandro Dalcin   vbinary->mfsub    = MPI_FILE_NULL;
1539e8a65b7dSLisandro Dalcin #endif
1540cc843e7aSLisandro Dalcin   vbinary->filename        = NULL;
15417e4fd573SVaclav Hapla   vbinary->filemode        = FILE_MODE_UNDEFINED;
154202c9f0b5SLisandro Dalcin   vbinary->fdes_info       = NULL;
15435c6c1daeSBarry Smith   vbinary->skipinfo        = PETSC_FALSE;
15445c6c1daeSBarry Smith   vbinary->skipoptions     = PETSC_TRUE;
15455c6c1daeSBarry Smith   vbinary->skipheader      = PETSC_FALSE;
15465c6c1daeSBarry Smith   vbinary->storecompressed = PETSC_FALSE;
1547f90597f1SStefano Zampini   vbinary->ogzfilename     = NULL;
15485c6c1daeSBarry Smith   vbinary->flowcontrol     = 256; /* seems a good number for Cray XT-5 */
15495c6c1daeSBarry Smith 
1550cc843e7aSLisandro Dalcin   vbinary->setfromoptionscalled = PETSC_FALSE;
1551cc843e7aSLisandro Dalcin 
15529566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetFlowControl_C", PetscViewerBinaryGetFlowControl_Binary));
15539566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetFlowControl_C", PetscViewerBinarySetFlowControl_Binary));
15549566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipHeader_C", PetscViewerBinaryGetSkipHeader_Binary));
15559566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipHeader_C", PetscViewerBinarySetSkipHeader_Binary));
15569566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipOptions_C", PetscViewerBinaryGetSkipOptions_Binary));
15579566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipOptions_C", PetscViewerBinarySetSkipOptions_Binary));
15589566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipInfo_C", PetscViewerBinaryGetSkipInfo_Binary));
15599566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipInfo_C", PetscViewerBinarySetSkipInfo_Binary));
15609566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetInfoPointer_C", PetscViewerBinaryGetInfoPointer_Binary));
15619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_Binary));
15629566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_Binary));
15639566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_Binary));
15649566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_Binary));
15655c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
15669566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetUseMPIIO_C", PetscViewerBinaryGetUseMPIIO_Binary));
15679566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetUseMPIIO_C", PetscViewerBinarySetUseMPIIO_Binary));
15685c6c1daeSBarry Smith #endif
15695c6c1daeSBarry Smith   PetscFunctionReturn(0);
15705c6c1daeSBarry Smith }
15715c6c1daeSBarry Smith 
15725c6c1daeSBarry Smith /* ---------------------------------------------------------------------*/
15735c6c1daeSBarry Smith /*
15745c6c1daeSBarry Smith     The variable Petsc_Viewer_Binary_keyval is used to indicate an MPI attribute that
15755c6c1daeSBarry Smith   is attached to a communicator, in this case the attribute is a PetscViewer.
15765c6c1daeSBarry Smith */
1577d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_Binary_keyval = MPI_KEYVAL_INVALID;
15785c6c1daeSBarry Smith 
15795c6c1daeSBarry Smith /*@C
1580811af0c4SBarry Smith      PETSC_VIEWER_BINARY_ - Creates a `PETSCVIEWERBINARY` `PetscViewer` shared by all processors
15815c6c1daeSBarry Smith                      in a communicator.
15825c6c1daeSBarry Smith 
1583d083f849SBarry Smith      Collective
15845c6c1daeSBarry Smith 
15855c6c1daeSBarry Smith      Input Parameter:
1586811af0c4SBarry Smith .    comm - the MPI communicator to share the `PETSCVIEWERBINARY`
15875c6c1daeSBarry Smith 
15885c6c1daeSBarry Smith      Level: intermediate
15895c6c1daeSBarry Smith 
15905c6c1daeSBarry Smith    Options Database Keys:
159110699b91SBarry Smith +    -viewer_binary_filename <name> - filename in which to store the binary data, defaults to binaryoutput
159210699b91SBarry Smith .    -viewer_binary_skip_info - true means do not create .info file for this viewer
159310699b91SBarry Smith .    -viewer_binary_skip_options - true means do not use the options database for this viewer
159410699b91SBarry Smith .    -viewer_binary_skip_header - true means do not store the usual header information in the binary file
159510699b91SBarry Smith -    -viewer_binary_mpiio - true means use the file via MPI-IO, maybe faster for large files and many MPI ranks
15965c6c1daeSBarry Smith 
1597811af0c4SBarry Smith    Environmental variable:
159810699b91SBarry Smith -   PETSC_VIEWER_BINARY_FILENAME - filename in which to store the binary data, defaults to binaryoutput
15995c6c1daeSBarry Smith 
1600811af0c4SBarry Smith      Note:
1601811af0c4SBarry Smith      Unlike almost all other PETSc routines, `PETSC_VIEWER_BINARY_` does not return
16025c6c1daeSBarry Smith      an error code.  The binary PetscViewer is usually used in the form
16035c6c1daeSBarry Smith $       XXXView(XXX object,PETSC_VIEWER_BINARY_(comm));
16045c6c1daeSBarry Smith 
1605811af0c4SBarry Smith .seealso: `PETSCVIEWERBINARY`, `PETSC_VIEWER_BINARY_WORLD`, `PETSC_VIEWER_BINARY_SELF`, `PetscViewerBinaryOpen()`, `PetscViewerCreate()`,
1606db781477SPatrick Sanan           `PetscViewerDestroy()`
16075c6c1daeSBarry Smith @*/
1608d71ae5a4SJacob Faibussowitsch PetscViewer PETSC_VIEWER_BINARY_(MPI_Comm comm)
1609d71ae5a4SJacob Faibussowitsch {
16105c6c1daeSBarry Smith   PetscErrorCode ierr;
16115c6c1daeSBarry Smith   PetscBool      flg;
16125c6c1daeSBarry Smith   PetscViewer    viewer;
16135c6c1daeSBarry Smith   char           fname[PETSC_MAX_PATH_LEN];
16145c6c1daeSBarry Smith   MPI_Comm       ncomm;
16155c6c1daeSBarry Smith 
16165c6c1daeSBarry Smith   PetscFunctionBegin;
16179371c9d4SSatish Balay   ierr = PetscCommDuplicate(comm, &ncomm, NULL);
16189371c9d4SSatish Balay   if (ierr) {
16199371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16209371c9d4SSatish Balay     PetscFunctionReturn(NULL);
16219371c9d4SSatish Balay   }
16225c6c1daeSBarry Smith   if (Petsc_Viewer_Binary_keyval == MPI_KEYVAL_INVALID) {
162302c9f0b5SLisandro Dalcin     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Viewer_Binary_keyval, NULL);
16249371c9d4SSatish Balay     if (ierr) {
16259371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16269371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16279371c9d4SSatish Balay     }
16285c6c1daeSBarry Smith   }
162947435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm, Petsc_Viewer_Binary_keyval, (void **)&viewer, (int *)&flg);
16309371c9d4SSatish Balay   if (ierr) {
16319371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16329371c9d4SSatish Balay     PetscFunctionReturn(NULL);
16339371c9d4SSatish Balay   }
16345c6c1daeSBarry Smith   if (!flg) { /* PetscViewer not yet created */
16355c6c1daeSBarry Smith     ierr = PetscOptionsGetenv(ncomm, "PETSC_VIEWER_BINARY_FILENAME", fname, PETSC_MAX_PATH_LEN, &flg);
16369371c9d4SSatish Balay     if (ierr) {
16379371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16389371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16399371c9d4SSatish Balay     }
16405c6c1daeSBarry Smith     if (!flg) {
16415c6c1daeSBarry Smith       ierr = PetscStrcpy(fname, "binaryoutput");
16429371c9d4SSatish Balay       if (ierr) {
16439371c9d4SSatish Balay         PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16449371c9d4SSatish Balay         PetscFunctionReturn(NULL);
16459371c9d4SSatish Balay       }
16465c6c1daeSBarry Smith     }
16475c6c1daeSBarry Smith     ierr = PetscViewerBinaryOpen(ncomm, fname, FILE_MODE_WRITE, &viewer);
16489371c9d4SSatish Balay     if (ierr) {
16499371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16509371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16519371c9d4SSatish Balay     }
16525c6c1daeSBarry Smith     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
16539371c9d4SSatish Balay     if (ierr) {
16549371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16559371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16569371c9d4SSatish Balay     }
165747435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm, Petsc_Viewer_Binary_keyval, (void *)viewer);
16589371c9d4SSatish Balay     if (ierr) {
16599371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16609371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16619371c9d4SSatish Balay     }
16625c6c1daeSBarry Smith   }
16635c6c1daeSBarry Smith   ierr = PetscCommDestroy(&ncomm);
16649371c9d4SSatish Balay   if (ierr) {
16659371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16669371c9d4SSatish Balay     PetscFunctionReturn(NULL);
16679371c9d4SSatish Balay   }
16685c6c1daeSBarry Smith   PetscFunctionReturn(viewer);
16695c6c1daeSBarry Smith }
1670