xref: /petsc/src/sys/classes/viewer/impls/binary/binv.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1af0996ceSBarry Smith #include <petsc/private/viewerimpl.h> /*I   "petscviewer.h"   I*/
25c6c1daeSBarry Smith 
35c6c1daeSBarry Smith typedef struct {
45c6c1daeSBarry Smith   int fdes; /* file descriptor, ignored if using MPI IO */
55c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
6bc196f7cSDave May   PetscBool  usempiio;
75c6c1daeSBarry Smith   MPI_File   mfdes; /* ignored unless using MPI IO */
8e8a65b7dSLisandro Dalcin   MPI_File   mfsub; /* subviewer support */
95c6c1daeSBarry Smith   MPI_Offset moff;
105c6c1daeSBarry Smith #endif
11cc843e7aSLisandro Dalcin   char         *filename;            /* file name */
12cc843e7aSLisandro Dalcin   PetscFileMode filemode;            /* read/write/append mode */
135c6c1daeSBarry Smith   FILE         *fdes_info;           /* optional file containing info on binary file*/
145c6c1daeSBarry Smith   PetscBool     storecompressed;     /* gzip the write binary file when closing it*/
15f90597f1SStefano Zampini   char         *ogzfilename;         /* gzip can be run after the filename has been updated */
165c6c1daeSBarry Smith   PetscBool     skipinfo;            /* Don't create info file for writing; don't use for reading */
175c6c1daeSBarry Smith   PetscBool     skipoptions;         /* don't use PETSc options database when loading */
185c6c1daeSBarry Smith   PetscInt      flowcontrol;         /* allow only <flowcontrol> messages outstanding at a time while doing IO */
195c6c1daeSBarry Smith   PetscBool     skipheader;          /* don't write header, only raw data */
20a261c58fSBarry Smith   PetscBool     matlabheaderwritten; /* if format is PETSC_VIEWER_BINARY_MATLAB has the MATLAB .info header been written yet */
21c98fd787SBarry Smith   PetscBool     setfromoptionscalled;
225c6c1daeSBarry Smith } PetscViewer_Binary;
235c6c1daeSBarry Smith 
249371c9d4SSatish Balay static PetscErrorCode PetscViewerBinaryClearFunctionList(PetscViewer v) {
252e956fe4SStefano Zampini   PetscFunctionBegin;
262e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetFlowControl_C", NULL));
272e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetFlowControl_C", NULL));
282e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipHeader_C", NULL));
292e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipHeader_C", NULL));
302e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipOptions_C", NULL));
312e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipOptions_C", NULL));
322e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipInfo_C", NULL));
332e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipInfo_C", NULL));
342e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetInfoPointer_C", NULL));
352e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", NULL));
362e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", NULL));
372e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", NULL));
382e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", NULL));
392e956fe4SStefano Zampini #if defined(PETSC_HAVE_MPIIO)
402e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetUseMPIIO_C", NULL));
412e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetUseMPIIO_C", NULL));
422e956fe4SStefano Zampini #endif
432e956fe4SStefano Zampini   PetscFunctionReturn(0);
442e956fe4SStefano Zampini }
452e956fe4SStefano Zampini 
46cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
479371c9d4SSatish Balay static PetscErrorCode PetscViewerBinarySyncMPIIO(PetscViewer viewer) {
48cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
49cc843e7aSLisandro Dalcin 
50cc843e7aSLisandro Dalcin   PetscFunctionBegin;
51cc843e7aSLisandro Dalcin   if (vbinary->filemode == FILE_MODE_READ) PetscFunctionReturn(0);
52*48a46eb9SPierre Jolivet   if (vbinary->mfsub != MPI_FILE_NULL) PetscCallMPI(MPI_File_sync(vbinary->mfsub));
53cc843e7aSLisandro Dalcin   if (vbinary->mfdes != MPI_FILE_NULL) {
549566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer)));
559566063dSJacob Faibussowitsch     PetscCallMPI(MPI_File_sync(vbinary->mfdes));
56cc843e7aSLisandro Dalcin   }
57cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
58cc843e7aSLisandro Dalcin }
59cc843e7aSLisandro Dalcin #endif
60cc843e7aSLisandro Dalcin 
619371c9d4SSatish Balay static PetscErrorCode PetscViewerGetSubViewer_Binary(PetscViewer viewer, MPI_Comm comm, PetscViewer *outviewer) {
62e8a65b7dSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
63cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
645c6c1daeSBarry Smith 
655c6c1daeSBarry Smith   PetscFunctionBegin;
669566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
67e8a65b7dSLisandro Dalcin 
68e8a65b7dSLisandro Dalcin   /* Return subviewer in process zero */
699566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
70dd400576SPatrick Sanan   if (rank == 0) {
71e5afcf28SBarry Smith     PetscMPIInt flg;
72e5afcf28SBarry Smith 
739566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_compare(PETSC_COMM_SELF, comm, &flg));
74cc73adaaSBarry Smith     PetscCheck(flg == MPI_IDENT || flg == MPI_CONGRUENT, PETSC_COMM_SELF, PETSC_ERR_SUP, "PetscViewerGetSubViewer() for PETSCVIEWERBINARY requires a singleton MPI_Comm");
759566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(comm, outviewer));
769566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(*outviewer, PETSCVIEWERBINARY));
779566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy((*outviewer)->data, vbinary, sizeof(PetscViewer_Binary)));
7812f4c3a9SLisandro Dalcin     (*outviewer)->setupcalled = PETSC_TRUE;
7996509d17SLisandro Dalcin   } else {
8096509d17SLisandro Dalcin     *outviewer = NULL;
8196509d17SLisandro Dalcin   }
82e8a65b7dSLisandro Dalcin 
83e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
84e8a65b7dSLisandro Dalcin   if (vbinary->usempiio && *outviewer) {
85e8a65b7dSLisandro Dalcin     PetscViewer_Binary *obinary = (PetscViewer_Binary *)(*outviewer)->data;
86e8a65b7dSLisandro Dalcin     /* Parent viewer opens a new MPI file handle on PETSC_COMM_SELF and keeps track of it for future reuse */
87e8a65b7dSLisandro Dalcin     if (vbinary->mfsub == MPI_FILE_NULL) {
88e8a65b7dSLisandro Dalcin       int amode;
89cc843e7aSLisandro Dalcin       switch (vbinary->filemode) {
90e8a65b7dSLisandro Dalcin       case FILE_MODE_READ: amode = MPI_MODE_RDONLY; break;
91e8a65b7dSLisandro Dalcin       case FILE_MODE_WRITE: amode = MPI_MODE_WRONLY; break;
9222a8f86cSLisandro Dalcin       case FILE_MODE_APPEND: amode = MPI_MODE_WRONLY; break;
9398921bdaSJacob Faibussowitsch       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[vbinary->filemode]);
94e8a65b7dSLisandro Dalcin       }
959566063dSJacob Faibussowitsch       PetscCallMPI(MPI_File_open(PETSC_COMM_SELF, vbinary->filename, amode, MPI_INFO_NULL, &vbinary->mfsub));
96e8a65b7dSLisandro Dalcin     }
97e8a65b7dSLisandro Dalcin     /* Subviewer gets the MPI file handle on PETSC_COMM_SELF */
98e8a65b7dSLisandro Dalcin     obinary->mfdes = vbinary->mfsub;
99e8a65b7dSLisandro Dalcin     obinary->mfsub = MPI_FILE_NULL;
100e8a65b7dSLisandro Dalcin     obinary->moff  = vbinary->moff;
101e8a65b7dSLisandro Dalcin   }
102e8a65b7dSLisandro Dalcin #endif
103cc843e7aSLisandro Dalcin 
104cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
1059566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinarySyncMPIIO(viewer));
106cc843e7aSLisandro Dalcin #endif
1075c6c1daeSBarry Smith   PetscFunctionReturn(0);
1085c6c1daeSBarry Smith }
1095c6c1daeSBarry Smith 
1109371c9d4SSatish Balay static PetscErrorCode PetscViewerRestoreSubViewer_Binary(PetscViewer viewer, MPI_Comm comm, PetscViewer *outviewer) {
111e8a65b7dSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
112cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
113e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
114e8a65b7dSLisandro Dalcin   MPI_Offset moff = 0;
115e8a65b7dSLisandro Dalcin #endif
1165c6c1daeSBarry Smith 
1175c6c1daeSBarry Smith   PetscFunctionBegin;
1189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
119c5853193SPierre Jolivet   PetscCheck(rank == 0 || !*outviewer, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Subviewer not obtained from viewer");
120e8a65b7dSLisandro Dalcin 
121e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
122e8a65b7dSLisandro Dalcin   if (vbinary->usempiio && *outviewer) {
123e8a65b7dSLisandro Dalcin     PetscViewer_Binary *obinary = (PetscViewer_Binary *)(*outviewer)->data;
12408401ef6SPierre Jolivet     PetscCheck(obinary->mfdes == vbinary->mfsub, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Subviewer not obtained from viewer");
1259566063dSJacob Faibussowitsch     if (obinary->mfsub != MPI_FILE_NULL) PetscCallMPI(MPI_File_close(&obinary->mfsub));
126e8a65b7dSLisandro Dalcin     moff = obinary->moff;
127e8a65b7dSLisandro Dalcin   }
128e8a65b7dSLisandro Dalcin #endif
129e8a65b7dSLisandro Dalcin 
130e8a65b7dSLisandro Dalcin   if (*outviewer) {
131e8a65b7dSLisandro Dalcin     PetscViewer_Binary *obinary = (PetscViewer_Binary *)(*outviewer)->data;
13208401ef6SPierre Jolivet     PetscCheck(obinary->fdes == vbinary->fdes, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Subviewer not obtained from viewer");
1339566063dSJacob Faibussowitsch     PetscCall(PetscFree((*outviewer)->data));
1342e956fe4SStefano Zampini     PetscCall(PetscViewerBinaryClearFunctionList(*outviewer));
1359566063dSJacob Faibussowitsch     PetscCall(PetscHeaderDestroy(outviewer));
1365c6c1daeSBarry Smith   }
137e8a65b7dSLisandro Dalcin 
138e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
139e8a65b7dSLisandro Dalcin   if (vbinary->usempiio) {
140e8a65b7dSLisandro Dalcin     PetscInt64 ioff = (PetscInt64)moff; /* We could use MPI_OFFSET datatype (requires MPI 2.2) */
1419566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&ioff, 1, MPIU_INT64, 0, PetscObjectComm((PetscObject)viewer)));
142e8a65b7dSLisandro Dalcin     vbinary->moff = (MPI_Offset)ioff;
143e8a65b7dSLisandro Dalcin   }
144e8a65b7dSLisandro Dalcin #endif
145cc843e7aSLisandro Dalcin 
146cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
1479566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinarySyncMPIIO(viewer));
148cc843e7aSLisandro Dalcin #endif
1495c6c1daeSBarry Smith   PetscFunctionReturn(0);
1505c6c1daeSBarry Smith }
1515c6c1daeSBarry Smith 
1525c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
1535c6c1daeSBarry Smith /*@C
15485b8072bSPatrick Sanan     PetscViewerBinaryGetMPIIOOffset - Gets the current global offset that should be passed to `MPI_File_set_view()` or `MPI_File_{write|read}_at[_all]()`
1555c6c1daeSBarry Smith 
1565c6c1daeSBarry Smith     Not Collective
1575c6c1daeSBarry Smith 
1585c6c1daeSBarry Smith     Input Parameter:
15985b8072bSPatrick Sanan .   viewer - PetscViewer context, obtained from `PetscViewerBinaryOpen()`
1605c6c1daeSBarry Smith 
1615c6c1daeSBarry Smith     Output Parameter:
16222a8f86cSLisandro Dalcin .   off - the current global offset
1635c6c1daeSBarry Smith 
1645c6c1daeSBarry Smith     Level: advanced
1655c6c1daeSBarry Smith 
1665c6c1daeSBarry Smith     Fortran Note:
1675c6c1daeSBarry Smith     This routine is not supported in Fortran.
1685c6c1daeSBarry Smith 
16985b8072bSPatrick Sanan     Use `PetscViewerBinaryAddMPIIOOffset()` to increase this value after you have written a view.
1705c6c1daeSBarry Smith 
171db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryAddMPIIOOffset()`
1725c6c1daeSBarry Smith @*/
1739371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryGetMPIIOOffset(PetscViewer viewer, MPI_Offset *off) {
17422a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
1755c6c1daeSBarry Smith 
1765c6c1daeSBarry Smith   PetscFunctionBegin;
17722a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
17822a8f86cSLisandro Dalcin   PetscValidPointer(off, 2);
17922a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
1805c6c1daeSBarry Smith   *off    = vbinary->moff;
1815c6c1daeSBarry Smith   PetscFunctionReturn(0);
1825c6c1daeSBarry Smith }
1835c6c1daeSBarry Smith 
1845c6c1daeSBarry Smith /*@C
18522a8f86cSLisandro Dalcin     PetscViewerBinaryAddMPIIOOffset - Adds to the current global offset
1865c6c1daeSBarry Smith 
18722a8f86cSLisandro Dalcin     Logically Collective
1885c6c1daeSBarry Smith 
1895c6c1daeSBarry Smith     Input Parameters:
1905c6c1daeSBarry Smith +   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
19122a8f86cSLisandro Dalcin -   off - the addition to the global offset
1925c6c1daeSBarry Smith 
1935c6c1daeSBarry Smith     Level: advanced
1945c6c1daeSBarry Smith 
1955c6c1daeSBarry Smith     Fortran Note:
1965c6c1daeSBarry Smith     This routine is not supported in Fortran.
1975c6c1daeSBarry Smith 
19885b8072bSPatrick Sanan     Use `PetscViewerBinaryGetMPIIOOffset()` to get the value that you should pass to `MPI_File_set_view()` or `MPI_File_{write|read}_at[_all]()`
1995c6c1daeSBarry Smith 
200c2e3fba1SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
2015c6c1daeSBarry Smith @*/
2029371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryAddMPIIOOffset(PetscViewer viewer, MPI_Offset off) {
20322a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
2045c6c1daeSBarry Smith 
2055c6c1daeSBarry Smith   PetscFunctionBegin;
20622a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
20722a8f86cSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, (PetscInt)off, 2);
20822a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
2095c6c1daeSBarry Smith   vbinary->moff += off;
2105c6c1daeSBarry Smith   PetscFunctionReturn(0);
2115c6c1daeSBarry Smith }
2125c6c1daeSBarry Smith 
2135c6c1daeSBarry Smith /*@C
2145c6c1daeSBarry Smith     PetscViewerBinaryGetMPIIODescriptor - Extracts the MPI IO file descriptor from a PetscViewer.
2155c6c1daeSBarry Smith 
2165c6c1daeSBarry Smith     Not Collective
2175c6c1daeSBarry Smith 
2185c6c1daeSBarry Smith     Input Parameter:
2195c6c1daeSBarry Smith .   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
2205c6c1daeSBarry Smith 
2215c6c1daeSBarry Smith     Output Parameter:
2225c6c1daeSBarry Smith .   fdes - file descriptor
2235c6c1daeSBarry Smith 
2245c6c1daeSBarry Smith     Level: advanced
2255c6c1daeSBarry Smith 
2265c6c1daeSBarry Smith     Fortran Note:
2275c6c1daeSBarry Smith     This routine is not supported in Fortran.
2285c6c1daeSBarry Smith 
229c2e3fba1SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
2305c6c1daeSBarry Smith @*/
2319371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryGetMPIIODescriptor(PetscViewer viewer, MPI_File *fdes) {
23222a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
2335c6c1daeSBarry Smith 
2345c6c1daeSBarry Smith   PetscFunctionBegin;
23522a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
23622a8f86cSLisandro Dalcin   PetscValidPointer(fdes, 2);
2379566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
23822a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
2395c6c1daeSBarry Smith   *fdes   = vbinary->mfdes;
2405c6c1daeSBarry Smith   PetscFunctionReturn(0);
2415c6c1daeSBarry Smith }
242cc843e7aSLisandro Dalcin #endif
2435c6c1daeSBarry Smith 
244cc843e7aSLisandro Dalcin /*@
245cc843e7aSLisandro Dalcin     PetscViewerBinarySetUseMPIIO - Sets a binary viewer to use MPI-IO for reading/writing. Must be called
246cc843e7aSLisandro Dalcin         before PetscViewerFileSetName()
247cc843e7aSLisandro Dalcin 
248cc843e7aSLisandro Dalcin     Logically Collective on PetscViewer
249cc843e7aSLisandro Dalcin 
250cc843e7aSLisandro Dalcin     Input Parameters:
251cc843e7aSLisandro Dalcin +   viewer - the PetscViewer; must be a binary
252cc843e7aSLisandro Dalcin -   use - PETSC_TRUE means MPI-IO will be used
253cc843e7aSLisandro Dalcin 
254cc843e7aSLisandro Dalcin     Options Database:
255cc843e7aSLisandro Dalcin     -viewer_binary_mpiio : Flag for using MPI-IO
256cc843e7aSLisandro Dalcin 
257cc843e7aSLisandro Dalcin     Level: advanced
258cc843e7aSLisandro Dalcin 
259db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`,
260db781477SPatrick Sanan           `PetscViewerBinaryGetUseMPIIO()`
261cc843e7aSLisandro Dalcin 
262cc843e7aSLisandro Dalcin @*/
2639371c9d4SSatish Balay PetscErrorCode PetscViewerBinarySetUseMPIIO(PetscViewer viewer, PetscBool use) {
264a8aae444SDave May   PetscFunctionBegin;
265cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
266cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, use, 2);
267cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetUseMPIIO_C", (PetscViewer, PetscBool), (viewer, use));
268cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
269cc843e7aSLisandro Dalcin }
270cc843e7aSLisandro Dalcin 
271cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
2729371c9d4SSatish Balay static PetscErrorCode PetscViewerBinarySetUseMPIIO_Binary(PetscViewer viewer, PetscBool use) {
273cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
274cc843e7aSLisandro Dalcin   PetscFunctionBegin;
275cc73adaaSBarry Smith   PetscCheck(!viewer->setupcalled || vbinary->usempiio == use, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Cannot change MPIIO to %s after setup", PetscBools[use]);
276cc843e7aSLisandro Dalcin   vbinary->usempiio = use;
277a8aae444SDave May   PetscFunctionReturn(0);
278a8aae444SDave May }
279a8aae444SDave May #endif
280a8aae444SDave May 
281cc843e7aSLisandro Dalcin /*@
282bc196f7cSDave May     PetscViewerBinaryGetUseMPIIO - Returns PETSC_TRUE if the binary viewer uses MPI-IO.
2835c6c1daeSBarry Smith 
2845c6c1daeSBarry Smith     Not Collective
2855c6c1daeSBarry Smith 
2865c6c1daeSBarry Smith     Input Parameter:
2875c6c1daeSBarry Smith .   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
2885c6c1daeSBarry Smith 
2895c6c1daeSBarry Smith     Output Parameter:
290d8d19677SJose E. Roman .   use - PETSC_TRUE if MPI-IO is being used
2915c6c1daeSBarry Smith 
2925c6c1daeSBarry Smith     Options Database:
2935c6c1daeSBarry Smith     -viewer_binary_mpiio : Flag for using MPI-IO
2945c6c1daeSBarry Smith 
2955c6c1daeSBarry Smith     Level: advanced
2965c6c1daeSBarry Smith 
297bc196f7cSDave May     Note:
298bc196f7cSDave May     If MPI-IO is not available, this function will always return PETSC_FALSE
299bc196f7cSDave May 
3005c6c1daeSBarry Smith     Fortran Note:
3015c6c1daeSBarry Smith     This routine is not supported in Fortran.
3025c6c1daeSBarry Smith 
303db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
3045c6c1daeSBarry Smith @*/
3059371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryGetUseMPIIO(PetscViewer viewer, PetscBool *use) {
3065c6c1daeSBarry Smith   PetscFunctionBegin;
307cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
308cc843e7aSLisandro Dalcin   PetscValidBoolPointer(use, 2);
309cc843e7aSLisandro Dalcin   *use = PETSC_FALSE;
310cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinaryGetUseMPIIO_C", (PetscViewer, PetscBool *), (viewer, use));
3115c6c1daeSBarry Smith   PetscFunctionReturn(0);
3125c6c1daeSBarry Smith }
3135c6c1daeSBarry Smith 
314cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
3159371c9d4SSatish Balay static PetscErrorCode PetscViewerBinaryGetUseMPIIO_Binary(PetscViewer viewer, PetscBool *use) {
3165c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
3175c6c1daeSBarry Smith 
3185c6c1daeSBarry Smith   PetscFunctionBegin;
319cc843e7aSLisandro Dalcin   *use = vbinary->usempiio;
320cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
321cc843e7aSLisandro Dalcin }
322cc843e7aSLisandro Dalcin #endif
323cc843e7aSLisandro Dalcin 
324cc843e7aSLisandro Dalcin /*@
325cc843e7aSLisandro Dalcin     PetscViewerBinarySetFlowControl - Sets how many messages are allowed to outstanding at the same time during parallel IO reads/writes
326cc843e7aSLisandro Dalcin 
327cc843e7aSLisandro Dalcin     Not Collective
328cc843e7aSLisandro Dalcin 
329d8d19677SJose E. Roman     Input Parameters:
330cc843e7aSLisandro Dalcin +   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
331cc843e7aSLisandro Dalcin -   fc - the number of messages, defaults to 256 if this function was not called
332cc843e7aSLisandro Dalcin 
333cc843e7aSLisandro Dalcin     Level: advanced
334cc843e7aSLisandro Dalcin 
335c2e3fba1SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetFlowControl()`
336cc843e7aSLisandro Dalcin 
337cc843e7aSLisandro Dalcin @*/
3389371c9d4SSatish Balay PetscErrorCode PetscViewerBinarySetFlowControl(PetscViewer viewer, PetscInt fc) {
339cc843e7aSLisandro Dalcin   PetscFunctionBegin;
340cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
341cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, fc, 2);
342cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetFlowControl_C", (PetscViewer, PetscInt), (viewer, fc));
3435c6c1daeSBarry Smith   PetscFunctionReturn(0);
3445c6c1daeSBarry Smith }
3455c6c1daeSBarry Smith 
3469371c9d4SSatish Balay static PetscErrorCode PetscViewerBinarySetFlowControl_Binary(PetscViewer viewer, PetscInt fc) {
347cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
348cc843e7aSLisandro Dalcin 
349cc843e7aSLisandro Dalcin   PetscFunctionBegin;
35008401ef6SPierre Jolivet   PetscCheck(fc > 1, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_OUTOFRANGE, "Flow control count must be greater than 1, %" PetscInt_FMT " was set", fc);
351cc843e7aSLisandro Dalcin   vbinary->flowcontrol = fc;
352cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
353cc843e7aSLisandro Dalcin }
354cc843e7aSLisandro Dalcin 
355cc843e7aSLisandro Dalcin /*@
3565c6c1daeSBarry Smith     PetscViewerBinaryGetFlowControl - Returns how many messages are allowed to outstanding at the same time during parallel IO reads/writes
3575c6c1daeSBarry Smith 
3585c6c1daeSBarry Smith     Not Collective
3595c6c1daeSBarry Smith 
3605c6c1daeSBarry Smith     Input Parameter:
3615c6c1daeSBarry Smith .   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
3625c6c1daeSBarry Smith 
3635c6c1daeSBarry Smith     Output Parameter:
3645c6c1daeSBarry Smith .   fc - the number of messages
3655c6c1daeSBarry Smith 
3665c6c1daeSBarry Smith     Level: advanced
3675c6c1daeSBarry Smith 
368c2e3fba1SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinarySetFlowControl()`
3695c6c1daeSBarry Smith 
3705c6c1daeSBarry Smith @*/
3719371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryGetFlowControl(PetscViewer viewer, PetscInt *fc) {
3725c6c1daeSBarry Smith   PetscFunctionBegin;
373cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
374cc843e7aSLisandro Dalcin   PetscValidIntPointer(fc, 2);
375cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetFlowControl_C", (PetscViewer, PetscInt *), (viewer, fc));
3765c6c1daeSBarry Smith   PetscFunctionReturn(0);
3775c6c1daeSBarry Smith }
3785c6c1daeSBarry Smith 
3799371c9d4SSatish Balay static PetscErrorCode PetscViewerBinaryGetFlowControl_Binary(PetscViewer viewer, PetscInt *fc) {
3805c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
3815c6c1daeSBarry Smith 
3825c6c1daeSBarry Smith   PetscFunctionBegin;
383cc843e7aSLisandro Dalcin   *fc = vbinary->flowcontrol;
3845c6c1daeSBarry Smith   PetscFunctionReturn(0);
3855c6c1daeSBarry Smith }
3865c6c1daeSBarry Smith 
3875c6c1daeSBarry Smith /*@C
3885c6c1daeSBarry Smith     PetscViewerBinaryGetDescriptor - Extracts the file descriptor from a PetscViewer.
3895c6c1daeSBarry Smith 
3905872f025SBarry Smith     Collective On PetscViewer
3915c6c1daeSBarry Smith 
3925c6c1daeSBarry Smith     Input Parameter:
3935c6c1daeSBarry Smith .   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
3945c6c1daeSBarry Smith 
3955c6c1daeSBarry Smith     Output Parameter:
3965c6c1daeSBarry Smith .   fdes - file descriptor
3975c6c1daeSBarry Smith 
3985c6c1daeSBarry Smith     Level: advanced
3995c6c1daeSBarry Smith 
4005c6c1daeSBarry Smith     Notes:
4015c6c1daeSBarry Smith       For writable binary PetscViewers, the descriptor will only be valid for the
4025c6c1daeSBarry Smith     first processor in the communicator that shares the PetscViewer. For readable
4035c6c1daeSBarry Smith     files it will only be valid on nodes that have the file. If node 0 does not
4045c6c1daeSBarry Smith     have the file it generates an error even if another node does have the file.
4055c6c1daeSBarry Smith 
4065c6c1daeSBarry Smith     Fortran Note:
4075c6c1daeSBarry Smith     This routine is not supported in Fortran.
4085c6c1daeSBarry Smith 
40995452b02SPatrick Sanan     Developer Notes:
41095452b02SPatrick Sanan     This must be called on all processes because Dave May changed
4115872f025SBarry Smith     the source code that this may be trigger a PetscViewerSetUp() call if it was not previously triggered.
4125872f025SBarry Smith 
413c2e3fba1SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`
4145c6c1daeSBarry Smith @*/
4159371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryGetDescriptor(PetscViewer viewer, int *fdes) {
41622a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
4175c6c1daeSBarry Smith 
4185c6c1daeSBarry Smith   PetscFunctionBegin;
41922a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
42022a8f86cSLisandro Dalcin   PetscValidPointer(fdes, 2);
4219566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
42222a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
4235c6c1daeSBarry Smith   *fdes   = vbinary->fdes;
4245c6c1daeSBarry Smith   PetscFunctionReturn(0);
4255c6c1daeSBarry Smith }
4265c6c1daeSBarry Smith 
4275c6c1daeSBarry Smith /*@
4285c6c1daeSBarry Smith     PetscViewerBinarySkipInfo - Binary file will not have .info file created with it
4295c6c1daeSBarry Smith 
4305c6c1daeSBarry Smith     Not Collective
4315c6c1daeSBarry Smith 
432fd292e60Sprj-     Input Parameter:
4335c6c1daeSBarry Smith .   viewer - PetscViewer context, obtained from PetscViewerCreate()
4345c6c1daeSBarry Smith 
4355c6c1daeSBarry Smith     Options Database Key:
43610699b91SBarry Smith .   -viewer_binary_skip_info - true indicates do not generate .info file
4375c6c1daeSBarry Smith 
4385c6c1daeSBarry Smith     Level: advanced
4395c6c1daeSBarry Smith 
44095452b02SPatrick Sanan     Notes:
44195452b02SPatrick Sanan     This must be called after PetscViewerSetType(). If you use PetscViewerBinaryOpen() then
4425c6c1daeSBarry Smith     you can only skip the info file with the -viewer_binary_skip_info flag. To use the function you must open the
443a2d7db39SDave May     viewer with PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinarySkipInfo().
4445c6c1daeSBarry Smith 
4455c6c1daeSBarry Smith     The .info contains meta information about the data in the binary file, for example the block size if it was
4465c6c1daeSBarry Smith     set for a vector or matrix.
4475c6c1daeSBarry Smith 
448db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySetSkipOptions()`,
449db781477SPatrick Sanan           `PetscViewerBinaryGetSkipOptions()`, `PetscViewerBinaryGetSkipInfo()`
4505c6c1daeSBarry Smith @*/
4519371c9d4SSatish Balay PetscErrorCode PetscViewerBinarySkipInfo(PetscViewer viewer) {
4525c6c1daeSBarry Smith   PetscFunctionBegin;
4539566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinarySetSkipInfo(viewer, PETSC_TRUE));
454807ea322SDave May   PetscFunctionReturn(0);
455807ea322SDave May }
456807ea322SDave May 
457807ea322SDave May /*@
458807ea322SDave May     PetscViewerBinarySetSkipInfo - Binary file will not have .info file created with it
459807ea322SDave May 
460807ea322SDave May     Not Collective
461807ea322SDave May 
462d8d19677SJose E. Roman     Input Parameters:
463cc843e7aSLisandro Dalcin +   viewer - PetscViewer context, obtained from PetscViewerCreate()
464cc843e7aSLisandro Dalcin -   skip - PETSC_TRUE implies the .info file will not be generated
465807ea322SDave May 
466807ea322SDave May     Options Database Key:
46710699b91SBarry Smith .   -viewer_binary_skip_info - true indicates do not generate .info file
468807ea322SDave May 
469807ea322SDave May     Level: advanced
470807ea322SDave May 
471db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySetSkipOptions()`,
472db781477SPatrick Sanan           `PetscViewerBinaryGetSkipOptions()`, `PetscViewerBinaryGetSkipInfo()`
473807ea322SDave May @*/
4749371c9d4SSatish Balay PetscErrorCode PetscViewerBinarySetSkipInfo(PetscViewer viewer, PetscBool skip) {
475807ea322SDave May   PetscFunctionBegin;
476cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
477cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, skip, 2);
478cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetSkipInfo_C", (PetscViewer, PetscBool), (viewer, skip));
479807ea322SDave May   PetscFunctionReturn(0);
480807ea322SDave May }
481807ea322SDave May 
4829371c9d4SSatish Balay static PetscErrorCode PetscViewerBinarySetSkipInfo_Binary(PetscViewer viewer, PetscBool skip) {
483807ea322SDave May   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
484807ea322SDave May 
485807ea322SDave May   PetscFunctionBegin;
486cc843e7aSLisandro Dalcin   vbinary->skipinfo = skip;
487807ea322SDave May   PetscFunctionReturn(0);
488807ea322SDave May }
489807ea322SDave May 
490807ea322SDave May /*@
491807ea322SDave May     PetscViewerBinaryGetSkipInfo - check if viewer wrote a .info file
492807ea322SDave May 
493807ea322SDave May     Not Collective
494807ea322SDave May 
495807ea322SDave May     Input Parameter:
496807ea322SDave May .   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
497807ea322SDave May 
498807ea322SDave May     Output Parameter:
499807ea322SDave May .   skip - PETSC_TRUE implies the .info file was not generated
500807ea322SDave May 
501807ea322SDave May     Level: advanced
502807ea322SDave May 
50395452b02SPatrick Sanan     Notes:
50495452b02SPatrick Sanan     This must be called after PetscViewerSetType()
505807ea322SDave May 
506db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
507db781477SPatrick Sanan           `PetscViewerBinarySetSkipOptions()`, `PetscViewerBinarySetSkipInfo()`
508807ea322SDave May @*/
5099371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryGetSkipInfo(PetscViewer viewer, PetscBool *skip) {
510807ea322SDave May   PetscFunctionBegin;
511cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
512cc843e7aSLisandro Dalcin   PetscValidBoolPointer(skip, 2);
513cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetSkipInfo_C", (PetscViewer, PetscBool *), (viewer, skip));
514807ea322SDave May   PetscFunctionReturn(0);
515807ea322SDave May }
516807ea322SDave May 
5179371c9d4SSatish Balay static PetscErrorCode PetscViewerBinaryGetSkipInfo_Binary(PetscViewer viewer, PetscBool *skip) {
518807ea322SDave May   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
519807ea322SDave May 
520807ea322SDave May   PetscFunctionBegin;
521cc843e7aSLisandro Dalcin   *skip = vbinary->skipinfo;
522807ea322SDave May   PetscFunctionReturn(0);
523807ea322SDave May }
524807ea322SDave May 
5255c6c1daeSBarry Smith /*@
5265c6c1daeSBarry Smith     PetscViewerBinarySetSkipOptions - do not use the PETSc options database when loading objects
5275c6c1daeSBarry Smith 
5285c6c1daeSBarry Smith     Not Collective
5295c6c1daeSBarry Smith 
5305c6c1daeSBarry Smith     Input Parameters:
5315c6c1daeSBarry Smith +   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
53210699b91SBarry Smith -   skip - PETSC_TRUE means do not use the options from the options database
5335c6c1daeSBarry Smith 
5345c6c1daeSBarry Smith     Options Database Key:
53510699b91SBarry Smith .   -viewer_binary_skip_options - true means do not use the options from the options database
5365c6c1daeSBarry Smith 
5375c6c1daeSBarry Smith     Level: advanced
5385c6c1daeSBarry Smith 
53995452b02SPatrick Sanan     Notes:
54095452b02SPatrick Sanan     This must be called after PetscViewerSetType()
5415c6c1daeSBarry Smith 
542db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
543db781477SPatrick Sanan           `PetscViewerBinaryGetSkipOptions()`
5445c6c1daeSBarry Smith @*/
5459371c9d4SSatish Balay PetscErrorCode PetscViewerBinarySetSkipOptions(PetscViewer viewer, PetscBool skip) {
5465c6c1daeSBarry Smith   PetscFunctionBegin;
547cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
548cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, skip, 2);
549cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetSkipOptions_C", (PetscViewer, PetscBool), (viewer, skip));
550807ea322SDave May   PetscFunctionReturn(0);
551807ea322SDave May }
552807ea322SDave May 
5539371c9d4SSatish Balay static PetscErrorCode PetscViewerBinarySetSkipOptions_Binary(PetscViewer viewer, PetscBool skip) {
554cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
555807ea322SDave May 
556807ea322SDave May   PetscFunctionBegin;
557cc843e7aSLisandro Dalcin   vbinary->skipoptions = skip;
5585c6c1daeSBarry Smith   PetscFunctionReturn(0);
5595c6c1daeSBarry Smith }
5605c6c1daeSBarry Smith 
5615c6c1daeSBarry Smith /*@
5625c6c1daeSBarry Smith     PetscViewerBinaryGetSkipOptions - checks if viewer uses the PETSc options database when loading objects
5635c6c1daeSBarry Smith 
5645c6c1daeSBarry Smith     Not Collective
5655c6c1daeSBarry Smith 
5665c6c1daeSBarry Smith     Input Parameter:
5675c6c1daeSBarry Smith .   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
5685c6c1daeSBarry Smith 
5695c6c1daeSBarry Smith     Output Parameter:
5705c6c1daeSBarry Smith .   skip - PETSC_TRUE means do not use
5715c6c1daeSBarry Smith 
5725c6c1daeSBarry Smith     Level: advanced
5735c6c1daeSBarry Smith 
57495452b02SPatrick Sanan     Notes:
57595452b02SPatrick Sanan     This must be called after PetscViewerSetType()
5765c6c1daeSBarry Smith 
577db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
578db781477SPatrick Sanan           `PetscViewerBinarySetSkipOptions()`
5795c6c1daeSBarry Smith @*/
5809371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryGetSkipOptions(PetscViewer viewer, PetscBool *skip) {
5815c6c1daeSBarry Smith   PetscFunctionBegin;
582cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
583cc843e7aSLisandro Dalcin   PetscValidBoolPointer(skip, 2);
584cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetSkipOptions_C", (PetscViewer, PetscBool *), (viewer, skip));
5855c6c1daeSBarry Smith   PetscFunctionReturn(0);
5865c6c1daeSBarry Smith }
5875c6c1daeSBarry Smith 
5889371c9d4SSatish Balay static PetscErrorCode PetscViewerBinaryGetSkipOptions_Binary(PetscViewer viewer, PetscBool *skip) {
589d21b9a37SPierre Jolivet   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
5905c6c1daeSBarry Smith 
5915c6c1daeSBarry Smith   PetscFunctionBegin;
592cc843e7aSLisandro Dalcin   *skip = vbinary->skipoptions;
5935c6c1daeSBarry Smith   PetscFunctionReturn(0);
5945c6c1daeSBarry Smith }
5955c6c1daeSBarry Smith 
5965c6c1daeSBarry Smith /*@
5975c6c1daeSBarry Smith     PetscViewerBinarySetSkipHeader - do not write a header with size information on output, just raw data
5985c6c1daeSBarry Smith 
5995c6c1daeSBarry Smith     Not Collective
6005c6c1daeSBarry Smith 
6015c6c1daeSBarry Smith     Input Parameters:
6025c6c1daeSBarry Smith +   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
6035c6c1daeSBarry Smith -   skip - PETSC_TRUE means do not write header
6045c6c1daeSBarry Smith 
6055c6c1daeSBarry Smith     Options Database Key:
60610699b91SBarry Smith .   -viewer_binary_skip_header - PETSC_TRUE means do not write header
6075c6c1daeSBarry Smith 
6085c6c1daeSBarry Smith     Level: advanced
6095c6c1daeSBarry Smith 
61095452b02SPatrick Sanan     Notes:
61195452b02SPatrick Sanan       This must be called after PetscViewerSetType()
6125c6c1daeSBarry Smith 
61310699b91SBarry Smith       Is ignored on anything but a binary viewer
6145c6c1daeSBarry Smith 
615db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
616db781477SPatrick Sanan           `PetscViewerBinaryGetSkipHeader()`
6175c6c1daeSBarry Smith @*/
6189371c9d4SSatish Balay PetscErrorCode PetscViewerBinarySetSkipHeader(PetscViewer viewer, PetscBool skip) {
6195c6c1daeSBarry Smith   PetscFunctionBegin;
620cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
621cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer, skip, 2);
622cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinarySetSkipHeader_C", (PetscViewer, PetscBool), (viewer, skip));
6235c6c1daeSBarry Smith   PetscFunctionReturn(0);
6245c6c1daeSBarry Smith }
6255c6c1daeSBarry Smith 
6269371c9d4SSatish Balay static PetscErrorCode PetscViewerBinarySetSkipHeader_Binary(PetscViewer viewer, PetscBool skip) {
6275c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
6285c6c1daeSBarry Smith 
6295c6c1daeSBarry Smith   PetscFunctionBegin;
630cc843e7aSLisandro Dalcin   vbinary->skipheader = skip;
6315c6c1daeSBarry Smith   PetscFunctionReturn(0);
6325c6c1daeSBarry Smith }
6335c6c1daeSBarry Smith 
6345c6c1daeSBarry Smith /*@
6355c6c1daeSBarry Smith     PetscViewerBinaryGetSkipHeader - checks whether to write a header with size information on output, or just raw data
6365c6c1daeSBarry Smith 
6375c6c1daeSBarry Smith     Not Collective
6385c6c1daeSBarry Smith 
6395c6c1daeSBarry Smith     Input Parameter:
6405c6c1daeSBarry Smith .   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
6415c6c1daeSBarry Smith 
6425c6c1daeSBarry Smith     Output Parameter:
6435c6c1daeSBarry Smith .   skip - PETSC_TRUE means do not write header
6445c6c1daeSBarry Smith 
6455c6c1daeSBarry Smith     Level: advanced
6465c6c1daeSBarry Smith 
64795452b02SPatrick Sanan     Notes:
64895452b02SPatrick Sanan     This must be called after PetscViewerSetType()
6495c6c1daeSBarry Smith 
6505c6c1daeSBarry Smith             Returns false for PETSCSOCKETVIEWER, you cannot skip the header for it.
6515c6c1daeSBarry Smith 
652db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
653db781477SPatrick Sanan           `PetscViewerBinarySetSkipHeader()`
6545c6c1daeSBarry Smith @*/
6559371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryGetSkipHeader(PetscViewer viewer, PetscBool *skip) {
6565c6c1daeSBarry Smith   PetscFunctionBegin;
657cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
658cc843e7aSLisandro Dalcin   PetscValidBoolPointer(skip, 2);
659cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerBinaryGetSkipHeader_C", (PetscViewer, PetscBool *), (viewer, skip));
6605c6c1daeSBarry Smith   PetscFunctionReturn(0);
6615c6c1daeSBarry Smith }
6625c6c1daeSBarry Smith 
6639371c9d4SSatish Balay static PetscErrorCode PetscViewerBinaryGetSkipHeader_Binary(PetscViewer viewer, PetscBool *skip) {
664f3b211e4SSatish Balay   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
6655c6c1daeSBarry Smith 
6665c6c1daeSBarry Smith   PetscFunctionBegin;
667cc843e7aSLisandro Dalcin   *skip = vbinary->skipheader;
6685c6c1daeSBarry Smith   PetscFunctionReturn(0);
6695c6c1daeSBarry Smith }
6705c6c1daeSBarry Smith 
6715c6c1daeSBarry Smith /*@C
6725c6c1daeSBarry Smith     PetscViewerBinaryGetInfoPointer - Extracts the file pointer for the ASCII
6735c6c1daeSBarry Smith           info file associated with a binary file.
6745c6c1daeSBarry Smith 
6755c6c1daeSBarry Smith     Not Collective
6765c6c1daeSBarry Smith 
6775c6c1daeSBarry Smith     Input Parameter:
6785c6c1daeSBarry Smith .   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
6795c6c1daeSBarry Smith 
6805c6c1daeSBarry Smith     Output Parameter:
6810298fd71SBarry Smith .   file - file pointer  Always returns NULL if not a binary viewer
6825c6c1daeSBarry Smith 
6835c6c1daeSBarry Smith     Level: advanced
6845c6c1daeSBarry Smith 
6855c6c1daeSBarry Smith     Notes:
6865c6c1daeSBarry Smith       For writable binary PetscViewers, the descriptor will only be valid for the
6875c6c1daeSBarry Smith     first processor in the communicator that shares the PetscViewer.
6885c6c1daeSBarry Smith 
6895c6c1daeSBarry Smith     Fortran Note:
6905c6c1daeSBarry Smith     This routine is not supported in Fortran.
6915c6c1daeSBarry Smith 
692c2e3fba1SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`
6935c6c1daeSBarry Smith @*/
6949371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryGetInfoPointer(PetscViewer viewer, FILE **file) {
6955c6c1daeSBarry Smith   PetscFunctionBegin;
696cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
697cc843e7aSLisandro Dalcin   PetscValidPointer(file, 2);
6980298fd71SBarry Smith   *file = NULL;
699cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerBinaryGetInfoPointer_C", (PetscViewer, FILE **), (viewer, file));
7005c6c1daeSBarry Smith   PetscFunctionReturn(0);
7015c6c1daeSBarry Smith }
7025c6c1daeSBarry Smith 
7039371c9d4SSatish Balay static PetscErrorCode PetscViewerBinaryGetInfoPointer_Binary(PetscViewer viewer, FILE **file) {
704cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
7055c6c1daeSBarry Smith 
7065c6c1daeSBarry Smith   PetscFunctionBegin;
7079566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
708cc843e7aSLisandro Dalcin   *file = vbinary->fdes_info;
709cc843e7aSLisandro Dalcin   if (viewer->format == PETSC_VIEWER_BINARY_MATLAB && !vbinary->matlabheaderwritten) {
7105c6c1daeSBarry Smith     if (vbinary->fdes_info) {
711cc843e7aSLisandro Dalcin       FILE *info = vbinary->fdes_info;
7129566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- begin code written by PetscViewerBinary for MATLAB format ---#\n"));
7139566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#$$ Set.filename = '%s';\n", vbinary->filename));
7149566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#$$ fd = PetscOpenFile(Set.filename);\n"));
7159566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n"));
716cc843e7aSLisandro Dalcin     }
717cc843e7aSLisandro Dalcin     vbinary->matlabheaderwritten = PETSC_TRUE;
7185c6c1daeSBarry Smith   }
7195c6c1daeSBarry Smith   PetscFunctionReturn(0);
7205c6c1daeSBarry Smith }
7215c6c1daeSBarry Smith 
722e0385b85SDave May #if defined(PETSC_HAVE_MPIIO)
7239371c9d4SSatish Balay static PetscErrorCode PetscViewerFileClose_BinaryMPIIO(PetscViewer v) {
724e0385b85SDave May   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
725e0385b85SDave May 
726e0385b85SDave May   PetscFunctionBegin;
727*48a46eb9SPierre Jolivet   if (vbinary->mfdes != MPI_FILE_NULL) PetscCallMPI(MPI_File_close(&vbinary->mfdes));
728*48a46eb9SPierre Jolivet   if (vbinary->mfsub != MPI_FILE_NULL) PetscCallMPI(MPI_File_close(&vbinary->mfsub));
729cc843e7aSLisandro Dalcin   vbinary->moff = 0;
730e0385b85SDave May   PetscFunctionReturn(0);
731e0385b85SDave May }
732e0385b85SDave May #endif
733e0385b85SDave May 
7349371c9d4SSatish Balay static PetscErrorCode PetscViewerFileClose_BinarySTDIO(PetscViewer v) {
735cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
736cc843e7aSLisandro Dalcin 
737cc843e7aSLisandro Dalcin   PetscFunctionBegin;
738cc843e7aSLisandro Dalcin   if (vbinary->fdes != -1) {
7399566063dSJacob Faibussowitsch     PetscCall(PetscBinaryClose(vbinary->fdes));
740cc843e7aSLisandro Dalcin     vbinary->fdes = -1;
741cc843e7aSLisandro Dalcin     if (vbinary->storecompressed) {
742cc843e7aSLisandro Dalcin       char        cmd[8 + PETSC_MAX_PATH_LEN], out[64 + PETSC_MAX_PATH_LEN] = "";
743cc843e7aSLisandro Dalcin       const char *gzfilename = vbinary->ogzfilename ? vbinary->ogzfilename : vbinary->filename;
744cc843e7aSLisandro Dalcin       /* compress the file */
7459566063dSJacob Faibussowitsch       PetscCall(PetscStrncpy(cmd, "gzip -f ", sizeof(cmd)));
7469566063dSJacob Faibussowitsch       PetscCall(PetscStrlcat(cmd, gzfilename, sizeof(cmd)));
747cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_POPEN)
748cc843e7aSLisandro Dalcin       {
749cc843e7aSLisandro Dalcin         FILE *fp;
7509566063dSJacob Faibussowitsch         PetscCall(PetscPOpen(PETSC_COMM_SELF, NULL, cmd, "r", &fp));
751cc73adaaSBarry Smith         PetscCheck(!fgets(out, (int)(sizeof(out) - 1), fp), PETSC_COMM_SELF, PETSC_ERR_LIB, "Error from command %s\n%s", cmd, out);
7529566063dSJacob Faibussowitsch         PetscCall(PetscPClose(PETSC_COMM_SELF, fp));
753cc843e7aSLisandro Dalcin       }
754cc843e7aSLisandro Dalcin #endif
755cc843e7aSLisandro Dalcin     }
756cc843e7aSLisandro Dalcin   }
7579566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary->ogzfilename));
758cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
759cc843e7aSLisandro Dalcin }
760cc843e7aSLisandro Dalcin 
7619371c9d4SSatish Balay static PetscErrorCode PetscViewerFileClose_BinaryInfo(PetscViewer v) {
762cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
763cc843e7aSLisandro Dalcin 
764cc843e7aSLisandro Dalcin   PetscFunctionBegin;
765cc843e7aSLisandro Dalcin   if (v->format == PETSC_VIEWER_BINARY_MATLAB && vbinary->matlabheaderwritten) {
766cc843e7aSLisandro Dalcin     if (vbinary->fdes_info) {
767cc843e7aSLisandro Dalcin       FILE *info = vbinary->fdes_info;
7689566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- begin code written by PetscViewerBinary for MATLAB format ---#\n"));
7699566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#$$ close(fd);\n"));
7709566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n"));
771cc843e7aSLisandro Dalcin     }
772cc843e7aSLisandro Dalcin   }
773cc843e7aSLisandro Dalcin   if (vbinary->fdes_info) {
774cc843e7aSLisandro Dalcin     FILE *info         = vbinary->fdes_info;
775cc843e7aSLisandro Dalcin     vbinary->fdes_info = NULL;
776cc73adaaSBarry Smith     PetscCheck(!fclose(info), PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file");
777cc843e7aSLisandro Dalcin   }
778cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
779cc843e7aSLisandro Dalcin }
780cc843e7aSLisandro Dalcin 
7819371c9d4SSatish Balay static PetscErrorCode PetscViewerFileClose_Binary(PetscViewer v) {
782cc843e7aSLisandro Dalcin   PetscFunctionBegin;
783cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
7849566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_BinaryMPIIO(v));
785cc843e7aSLisandro Dalcin #endif
7869566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_BinarySTDIO(v));
7879566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_BinaryInfo(v));
788cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
789cc843e7aSLisandro Dalcin }
790cc843e7aSLisandro Dalcin 
7919371c9d4SSatish Balay static PetscErrorCode PetscViewerDestroy_Binary(PetscViewer v) {
7925c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
7935c6c1daeSBarry Smith 
7945c6c1daeSBarry Smith   PetscFunctionBegin;
7959566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_Binary(v));
7969566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary->filename));
7979566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary));
7982e956fe4SStefano Zampini   PetscCall(PetscViewerBinaryClearFunctionList(v));
799e0385b85SDave May   PetscFunctionReturn(0);
800e0385b85SDave May }
8015c6c1daeSBarry Smith 
8025c6c1daeSBarry Smith /*@C
8035c6c1daeSBarry Smith    PetscViewerBinaryOpen - Opens a file for binary input/output.
8045c6c1daeSBarry Smith 
805d083f849SBarry Smith    Collective
8065c6c1daeSBarry Smith 
8075c6c1daeSBarry Smith    Input Parameters:
8085c6c1daeSBarry Smith +  comm - MPI communicator
8095c6c1daeSBarry Smith .  name - name of file
810cc843e7aSLisandro Dalcin -  mode - open mode of file
8115c6c1daeSBarry Smith $    FILE_MODE_WRITE - create new file for binary output
8125c6c1daeSBarry Smith $    FILE_MODE_READ - open existing file for binary input
8135c6c1daeSBarry Smith $    FILE_MODE_APPEND - open existing file for binary output
8145c6c1daeSBarry Smith 
8155c6c1daeSBarry Smith    Output Parameter:
816cc843e7aSLisandro Dalcin .  viewer - PetscViewer for binary input/output to use with the specified file
8175c6c1daeSBarry Smith 
8185c6c1daeSBarry Smith     Options Database Keys:
81963c55180SPatrick Sanan +    -viewer_binary_filename <name> -
82063c55180SPatrick Sanan .    -viewer_binary_skip_info -
82163c55180SPatrick Sanan .    -viewer_binary_skip_options -
82263c55180SPatrick Sanan .    -viewer_binary_skip_header -
82363c55180SPatrick Sanan -    -viewer_binary_mpiio -
8245c6c1daeSBarry Smith 
8255c6c1daeSBarry Smith    Level: beginner
8265c6c1daeSBarry Smith 
8275c6c1daeSBarry Smith    Note:
8285c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
8295c6c1daeSBarry Smith 
8305c6c1daeSBarry Smith     For reading files, the filename may begin with ftp:// or http:// and/or
8315c6c1daeSBarry Smith     end with .gz; in this case file is brought over and uncompressed.
8325c6c1daeSBarry Smith 
8335c6c1daeSBarry Smith     For creating files, if the file name ends with .gz it is automatically
8345c6c1daeSBarry Smith     compressed when closed.
8355c6c1daeSBarry Smith 
836db781477SPatrick Sanan .seealso: `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
837db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
838db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`, `PetscViewerBinarySetUseMPIIO()`,
839db781477SPatrick Sanan           `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
8405c6c1daeSBarry Smith @*/
8419371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryOpen(MPI_Comm comm, const char name[], PetscFileMode mode, PetscViewer *viewer) {
8425c6c1daeSBarry Smith   PetscFunctionBegin;
8439566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(comm, viewer));
8449566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(*viewer, PETSCVIEWERBINARY));
8459566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(*viewer, mode));
8469566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetName(*viewer, name));
8479566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions(*viewer));
8485c6c1daeSBarry Smith   PetscFunctionReturn(0);
8495c6c1daeSBarry Smith }
8505c6c1daeSBarry Smith 
8515c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
8529371c9d4SSatish Balay static PetscErrorCode PetscViewerBinaryWriteReadMPIIO(PetscViewer viewer, void *data, PetscInt num, PetscInt *count, PetscDataType dtype, PetscBool write) {
85322a8f86cSLisandro Dalcin   MPI_Comm            comm    = PetscObjectComm((PetscObject)viewer);
8545c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
85522a8f86cSLisandro Dalcin   MPI_File            mfdes   = vbinary->mfdes;
8565c6c1daeSBarry Smith   MPI_Datatype        mdtype;
85769e10bbaSLisandro Dalcin   PetscMPIInt         rank, cnt;
8585c6c1daeSBarry Smith   MPI_Status          status;
8595c6c1daeSBarry Smith   MPI_Aint            ul, dsize;
8605c6c1daeSBarry Smith 
8615c6c1daeSBarry Smith   PetscFunctionBegin;
8629566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
8639566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(num, &cnt));
8649566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToMPIDataType(dtype, &mdtype));
8655c6c1daeSBarry Smith   if (write) {
866*48a46eb9SPierre Jolivet     if (rank == 0) PetscCall(MPIU_File_write_at(mfdes, vbinary->moff, data, cnt, mdtype, &status));
8675c6c1daeSBarry Smith   } else {
868dd400576SPatrick Sanan     if (rank == 0) {
8699566063dSJacob Faibussowitsch       PetscCall(MPIU_File_read_at(mfdes, vbinary->moff, data, cnt, mdtype, &status));
8709566063dSJacob Faibussowitsch       if (cnt > 0) PetscCallMPI(MPI_Get_count(&status, mdtype, &cnt));
8715c6c1daeSBarry Smith     }
8729566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&cnt, 1, MPI_INT, 0, comm));
8739566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(data, cnt, mdtype, 0, comm));
87469e10bbaSLisandro Dalcin   }
8759566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_get_extent(mdtype, &ul, &dsize));
8765c6c1daeSBarry Smith   vbinary->moff += dsize * cnt;
8779860990eSLisandro Dalcin   if (count) *count = cnt;
8785c6c1daeSBarry Smith   PetscFunctionReturn(0);
8795c6c1daeSBarry Smith }
8805c6c1daeSBarry Smith #endif
8815c6c1daeSBarry Smith 
8825c6c1daeSBarry Smith /*@C
8835c6c1daeSBarry Smith    PetscViewerBinaryRead - Reads from a binary file, all processors get the same result
8845c6c1daeSBarry Smith 
885d083f849SBarry Smith    Collective
8865c6c1daeSBarry Smith 
8875c6c1daeSBarry Smith    Input Parameters:
8885c6c1daeSBarry Smith +  viewer - the binary viewer
889907376e6SBarry Smith .  data - location of the data to be written
890060da220SMatthew G. Knepley .  num - number of items of data to read
891907376e6SBarry Smith -  dtype - type of data to read
8925c6c1daeSBarry Smith 
893f8e4bde8SMatthew G. Knepley    Output Parameters:
8945972f5f3SLisandro Dalcin .  count - number of items of data actually read, or NULL.
895f8e4bde8SMatthew G. Knepley 
8965c6c1daeSBarry Smith    Level: beginner
8975c6c1daeSBarry Smith 
898db781477SPatrick Sanan .seealso: `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
899db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
900db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
9015c6c1daeSBarry Smith @*/
9029371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryRead(PetscViewer viewer, void *data, PetscInt num, PetscInt *count, PetscDataType dtype) {
90322a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
9045c6c1daeSBarry Smith 
90522a8f86cSLisandro Dalcin   PetscFunctionBegin;
90622a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
90722a8f86cSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, num, 3);
9089566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
90922a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
9105c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
911bc196f7cSDave May   if (vbinary->usempiio) {
9129566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWriteReadMPIIO(viewer, data, num, count, dtype, PETSC_FALSE));
9135c6c1daeSBarry Smith   } else {
9145c6c1daeSBarry Smith #endif
9159566063dSJacob Faibussowitsch     PetscCall(PetscBinarySynchronizedRead(PetscObjectComm((PetscObject)viewer), vbinary->fdes, data, num, count, dtype));
9165c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
9175c6c1daeSBarry Smith   }
9185c6c1daeSBarry Smith #endif
9195c6c1daeSBarry Smith   PetscFunctionReturn(0);
9205c6c1daeSBarry Smith }
9215c6c1daeSBarry Smith 
9225c6c1daeSBarry Smith /*@C
9235c6c1daeSBarry Smith    PetscViewerBinaryWrite - writes to a binary file, only from the first process
9245c6c1daeSBarry Smith 
925d083f849SBarry Smith    Collective
9265c6c1daeSBarry Smith 
9275c6c1daeSBarry Smith    Input Parameters:
9285c6c1daeSBarry Smith +  viewer - the binary viewer
9295c6c1daeSBarry Smith .  data - location of data
9305824c9d0SPatrick Sanan .  count - number of items of data to write
931f253e43cSLisandro Dalcin -  dtype - type of data to write
9325c6c1daeSBarry Smith 
9335c6c1daeSBarry Smith    Level: beginner
9345c6c1daeSBarry Smith 
935db781477SPatrick Sanan .seealso: `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
936db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`, `PetscDataType`
937db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
9385c6c1daeSBarry Smith @*/
9399371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryWrite(PetscViewer viewer, const void *data, PetscInt count, PetscDataType dtype) {
94022a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
9415c6c1daeSBarry Smith 
9425c6c1daeSBarry Smith   PetscFunctionBegin;
94322a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
94422a8f86cSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer, count, 3);
9459566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
94622a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary *)viewer->data;
9475c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
948bc196f7cSDave May   if (vbinary->usempiio) {
9499566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWriteReadMPIIO(viewer, (void *)data, count, NULL, dtype, PETSC_TRUE));
9505c6c1daeSBarry Smith   } else {
9515c6c1daeSBarry Smith #endif
9529566063dSJacob Faibussowitsch     PetscCall(PetscBinarySynchronizedWrite(PetscObjectComm((PetscObject)viewer), vbinary->fdes, data, count, dtype));
9535c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
9545c6c1daeSBarry Smith   }
9555c6c1daeSBarry Smith #endif
9565c6c1daeSBarry Smith   PetscFunctionReturn(0);
9575c6c1daeSBarry Smith }
9585c6c1daeSBarry Smith 
9599371c9d4SSatish Balay static PetscErrorCode PetscViewerBinaryWriteReadAll(PetscViewer viewer, PetscBool write, void *data, PetscInt count, PetscInt start, PetscInt total, PetscDataType dtype) {
9605972f5f3SLisandro Dalcin   MPI_Comm              comm = PetscObjectComm((PetscObject)viewer);
9615972f5f3SLisandro Dalcin   PetscMPIInt           size, rank;
9625972f5f3SLisandro Dalcin   MPI_Datatype          mdtype;
963ec4bef21SJose E. Roman   PETSC_UNUSED MPI_Aint lb;
964ec4bef21SJose E. Roman   MPI_Aint              dsize;
9655972f5f3SLisandro Dalcin   PetscBool             useMPIIO;
9665972f5f3SLisandro Dalcin 
9675972f5f3SLisandro Dalcin   PetscFunctionBegin;
9685972f5f3SLisandro Dalcin   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
969064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveBool(viewer, ((start >= 0) || (start == PETSC_DETERMINE)), 5);
970064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveBool(viewer, ((total >= 0) || (total == PETSC_DETERMINE)), 6);
971064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveInt(viewer, total, 6);
9729566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
9735972f5f3SLisandro Dalcin 
9749566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToMPIDataType(dtype, &mdtype));
9759566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_get_extent(mdtype, &lb, &dsize));
9769566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
9779566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
9785972f5f3SLisandro Dalcin 
9799566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetUseMPIIO(viewer, &useMPIIO));
9805972f5f3SLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
9815972f5f3SLisandro Dalcin   if (useMPIIO) {
9825972f5f3SLisandro Dalcin     MPI_File    mfdes;
9835972f5f3SLisandro Dalcin     MPI_Offset  off;
9845972f5f3SLisandro Dalcin     PetscMPIInt cnt;
9855972f5f3SLisandro Dalcin 
9865972f5f3SLisandro Dalcin     if (start == PETSC_DETERMINE) {
9879566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Scan(&count, &start, 1, MPIU_INT, MPI_SUM, comm));
9885972f5f3SLisandro Dalcin       start -= count;
9895972f5f3SLisandro Dalcin     }
9905972f5f3SLisandro Dalcin     if (total == PETSC_DETERMINE) {
9915972f5f3SLisandro Dalcin       total = start + count;
9929566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Bcast(&total, 1, MPIU_INT, size - 1, comm));
9935972f5f3SLisandro Dalcin     }
9949566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(count, &cnt));
9959566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryGetMPIIODescriptor(viewer, &mfdes));
9969566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryGetMPIIOOffset(viewer, &off));
9975972f5f3SLisandro Dalcin     off += (MPI_Offset)(start * dsize);
9985972f5f3SLisandro Dalcin     if (write) {
9999566063dSJacob Faibussowitsch       PetscCall(MPIU_File_write_at_all(mfdes, off, data, cnt, mdtype, MPI_STATUS_IGNORE));
10005972f5f3SLisandro Dalcin     } else {
10019566063dSJacob Faibussowitsch       PetscCall(MPIU_File_read_at_all(mfdes, off, data, cnt, mdtype, MPI_STATUS_IGNORE));
10025972f5f3SLisandro Dalcin     }
10035972f5f3SLisandro Dalcin     off = (MPI_Offset)(total * dsize);
10049566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryAddMPIIOOffset(viewer, off));
10055972f5f3SLisandro Dalcin     PetscFunctionReturn(0);
10065972f5f3SLisandro Dalcin   }
10075972f5f3SLisandro Dalcin #endif
10085972f5f3SLisandro Dalcin   {
10095972f5f3SLisandro Dalcin     int         fdes;
10105972f5f3SLisandro Dalcin     char       *workbuf = NULL;
1011dd400576SPatrick Sanan     PetscInt    tcount = rank == 0 ? 0 : count, maxcount = 0, message_count, flowcontrolcount;
10125972f5f3SLisandro Dalcin     PetscMPIInt tag, cnt, maxcnt, scnt = 0, rcnt = 0, j;
10135972f5f3SLisandro Dalcin     MPI_Status  status;
10145972f5f3SLisandro Dalcin 
10159566063dSJacob Faibussowitsch     PetscCall(PetscCommGetNewTag(comm, &tag));
10169566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Reduce(&tcount, &maxcount, 1, MPIU_INT, MPI_MAX, 0, comm));
10179566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(maxcount, &maxcnt));
10189566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(count, &cnt));
10195972f5f3SLisandro Dalcin 
10209566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryGetDescriptor(viewer, &fdes));
10219566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlowControlStart(viewer, &message_count, &flowcontrolcount));
1022dd400576SPatrick Sanan     if (rank == 0) {
10239566063dSJacob Faibussowitsch       PetscCall(PetscMalloc(maxcnt * dsize, &workbuf));
10245972f5f3SLisandro Dalcin       if (write) {
10259566063dSJacob Faibussowitsch         PetscCall(PetscBinaryWrite(fdes, data, cnt, dtype));
10265972f5f3SLisandro Dalcin       } else {
10279566063dSJacob Faibussowitsch         PetscCall(PetscBinaryRead(fdes, data, cnt, NULL, dtype));
10285972f5f3SLisandro Dalcin       }
10295972f5f3SLisandro Dalcin       for (j = 1; j < size; j++) {
10309566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlowControlStepMain(viewer, j, &message_count, flowcontrolcount));
10315972f5f3SLisandro Dalcin         if (write) {
10329566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Recv(workbuf, maxcnt, mdtype, j, tag, comm, &status));
10339566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Get_count(&status, mdtype, &rcnt));
10349566063dSJacob Faibussowitsch           PetscCall(PetscBinaryWrite(fdes, workbuf, rcnt, dtype));
10355972f5f3SLisandro Dalcin         } else {
10369566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Recv(&scnt, 1, MPI_INT, j, tag, comm, MPI_STATUS_IGNORE));
10379566063dSJacob Faibussowitsch           PetscCall(PetscBinaryRead(fdes, workbuf, scnt, NULL, dtype));
10389566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Send(workbuf, scnt, mdtype, j, tag, comm));
10395972f5f3SLisandro Dalcin         }
10405972f5f3SLisandro Dalcin       }
10419566063dSJacob Faibussowitsch       PetscCall(PetscFree(workbuf));
10429566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlowControlEndMain(viewer, &message_count));
10435972f5f3SLisandro Dalcin     } else {
10449566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlowControlStepWorker(viewer, rank, &message_count));
10455972f5f3SLisandro Dalcin       if (write) {
10469566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(data, cnt, mdtype, 0, tag, comm));
10475972f5f3SLisandro Dalcin       } else {
10489566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(&cnt, 1, MPI_INT, 0, tag, comm));
10499566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Recv(data, cnt, mdtype, 0, tag, comm, MPI_STATUS_IGNORE));
10505972f5f3SLisandro Dalcin       }
10519566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlowControlEndWorker(viewer, &message_count));
10525972f5f3SLisandro Dalcin     }
10535972f5f3SLisandro Dalcin   }
10545972f5f3SLisandro Dalcin   PetscFunctionReturn(0);
10555972f5f3SLisandro Dalcin }
10565972f5f3SLisandro Dalcin 
10575972f5f3SLisandro Dalcin /*@C
10585972f5f3SLisandro Dalcin    PetscViewerBinaryReadAll - reads from a binary file from all processes
10595972f5f3SLisandro Dalcin 
10605972f5f3SLisandro Dalcin    Collective
10615972f5f3SLisandro Dalcin 
10625972f5f3SLisandro Dalcin    Input Parameters:
10635972f5f3SLisandro Dalcin +  viewer - the binary viewer
10645972f5f3SLisandro Dalcin .  data - location of data
10655972f5f3SLisandro Dalcin .  count - local number of items of data to read
10665972f5f3SLisandro Dalcin .  start - local start, can be PETSC_DETERMINE
10675972f5f3SLisandro Dalcin .  total - global number of items of data to read, can be PETSC_DETERMINE
10685972f5f3SLisandro Dalcin -  dtype - type of data to read
10695972f5f3SLisandro Dalcin 
10705972f5f3SLisandro Dalcin    Level: advanced
10715972f5f3SLisandro Dalcin 
1072db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryRead()`, `PetscViewerBinaryWriteAll()`
10735972f5f3SLisandro Dalcin @*/
10749371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryReadAll(PetscViewer viewer, void *data, PetscInt count, PetscInt start, PetscInt total, PetscDataType dtype) {
10755972f5f3SLisandro Dalcin   PetscFunctionBegin;
10769566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWriteReadAll(viewer, PETSC_FALSE, data, count, start, total, dtype));
10775972f5f3SLisandro Dalcin   PetscFunctionReturn(0);
10785972f5f3SLisandro Dalcin }
10795972f5f3SLisandro Dalcin 
10805972f5f3SLisandro Dalcin /*@C
10815972f5f3SLisandro Dalcin    PetscViewerBinaryWriteAll - writes to a binary file from all processes
10825972f5f3SLisandro Dalcin 
10835972f5f3SLisandro Dalcin    Collective
10845972f5f3SLisandro Dalcin 
10855972f5f3SLisandro Dalcin    Input Parameters:
10865972f5f3SLisandro Dalcin +  viewer - the binary viewer
10875972f5f3SLisandro Dalcin .  data - location of data
10885972f5f3SLisandro Dalcin .  count - local number of items of data to write
10895972f5f3SLisandro Dalcin .  start - local start, can be PETSC_DETERMINE
10905972f5f3SLisandro Dalcin .  total - global number of items of data to write, can be PETSC_DETERMINE
10915972f5f3SLisandro Dalcin -  dtype - type of data to write
10925972f5f3SLisandro Dalcin 
10935972f5f3SLisandro Dalcin    Level: advanced
10945972f5f3SLisandro Dalcin 
1095db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryWriteAll()`, `PetscViewerBinaryReadAll()`
10965972f5f3SLisandro Dalcin @*/
10979371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryWriteAll(PetscViewer viewer, const void *data, PetscInt count, PetscInt start, PetscInt total, PetscDataType dtype) {
10985972f5f3SLisandro Dalcin   PetscFunctionBegin;
10999566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWriteReadAll(viewer, PETSC_TRUE, (void *)data, count, start, total, dtype));
11005972f5f3SLisandro Dalcin   PetscFunctionReturn(0);
11015972f5f3SLisandro Dalcin }
11025972f5f3SLisandro Dalcin 
11035c6c1daeSBarry Smith /*@C
11045c6c1daeSBarry Smith    PetscViewerBinaryWriteStringArray - writes to a binary file, only from the first process an array of strings
11055c6c1daeSBarry Smith 
1106d083f849SBarry Smith    Collective
11075c6c1daeSBarry Smith 
11085c6c1daeSBarry Smith    Input Parameters:
11095c6c1daeSBarry Smith +  viewer - the binary viewer
11105c6c1daeSBarry Smith -  data - location of the array of strings
11115c6c1daeSBarry Smith 
11125c6c1daeSBarry Smith    Level: intermediate
11135c6c1daeSBarry Smith 
111495452b02SPatrick Sanan     Notes:
111595452b02SPatrick Sanan     array of strings is null terminated
11165c6c1daeSBarry Smith 
1117db781477SPatrick Sanan .seealso: `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
1118db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
1119db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
11205c6c1daeSBarry Smith @*/
11219371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryWriteStringArray(PetscViewer viewer, const char *const *data) {
11225c6c1daeSBarry Smith   PetscInt i, n = 0, *sizes;
1123cc843e7aSLisandro Dalcin   size_t   len;
11245c6c1daeSBarry Smith 
112522a8f86cSLisandro Dalcin   PetscFunctionBegin;
11269566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
11275c6c1daeSBarry Smith   /* count number of strings */
11289371c9d4SSatish Balay   while (data[n++])
11299371c9d4SSatish Balay     ;
11305c6c1daeSBarry Smith   n--;
11319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n + 1, &sizes));
11325c6c1daeSBarry Smith   sizes[0] = n;
11335c6c1daeSBarry Smith   for (i = 0; i < n; i++) {
11349566063dSJacob Faibussowitsch     PetscCall(PetscStrlen(data[i], &len));
1135cc843e7aSLisandro Dalcin     sizes[i + 1] = (PetscInt)len + 1; /* size includes space for the null terminator */
11365c6c1daeSBarry Smith   }
11379566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWrite(viewer, sizes, n + 1, PETSC_INT));
1138*48a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall(PetscViewerBinaryWrite(viewer, (void *)data[i], sizes[i + 1], PETSC_CHAR));
11399566063dSJacob Faibussowitsch   PetscCall(PetscFree(sizes));
11405c6c1daeSBarry Smith   PetscFunctionReturn(0);
11415c6c1daeSBarry Smith }
11425c6c1daeSBarry Smith 
11435c6c1daeSBarry Smith /*@C
11445c6c1daeSBarry Smith    PetscViewerBinaryReadStringArray - reads a binary file an array of strings
11455c6c1daeSBarry Smith 
1146d083f849SBarry Smith    Collective
11475c6c1daeSBarry Smith 
11485c6c1daeSBarry Smith    Input Parameter:
11495c6c1daeSBarry Smith .  viewer - the binary viewer
11505c6c1daeSBarry Smith 
11515c6c1daeSBarry Smith    Output Parameter:
11525c6c1daeSBarry Smith .  data - location of the array of strings
11535c6c1daeSBarry Smith 
11545c6c1daeSBarry Smith    Level: intermediate
11555c6c1daeSBarry Smith 
115695452b02SPatrick Sanan     Notes:
115795452b02SPatrick Sanan     array of strings is null terminated
11585c6c1daeSBarry Smith 
1159db781477SPatrick Sanan .seealso: `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
1160db781477SPatrick Sanan           `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
1161db781477SPatrick Sanan           `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
11625c6c1daeSBarry Smith @*/
11639371c9d4SSatish Balay PetscErrorCode PetscViewerBinaryReadStringArray(PetscViewer viewer, char ***data) {
1164060da220SMatthew G. Knepley   PetscInt i, n, *sizes, N = 0;
11655c6c1daeSBarry Smith 
116622a8f86cSLisandro Dalcin   PetscFunctionBegin;
11679566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
11685c6c1daeSBarry Smith   /* count number of strings */
11699566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &n, 1, NULL, PETSC_INT));
11709566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &sizes));
11719566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, sizes, n, NULL, PETSC_INT));
1172a297a907SKarl Rupp   for (i = 0; i < n; i++) N += sizes[i];
11739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc((n + 1) * sizeof(char *) + N * sizeof(char), data));
11745c6c1daeSBarry Smith   (*data)[0] = (char *)((*data) + n + 1);
1175a297a907SKarl Rupp   for (i = 1; i < n; i++) (*data)[i] = (*data)[i - 1] + sizes[i - 1];
11769566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, (*data)[0], N, NULL, PETSC_CHAR));
117702c9f0b5SLisandro Dalcin   (*data)[n] = NULL;
11789566063dSJacob Faibussowitsch   PetscCall(PetscFree(sizes));
11795c6c1daeSBarry Smith   PetscFunctionReturn(0);
11805c6c1daeSBarry Smith }
11815c6c1daeSBarry Smith 
1182cc843e7aSLisandro Dalcin /*@C
1183cc843e7aSLisandro Dalcin      PetscViewerFileSetMode - Sets the open mode of file
1184cc843e7aSLisandro Dalcin 
1185cc843e7aSLisandro Dalcin     Logically Collective on PetscViewer
1186cc843e7aSLisandro Dalcin 
1187cc843e7aSLisandro Dalcin   Input Parameters:
1188cc843e7aSLisandro Dalcin +  viewer - the PetscViewer; must be a a PETSCVIEWERBINARY, PETSCVIEWERMATLAB, PETSCVIEWERHDF5, or PETSCVIEWERASCII  PetscViewer
1189cc843e7aSLisandro Dalcin -  mode - open mode of file
1190cc843e7aSLisandro Dalcin $    FILE_MODE_WRITE - create new file for output
1191cc843e7aSLisandro Dalcin $    FILE_MODE_READ - open existing file for input
1192cc843e7aSLisandro Dalcin $    FILE_MODE_APPEND - open existing file for output
1193cc843e7aSLisandro Dalcin 
1194cc843e7aSLisandro Dalcin   Level: advanced
1195cc843e7aSLisandro Dalcin 
1196db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
1197cc843e7aSLisandro Dalcin 
1198cc843e7aSLisandro Dalcin @*/
11999371c9d4SSatish Balay PetscErrorCode PetscViewerFileSetMode(PetscViewer viewer, PetscFileMode mode) {
1200cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1201cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1202cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveEnum(viewer, mode, 2);
120308401ef6SPierre Jolivet   PetscCheck(mode != FILE_MODE_UNDEFINED, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Cannot set FILE_MODE_UNDEFINED");
1204f7d195e4SLawrence Mitchell   PetscCheck(mode >= FILE_MODE_UNDEFINED && mode <= FILE_MODE_APPEND_UPDATE, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_OUTOFRANGE, "Invalid file mode %d", (int)mode);
1205cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerFileSetMode_C", (PetscViewer, PetscFileMode), (viewer, mode));
1206cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1207cc843e7aSLisandro Dalcin }
1208cc843e7aSLisandro Dalcin 
12099371c9d4SSatish Balay static PetscErrorCode PetscViewerFileSetMode_Binary(PetscViewer viewer, PetscFileMode mode) {
1210cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1211cc843e7aSLisandro Dalcin 
1212cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1213cc73adaaSBarry Smith   PetscCheck(!viewer->setupcalled || vbinary->filemode == mode, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Cannot change mode to %s after setup", PetscFileModes[mode]);
1214cc843e7aSLisandro Dalcin   vbinary->filemode = mode;
1215cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1216cc843e7aSLisandro Dalcin }
1217cc843e7aSLisandro Dalcin 
1218cc843e7aSLisandro Dalcin /*@C
1219cc843e7aSLisandro Dalcin      PetscViewerFileGetMode - Gets the open mode of file
1220cc843e7aSLisandro Dalcin 
1221cc843e7aSLisandro Dalcin     Not Collective
1222cc843e7aSLisandro Dalcin 
1223cc843e7aSLisandro Dalcin   Input Parameter:
1224cc843e7aSLisandro Dalcin .  viewer - the PetscViewer; must be a PETSCVIEWERBINARY, PETSCVIEWERMATLAB, PETSCVIEWERHDF5, or PETSCVIEWERASCII  PetscViewer
1225cc843e7aSLisandro Dalcin 
1226cc843e7aSLisandro Dalcin   Output Parameter:
1227cc843e7aSLisandro Dalcin .  mode - open mode of file
1228cc843e7aSLisandro Dalcin $    FILE_MODE_WRITE - create new file for binary output
1229cc843e7aSLisandro Dalcin $    FILE_MODE_READ - open existing file for binary input
1230cc843e7aSLisandro Dalcin $    FILE_MODE_APPEND - open existing file for binary output
1231cc843e7aSLisandro Dalcin 
1232cc843e7aSLisandro Dalcin   Level: advanced
1233cc843e7aSLisandro Dalcin 
1234db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
1235cc843e7aSLisandro Dalcin 
1236cc843e7aSLisandro Dalcin @*/
12379371c9d4SSatish Balay PetscErrorCode PetscViewerFileGetMode(PetscViewer viewer, PetscFileMode *mode) {
1238cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1239cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1240cc843e7aSLisandro Dalcin   PetscValidPointer(mode, 2);
1241cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerFileGetMode_C", (PetscViewer, PetscFileMode *), (viewer, mode));
1242cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1243cc843e7aSLisandro Dalcin }
1244cc843e7aSLisandro Dalcin 
12459371c9d4SSatish Balay static PetscErrorCode PetscViewerFileGetMode_Binary(PetscViewer viewer, PetscFileMode *mode) {
1246cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1247cc843e7aSLisandro Dalcin 
1248cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1249cc843e7aSLisandro Dalcin   *mode = vbinary->filemode;
1250cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1251cc843e7aSLisandro Dalcin }
1252cc843e7aSLisandro Dalcin 
12539371c9d4SSatish Balay static PetscErrorCode PetscViewerFileSetName_Binary(PetscViewer viewer, const char name[]) {
1254cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1255cc843e7aSLisandro Dalcin 
1256cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1257cc843e7aSLisandro Dalcin   if (viewer->setupcalled && vbinary->filename) {
1258cc843e7aSLisandro Dalcin     /* gzip can be run after the file with the previous filename has been closed */
12599566063dSJacob Faibussowitsch     PetscCall(PetscFree(vbinary->ogzfilename));
12609566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(vbinary->filename, &vbinary->ogzfilename));
1261cc843e7aSLisandro Dalcin   }
12629566063dSJacob Faibussowitsch   PetscCall(PetscFree(vbinary->filename));
12639566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &vbinary->filename));
1264cc843e7aSLisandro Dalcin   viewer->setupcalled = PETSC_FALSE;
1265cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1266cc843e7aSLisandro Dalcin }
1267cc843e7aSLisandro Dalcin 
12689371c9d4SSatish Balay static PetscErrorCode PetscViewerFileGetName_Binary(PetscViewer viewer, const char **name) {
12695c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
12705c6c1daeSBarry Smith 
12715c6c1daeSBarry Smith   PetscFunctionBegin;
12725c6c1daeSBarry Smith   *name = vbinary->filename;
12735c6c1daeSBarry Smith   PetscFunctionReturn(0);
12745c6c1daeSBarry Smith }
12755c6c1daeSBarry Smith 
12765c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
12779371c9d4SSatish Balay static PetscErrorCode PetscViewerFileSetUp_BinaryMPIIO(PetscViewer viewer) {
12785c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1279e8a65b7dSLisandro Dalcin   int                 amode;
12805c6c1daeSBarry Smith 
12815c6c1daeSBarry Smith   PetscFunctionBegin;
1282a297a907SKarl Rupp   vbinary->storecompressed = PETSC_FALSE;
12835c6c1daeSBarry Smith 
1284cc843e7aSLisandro Dalcin   vbinary->moff = 0;
1285cc843e7aSLisandro Dalcin   switch (vbinary->filemode) {
1286e8a65b7dSLisandro Dalcin   case FILE_MODE_READ: amode = MPI_MODE_RDONLY; break;
1287e8a65b7dSLisandro Dalcin   case FILE_MODE_WRITE: amode = MPI_MODE_WRONLY | MPI_MODE_CREATE; break;
128822a8f86cSLisandro Dalcin   case FILE_MODE_APPEND: amode = MPI_MODE_WRONLY | MPI_MODE_CREATE | MPI_MODE_APPEND; break;
12897e4fd573SVaclav Hapla   case FILE_MODE_UNDEFINED: SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerSetUp()");
129098921bdaSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[vbinary->filemode]);
12915c6c1daeSBarry Smith   }
12929566063dSJacob Faibussowitsch   PetscCallMPI(MPI_File_open(PetscObjectComm((PetscObject)viewer), vbinary->filename, amode, MPI_INFO_NULL, &vbinary->mfdes));
129322a8f86cSLisandro Dalcin   /*
129422a8f86cSLisandro Dalcin       The MPI standard does not have MPI_MODE_TRUNCATE. We emulate this behavior by setting the file size to zero.
129522a8f86cSLisandro Dalcin   */
12969566063dSJacob Faibussowitsch   if (vbinary->filemode == FILE_MODE_WRITE) PetscCallMPI(MPI_File_set_size(vbinary->mfdes, 0));
129722a8f86cSLisandro Dalcin   /*
129822a8f86cSLisandro Dalcin       Initially, all processes view the file as a linear byte stream. Therefore, for files opened with MPI_MODE_APPEND,
129922a8f86cSLisandro Dalcin       MPI_File_get_position[_shared](fh, &offset) returns the absolute byte position at the end of file.
130022a8f86cSLisandro Dalcin       Otherwise, we would need to call MPI_File_get_byte_offset(fh, offset, &byte_offset) to convert
130122a8f86cSLisandro Dalcin       the offset in etype units to an absolute byte position.
130222a8f86cSLisandro Dalcin    */
13039566063dSJacob Faibussowitsch   if (vbinary->filemode == FILE_MODE_APPEND) PetscCallMPI(MPI_File_get_position(vbinary->mfdes, &vbinary->moff));
1304cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1305cc843e7aSLisandro Dalcin }
1306cc843e7aSLisandro Dalcin #endif
13075c6c1daeSBarry Smith 
13089371c9d4SSatish Balay static PetscErrorCode PetscViewerFileSetUp_BinarySTDIO(PetscViewer viewer) {
1309cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1310cc843e7aSLisandro Dalcin   const char         *fname;
1311cc843e7aSLisandro Dalcin   char                bname[PETSC_MAX_PATH_LEN], *gz;
1312cc843e7aSLisandro Dalcin   PetscBool           found;
1313cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
13145c6c1daeSBarry Smith 
1315cc843e7aSLisandro Dalcin   PetscFunctionBegin;
13169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
13175c6c1daeSBarry Smith 
1318cc843e7aSLisandro Dalcin   /* if file name ends in .gz strip that off and note user wants file compressed */
1319cc843e7aSLisandro Dalcin   vbinary->storecompressed = PETSC_FALSE;
1320cc843e7aSLisandro Dalcin   if (vbinary->filemode == FILE_MODE_WRITE) {
13219566063dSJacob Faibussowitsch     PetscCall(PetscStrstr(vbinary->filename, ".gz", &gz));
13229371c9d4SSatish Balay     if (gz && gz[3] == 0) {
13239371c9d4SSatish Balay       *gz                      = 0;
13249371c9d4SSatish Balay       vbinary->storecompressed = PETSC_TRUE;
13259371c9d4SSatish Balay     }
1326cc843e7aSLisandro Dalcin   }
1327cc843e7aSLisandro Dalcin #if !defined(PETSC_HAVE_POPEN)
132828b400f6SJacob Faibussowitsch   PetscCheck(!vbinary->storecompressed, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP_SYS, "Cannot run gzip on this machine");
1329cc843e7aSLisandro Dalcin #endif
1330cc843e7aSLisandro Dalcin 
1331cc843e7aSLisandro Dalcin   fname = vbinary->filename;
1332cc843e7aSLisandro Dalcin   if (vbinary->filemode == FILE_MODE_READ) { /* possibly get the file from remote site or compressed file */
13339566063dSJacob Faibussowitsch     PetscCall(PetscFileRetrieve(PetscObjectComm((PetscObject)viewer), fname, bname, PETSC_MAX_PATH_LEN, &found));
133428b400f6SJacob Faibussowitsch     PetscCheck(found, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_OPEN, "Cannot locate file: %s", fname);
1335cc843e7aSLisandro Dalcin     fname = bname;
13365c6c1daeSBarry Smith   }
13375c6c1daeSBarry Smith 
1338cc843e7aSLisandro Dalcin   vbinary->fdes = -1;
1339dd400576SPatrick Sanan   if (rank == 0) { /* only first processor opens file*/
1340cc843e7aSLisandro Dalcin     PetscFileMode mode = vbinary->filemode;
1341cc843e7aSLisandro Dalcin     if (mode == FILE_MODE_APPEND) {
1342cc843e7aSLisandro Dalcin       /* check if asked to append to a non-existing file */
13439566063dSJacob Faibussowitsch       PetscCall(PetscTestFile(fname, '\0', &found));
1344cc843e7aSLisandro Dalcin       if (!found) mode = FILE_MODE_WRITE;
1345cc843e7aSLisandro Dalcin     }
13469566063dSJacob Faibussowitsch     PetscCall(PetscBinaryOpen(fname, mode, &vbinary->fdes));
1347cc843e7aSLisandro Dalcin   }
1348cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1349cc843e7aSLisandro Dalcin }
1350cc843e7aSLisandro Dalcin 
13519371c9d4SSatish Balay static PetscErrorCode PetscViewerFileSetUp_BinaryInfo(PetscViewer viewer) {
1352cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1353cc843e7aSLisandro Dalcin   PetscMPIInt         rank;
1354cc843e7aSLisandro Dalcin   PetscBool           found;
1355cc843e7aSLisandro Dalcin 
1356cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1357cc843e7aSLisandro Dalcin   vbinary->fdes_info = NULL;
13589566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
1359dd400576SPatrick Sanan   if (!vbinary->skipinfo && (vbinary->filemode == FILE_MODE_READ || rank == 0)) {
1360cc843e7aSLisandro Dalcin     char infoname[PETSC_MAX_PATH_LEN], iname[PETSC_MAX_PATH_LEN], *gz;
1361cc843e7aSLisandro Dalcin 
13629566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(infoname, vbinary->filename, sizeof(infoname)));
1363cc843e7aSLisandro Dalcin     /* remove .gz if it ends file name */
13649566063dSJacob Faibussowitsch     PetscCall(PetscStrstr(infoname, ".gz", &gz));
1365cc843e7aSLisandro Dalcin     if (gz && gz[3] == 0) *gz = 0;
1366cc843e7aSLisandro Dalcin 
13679566063dSJacob Faibussowitsch     PetscCall(PetscStrlcat(infoname, ".info", sizeof(infoname)));
1368cc843e7aSLisandro Dalcin     if (vbinary->filemode == FILE_MODE_READ) {
13699566063dSJacob Faibussowitsch       PetscCall(PetscFixFilename(infoname, iname));
13709566063dSJacob Faibussowitsch       PetscCall(PetscFileRetrieve(PetscObjectComm((PetscObject)viewer), iname, infoname, PETSC_MAX_PATH_LEN, &found));
13719566063dSJacob Faibussowitsch       if (found) PetscCall(PetscOptionsInsertFile(PetscObjectComm((PetscObject)viewer), ((PetscObject)viewer)->options, infoname, PETSC_FALSE));
1372dd400576SPatrick Sanan     } else if (rank == 0) { /* write or append */
1373cc843e7aSLisandro Dalcin       const char *omode  = (vbinary->filemode == FILE_MODE_APPEND) ? "a" : "w";
1374cc843e7aSLisandro Dalcin       vbinary->fdes_info = fopen(infoname, omode);
137528b400f6SJacob Faibussowitsch       PetscCheck(vbinary->fdes_info, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open .info file %s for writing", infoname);
13765c6c1daeSBarry Smith     }
13775c6c1daeSBarry Smith   }
13785c6c1daeSBarry Smith   PetscFunctionReturn(0);
13795c6c1daeSBarry Smith }
13805c6c1daeSBarry Smith 
13819371c9d4SSatish Balay static PetscErrorCode PetscViewerSetUp_Binary(PetscViewer viewer) {
13825c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1383cc843e7aSLisandro Dalcin   PetscBool           usempiio;
1384cc843e7aSLisandro Dalcin 
13855c6c1daeSBarry Smith   PetscFunctionBegin;
13869566063dSJacob Faibussowitsch   if (!vbinary->setfromoptionscalled) PetscCall(PetscViewerSetFromOptions(viewer));
138728b400f6SJacob Faibussowitsch   PetscCheck(vbinary->filename, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetName()");
138808401ef6SPierre Jolivet   PetscCheck(vbinary->filemode != (PetscFileMode)-1, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode()");
13899566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_Binary(viewer));
1390cc843e7aSLisandro Dalcin 
13919566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetUseMPIIO(viewer, &usempiio));
1392cc843e7aSLisandro Dalcin   if (usempiio) {
1393cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
13949566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetUp_BinaryMPIIO(viewer));
1395cc843e7aSLisandro Dalcin #endif
1396cc843e7aSLisandro Dalcin   } else {
13979566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetUp_BinarySTDIO(viewer));
1398cc843e7aSLisandro Dalcin   }
13999566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetUp_BinaryInfo(viewer));
1400cc843e7aSLisandro Dalcin 
14019566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectState((PetscObject)viewer, "File: %s", vbinary->filename));
14025c6c1daeSBarry Smith   PetscFunctionReturn(0);
14035c6c1daeSBarry Smith }
14045c6c1daeSBarry Smith 
14059371c9d4SSatish Balay static PetscErrorCode PetscViewerView_Binary(PetscViewer v, PetscViewer viewer) {
1406cb6ad94fSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
1407cb6ad94fSLisandro Dalcin   const char         *fname   = vbinary->filename ? vbinary->filename : "not yet set";
1408cc843e7aSLisandro Dalcin   const char         *fmode   = vbinary->filemode != (PetscFileMode)-1 ? PetscFileModes[vbinary->filemode] : "not yet set";
1409cb6ad94fSLisandro Dalcin   PetscBool           usempiio;
14102bf49c77SBarry Smith 
14112bf49c77SBarry Smith   PetscFunctionBegin;
14129566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetUseMPIIO(v, &usempiio));
14139566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Filename: %s\n", fname));
14149566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Mode: %s (%s)\n", fmode, usempiio ? "mpiio" : "stdio"));
14152bf49c77SBarry Smith   PetscFunctionReturn(0);
14162bf49c77SBarry Smith }
14172bf49c77SBarry Smith 
14189371c9d4SSatish Balay static PetscErrorCode PetscViewerSetFromOptions_Binary(PetscViewer viewer, PetscOptionItems *PetscOptionsObject) {
141922a8f86cSLisandro Dalcin   PetscViewer_Binary *binary = (PetscViewer_Binary *)viewer->data;
1420d22fd6bcSDave May   char                defaultname[PETSC_MAX_PATH_LEN];
142103a1d59fSDave May   PetscBool           flg;
1422e0385b85SDave May 
142303a1d59fSDave May   PetscFunctionBegin;
142422a8f86cSLisandro Dalcin   if (viewer->setupcalled) PetscFunctionReturn(0);
1425d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Binary PetscViewer Options");
14269566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(defaultname, PETSC_MAX_PATH_LEN - 1, "binaryoutput"));
14279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-viewer_binary_filename", "Specify filename", "PetscViewerFileSetName", defaultname, defaultname, sizeof(defaultname), &flg));
14289566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscViewerFileSetName_Binary(viewer, defaultname));
14299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_skip_info", "Skip writing/reading .info file", "PetscViewerBinarySetSkipInfo", binary->skipinfo, &binary->skipinfo, NULL));
14309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_skip_options", "Skip parsing Vec/Mat load options", "PetscViewerBinarySetSkipOptions", binary->skipoptions, &binary->skipoptions, NULL));
14319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_skip_header", "Skip writing/reading header information", "PetscViewerBinarySetSkipHeader", binary->skipheader, &binary->skipheader, NULL));
143203a1d59fSDave May #if defined(PETSC_HAVE_MPIIO)
14339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_mpiio", "Use MPI-IO functionality to write/read binary file", "PetscViewerBinarySetUseMPIIO", binary->usempiio, &binary->usempiio, NULL));
14345972f5f3SLisandro Dalcin #else
14359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_binary_mpiio", "Use MPI-IO functionality to write/read binary file (NOT AVAILABLE)", "PetscViewerBinarySetUseMPIIO", PETSC_FALSE, NULL, NULL));
143603a1d59fSDave May #endif
1437d0609cedSBarry Smith   PetscOptionsHeadEnd();
1438bc196f7cSDave May   binary->setfromoptionscalled = PETSC_TRUE;
143903a1d59fSDave May   PetscFunctionReturn(0);
144003a1d59fSDave May }
144103a1d59fSDave May 
14428556b5ebSBarry Smith /*MC
14438556b5ebSBarry Smith    PETSCVIEWERBINARY - A viewer that saves to binary files
14448556b5ebSBarry Smith 
1445c2e3fba1SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PETSC_VIEWER_STDOUT_()`, `PETSC_VIEWER_STDOUT_SELF`, `PETSC_VIEWER_STDOUT_WORLD`, `PetscViewerCreate()`, `PetscViewerASCIIOpen()`,
1446db781477SPatrick Sanan           `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`, `PETSCVIEWERDRAW`,
1447db781477SPatrick Sanan           `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()`,
1448db781477SPatrick Sanan           `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`
14498556b5ebSBarry Smith 
14501b266c99SBarry Smith   Level: beginner
14511b266c99SBarry Smith 
14528556b5ebSBarry Smith M*/
14538556b5ebSBarry Smith 
14549371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscViewerCreate_Binary(PetscViewer v) {
14555c6c1daeSBarry Smith   PetscViewer_Binary *vbinary;
14565c6c1daeSBarry Smith 
14575c6c1daeSBarry Smith   PetscFunctionBegin;
14589566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(v, &vbinary));
14595c6c1daeSBarry Smith   v->data = (void *)vbinary;
1460cc843e7aSLisandro Dalcin 
146103a1d59fSDave May   v->ops->setfromoptions   = PetscViewerSetFromOptions_Binary;
14625c6c1daeSBarry Smith   v->ops->destroy          = PetscViewerDestroy_Binary;
14632bf49c77SBarry Smith   v->ops->view             = PetscViewerView_Binary;
1464c98fd787SBarry Smith   v->ops->setup            = PetscViewerSetUp_Binary;
1465cc843e7aSLisandro Dalcin   v->ops->flush            = NULL; /* Should we support Flush() ? */
1466cc843e7aSLisandro Dalcin   v->ops->getsubviewer     = PetscViewerGetSubViewer_Binary;
1467cc843e7aSLisandro Dalcin   v->ops->restoresubviewer = PetscViewerRestoreSubViewer_Binary;
1468cc843e7aSLisandro Dalcin   v->ops->read             = PetscViewerBinaryRead;
1469cc843e7aSLisandro Dalcin 
1470cc843e7aSLisandro Dalcin   vbinary->fdes = -1;
1471e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
1472cc843e7aSLisandro Dalcin   vbinary->usempiio = PETSC_FALSE;
1473e8a65b7dSLisandro Dalcin   vbinary->mfdes    = MPI_FILE_NULL;
1474e8a65b7dSLisandro Dalcin   vbinary->mfsub    = MPI_FILE_NULL;
1475e8a65b7dSLisandro Dalcin #endif
1476cc843e7aSLisandro Dalcin   vbinary->filename        = NULL;
14777e4fd573SVaclav Hapla   vbinary->filemode        = FILE_MODE_UNDEFINED;
147802c9f0b5SLisandro Dalcin   vbinary->fdes_info       = NULL;
14795c6c1daeSBarry Smith   vbinary->skipinfo        = PETSC_FALSE;
14805c6c1daeSBarry Smith   vbinary->skipoptions     = PETSC_TRUE;
14815c6c1daeSBarry Smith   vbinary->skipheader      = PETSC_FALSE;
14825c6c1daeSBarry Smith   vbinary->storecompressed = PETSC_FALSE;
1483f90597f1SStefano Zampini   vbinary->ogzfilename     = NULL;
14845c6c1daeSBarry Smith   vbinary->flowcontrol     = 256; /* seems a good number for Cray XT-5 */
14855c6c1daeSBarry Smith 
1486cc843e7aSLisandro Dalcin   vbinary->setfromoptionscalled = PETSC_FALSE;
1487cc843e7aSLisandro Dalcin 
14889566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetFlowControl_C", PetscViewerBinaryGetFlowControl_Binary));
14899566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetFlowControl_C", PetscViewerBinarySetFlowControl_Binary));
14909566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipHeader_C", PetscViewerBinaryGetSkipHeader_Binary));
14919566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipHeader_C", PetscViewerBinarySetSkipHeader_Binary));
14929566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipOptions_C", PetscViewerBinaryGetSkipOptions_Binary));
14939566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipOptions_C", PetscViewerBinarySetSkipOptions_Binary));
14949566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipInfo_C", PetscViewerBinaryGetSkipInfo_Binary));
14959566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipInfo_C", PetscViewerBinarySetSkipInfo_Binary));
14969566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetInfoPointer_C", PetscViewerBinaryGetInfoPointer_Binary));
14979566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_Binary));
14989566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_Binary));
14999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_Binary));
15009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_Binary));
15015c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
15029566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetUseMPIIO_C", PetscViewerBinaryGetUseMPIIO_Binary));
15039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetUseMPIIO_C", PetscViewerBinarySetUseMPIIO_Binary));
15045c6c1daeSBarry Smith #endif
15055c6c1daeSBarry Smith   PetscFunctionReturn(0);
15065c6c1daeSBarry Smith }
15075c6c1daeSBarry Smith 
15085c6c1daeSBarry Smith /* ---------------------------------------------------------------------*/
15095c6c1daeSBarry Smith /*
15105c6c1daeSBarry Smith     The variable Petsc_Viewer_Binary_keyval is used to indicate an MPI attribute that
15115c6c1daeSBarry Smith   is attached to a communicator, in this case the attribute is a PetscViewer.
15125c6c1daeSBarry Smith */
1513d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_Binary_keyval = MPI_KEYVAL_INVALID;
15145c6c1daeSBarry Smith 
15155c6c1daeSBarry Smith /*@C
15165c6c1daeSBarry Smith      PETSC_VIEWER_BINARY_ - Creates a binary PetscViewer shared by all processors
15175c6c1daeSBarry Smith                      in a communicator.
15185c6c1daeSBarry Smith 
1519d083f849SBarry Smith      Collective
15205c6c1daeSBarry Smith 
15215c6c1daeSBarry Smith      Input Parameter:
15225c6c1daeSBarry Smith .    comm - the MPI communicator to share the binary PetscViewer
15235c6c1daeSBarry Smith 
15245c6c1daeSBarry Smith      Level: intermediate
15255c6c1daeSBarry Smith 
15265c6c1daeSBarry Smith    Options Database Keys:
152710699b91SBarry Smith +    -viewer_binary_filename <name> - filename in which to store the binary data, defaults to binaryoutput
152810699b91SBarry Smith .    -viewer_binary_skip_info - true means do not create .info file for this viewer
152910699b91SBarry Smith .    -viewer_binary_skip_options - true means do not use the options database for this viewer
153010699b91SBarry Smith .    -viewer_binary_skip_header - true means do not store the usual header information in the binary file
153110699b91SBarry Smith -    -viewer_binary_mpiio - true means use the file via MPI-IO, maybe faster for large files and many MPI ranks
15325c6c1daeSBarry Smith 
15335c6c1daeSBarry Smith    Environmental variables:
153410699b91SBarry Smith -   PETSC_VIEWER_BINARY_FILENAME - filename in which to store the binary data, defaults to binaryoutput
15355c6c1daeSBarry Smith 
15365c6c1daeSBarry Smith      Notes:
15375c6c1daeSBarry Smith      Unlike almost all other PETSc routines, PETSC_VIEWER_BINARY_ does not return
15385c6c1daeSBarry Smith      an error code.  The binary PetscViewer is usually used in the form
15395c6c1daeSBarry Smith $       XXXView(XXX object,PETSC_VIEWER_BINARY_(comm));
15405c6c1daeSBarry Smith 
1541db781477SPatrick Sanan .seealso: `PETSC_VIEWER_BINARY_WORLD`, `PETSC_VIEWER_BINARY_SELF`, `PetscViewerBinaryOpen()`, `PetscViewerCreate()`,
1542db781477SPatrick Sanan           `PetscViewerDestroy()`
15435c6c1daeSBarry Smith @*/
15449371c9d4SSatish Balay PetscViewer PETSC_VIEWER_BINARY_(MPI_Comm comm) {
15455c6c1daeSBarry Smith   PetscErrorCode ierr;
15465c6c1daeSBarry Smith   PetscBool      flg;
15475c6c1daeSBarry Smith   PetscViewer    viewer;
15485c6c1daeSBarry Smith   char           fname[PETSC_MAX_PATH_LEN];
15495c6c1daeSBarry Smith   MPI_Comm       ncomm;
15505c6c1daeSBarry Smith 
15515c6c1daeSBarry Smith   PetscFunctionBegin;
15529371c9d4SSatish Balay   ierr = PetscCommDuplicate(comm, &ncomm, NULL);
15539371c9d4SSatish Balay   if (ierr) {
15549371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
15559371c9d4SSatish Balay     PetscFunctionReturn(NULL);
15569371c9d4SSatish Balay   }
15575c6c1daeSBarry Smith   if (Petsc_Viewer_Binary_keyval == MPI_KEYVAL_INVALID) {
155802c9f0b5SLisandro Dalcin     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Viewer_Binary_keyval, NULL);
15599371c9d4SSatish Balay     if (ierr) {
15609371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
15619371c9d4SSatish Balay       PetscFunctionReturn(NULL);
15629371c9d4SSatish Balay     }
15635c6c1daeSBarry Smith   }
156447435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm, Petsc_Viewer_Binary_keyval, (void **)&viewer, (int *)&flg);
15659371c9d4SSatish Balay   if (ierr) {
15669371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
15679371c9d4SSatish Balay     PetscFunctionReturn(NULL);
15689371c9d4SSatish Balay   }
15695c6c1daeSBarry Smith   if (!flg) { /* PetscViewer not yet created */
15705c6c1daeSBarry Smith     ierr = PetscOptionsGetenv(ncomm, "PETSC_VIEWER_BINARY_FILENAME", fname, PETSC_MAX_PATH_LEN, &flg);
15719371c9d4SSatish Balay     if (ierr) {
15729371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
15739371c9d4SSatish Balay       PetscFunctionReturn(NULL);
15749371c9d4SSatish Balay     }
15755c6c1daeSBarry Smith     if (!flg) {
15765c6c1daeSBarry Smith       ierr = PetscStrcpy(fname, "binaryoutput");
15779371c9d4SSatish Balay       if (ierr) {
15789371c9d4SSatish Balay         PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
15799371c9d4SSatish Balay         PetscFunctionReturn(NULL);
15809371c9d4SSatish Balay       }
15815c6c1daeSBarry Smith     }
15825c6c1daeSBarry Smith     ierr = PetscViewerBinaryOpen(ncomm, fname, FILE_MODE_WRITE, &viewer);
15839371c9d4SSatish Balay     if (ierr) {
15849371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
15859371c9d4SSatish Balay       PetscFunctionReturn(NULL);
15869371c9d4SSatish Balay     }
15875c6c1daeSBarry Smith     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
15889371c9d4SSatish Balay     if (ierr) {
15899371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
15909371c9d4SSatish Balay       PetscFunctionReturn(NULL);
15919371c9d4SSatish Balay     }
159247435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm, Petsc_Viewer_Binary_keyval, (void *)viewer);
15939371c9d4SSatish Balay     if (ierr) {
15949371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
15959371c9d4SSatish Balay       PetscFunctionReturn(NULL);
15969371c9d4SSatish Balay     }
15975c6c1daeSBarry Smith   }
15985c6c1daeSBarry Smith   ierr = PetscCommDestroy(&ncomm);
15999371c9d4SSatish Balay   if (ierr) {
16009371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_BINARY_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16019371c9d4SSatish Balay     PetscFunctionReturn(NULL);
16029371c9d4SSatish Balay   }
16035c6c1daeSBarry Smith   PetscFunctionReturn(viewer);
16045c6c1daeSBarry Smith }
1605