1af0996ceSBarry Smith #include <petsc/private/viewerimpl.h> /*I "petscviewer.h" I*/
25c6c1daeSBarry Smith
376667918SBarry Smith /*
476667918SBarry Smith This needs to start the same as PetscViewer_Socket.
576667918SBarry Smith */
65c6c1daeSBarry Smith typedef struct {
75c6c1daeSBarry Smith int fdes; /* file descriptor, ignored if using MPI IO */
876667918SBarry Smith PetscInt flowcontrol; /* allow only <flowcontrol> messages outstanding at a time while doing IO */
976667918SBarry Smith PetscBool skipheader; /* don't write header, only raw data */
105c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
11bc196f7cSDave May PetscBool usempiio;
125c6c1daeSBarry Smith MPI_File mfdes; /* ignored unless using MPI IO */
13e8a65b7dSLisandro Dalcin MPI_File mfsub; /* subviewer support */
145c6c1daeSBarry Smith MPI_Offset moff;
155c6c1daeSBarry Smith #endif
16cc843e7aSLisandro Dalcin char *filename; /* file name */
17cc843e7aSLisandro Dalcin PetscFileMode filemode; /* read/write/append mode */
185c6c1daeSBarry Smith FILE *fdes_info; /* optional file containing info on binary file*/
195c6c1daeSBarry Smith PetscBool storecompressed; /* gzip the write binary file when closing it*/
20f90597f1SStefano Zampini char *ogzfilename; /* gzip can be run after the filename has been updated */
215c6c1daeSBarry Smith PetscBool skipinfo; /* Don't create info file for writing; don't use for reading */
225c6c1daeSBarry Smith PetscBool skipoptions; /* don't use PETSc options database when loading */
23a261c58fSBarry Smith PetscBool matlabheaderwritten; /* if format is PETSC_VIEWER_BINARY_MATLAB has the MATLAB .info header been written yet */
24c98fd787SBarry Smith PetscBool setfromoptionscalled;
255c6c1daeSBarry Smith } PetscViewer_Binary;
265c6c1daeSBarry Smith
PetscViewerBinaryClearFunctionList(PetscViewer v)27d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryClearFunctionList(PetscViewer v)
28d71ae5a4SJacob Faibussowitsch {
292e956fe4SStefano Zampini PetscFunctionBegin;
302e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetFlowControl_C", NULL));
312e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetFlowControl_C", NULL));
322e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipHeader_C", NULL));
332e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipHeader_C", NULL));
342e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipOptions_C", NULL));
352e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipOptions_C", NULL));
362e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipInfo_C", NULL));
372e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipInfo_C", NULL));
382e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetInfoPointer_C", NULL));
392e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", NULL));
402e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", NULL));
412e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", NULL));
422e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", NULL));
432e956fe4SStefano Zampini #if defined(PETSC_HAVE_MPIIO)
442e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetUseMPIIO_C", NULL));
452e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetUseMPIIO_C", NULL));
462e956fe4SStefano Zampini #endif
473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
482e956fe4SStefano Zampini }
492e956fe4SStefano Zampini
50cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
PetscViewerBinarySyncMPIIO(PetscViewer viewer)51d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySyncMPIIO(PetscViewer viewer)
52d71ae5a4SJacob Faibussowitsch {
53cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
54cc843e7aSLisandro Dalcin
55cc843e7aSLisandro Dalcin PetscFunctionBegin;
563ba16761SJacob Faibussowitsch if (vbinary->filemode == FILE_MODE_READ) PetscFunctionReturn(PETSC_SUCCESS);
5748a46eb9SPierre Jolivet if (vbinary->mfsub != MPI_FILE_NULL) PetscCallMPI(MPI_File_sync(vbinary->mfsub));
58cc843e7aSLisandro Dalcin if (vbinary->mfdes != MPI_FILE_NULL) {
599566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer)));
609566063dSJacob Faibussowitsch PetscCallMPI(MPI_File_sync(vbinary->mfdes));
61cc843e7aSLisandro Dalcin }
623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
63cc843e7aSLisandro Dalcin }
64cc843e7aSLisandro Dalcin #endif
65cc843e7aSLisandro Dalcin
PetscViewerGetSubViewer_Binary(PetscViewer viewer,MPI_Comm comm,PetscViewer * outviewer)66d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerGetSubViewer_Binary(PetscViewer viewer, MPI_Comm comm, PetscViewer *outviewer)
67d71ae5a4SJacob Faibussowitsch {
68e8a65b7dSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
69cc843e7aSLisandro Dalcin PetscMPIInt rank;
705c6c1daeSBarry Smith
715c6c1daeSBarry Smith PetscFunctionBegin;
729566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer));
73e8a65b7dSLisandro Dalcin
74e8a65b7dSLisandro Dalcin /* Return subviewer in process zero */
759566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
76dd400576SPatrick Sanan if (rank == 0) {
77e5afcf28SBarry Smith PetscMPIInt flg;
78e5afcf28SBarry Smith
799566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(PETSC_COMM_SELF, comm, &flg));
80cc73adaaSBarry Smith PetscCheck(flg == MPI_IDENT || flg == MPI_CONGRUENT, PETSC_COMM_SELF, PETSC_ERR_SUP, "PetscViewerGetSubViewer() for PETSCVIEWERBINARY requires a singleton MPI_Comm");
819566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(comm, outviewer));
829566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(*outviewer, PETSCVIEWERBINARY));
839566063dSJacob Faibussowitsch PetscCall(PetscMemcpy((*outviewer)->data, vbinary, sizeof(PetscViewer_Binary)));
8412f4c3a9SLisandro Dalcin (*outviewer)->setupcalled = PETSC_TRUE;
8596509d17SLisandro Dalcin } else {
8696509d17SLisandro Dalcin *outviewer = NULL;
8796509d17SLisandro Dalcin }
88e8a65b7dSLisandro Dalcin
89e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
90e8a65b7dSLisandro Dalcin if (vbinary->usempiio && *outviewer) {
91e8a65b7dSLisandro Dalcin PetscViewer_Binary *obinary = (PetscViewer_Binary *)(*outviewer)->data;
92e8a65b7dSLisandro Dalcin /* Parent viewer opens a new MPI file handle on PETSC_COMM_SELF and keeps track of it for future reuse */
93e8a65b7dSLisandro Dalcin if (vbinary->mfsub == MPI_FILE_NULL) {
94e8a65b7dSLisandro Dalcin int amode;
95cc843e7aSLisandro Dalcin switch (vbinary->filemode) {
96d71ae5a4SJacob Faibussowitsch case FILE_MODE_READ:
97d71ae5a4SJacob Faibussowitsch amode = MPI_MODE_RDONLY;
98d71ae5a4SJacob Faibussowitsch break;
99d71ae5a4SJacob Faibussowitsch case FILE_MODE_WRITE:
100d71ae5a4SJacob Faibussowitsch amode = MPI_MODE_WRONLY;
101d71ae5a4SJacob Faibussowitsch break;
102d71ae5a4SJacob Faibussowitsch case FILE_MODE_APPEND:
103d71ae5a4SJacob Faibussowitsch amode = MPI_MODE_WRONLY;
104d71ae5a4SJacob Faibussowitsch break;
105d71ae5a4SJacob Faibussowitsch default:
106d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[vbinary->filemode]);
107e8a65b7dSLisandro Dalcin }
1089566063dSJacob Faibussowitsch PetscCallMPI(MPI_File_open(PETSC_COMM_SELF, vbinary->filename, amode, MPI_INFO_NULL, &vbinary->mfsub));
109e8a65b7dSLisandro Dalcin }
110e8a65b7dSLisandro Dalcin /* Subviewer gets the MPI file handle on PETSC_COMM_SELF */
111e8a65b7dSLisandro Dalcin obinary->mfdes = vbinary->mfsub;
112e8a65b7dSLisandro Dalcin obinary->mfsub = MPI_FILE_NULL;
113e8a65b7dSLisandro Dalcin obinary->moff = vbinary->moff;
114e8a65b7dSLisandro Dalcin }
115e8a65b7dSLisandro Dalcin #endif
116cc843e7aSLisandro Dalcin
117cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
1189566063dSJacob Faibussowitsch PetscCall(PetscViewerBinarySyncMPIIO(viewer));
119cc843e7aSLisandro Dalcin #endif
1203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1215c6c1daeSBarry Smith }
1225c6c1daeSBarry Smith
PetscViewerRestoreSubViewer_Binary(PetscViewer viewer,MPI_Comm comm,PetscViewer * outviewer)123d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerRestoreSubViewer_Binary(PetscViewer viewer, MPI_Comm comm, PetscViewer *outviewer)
124d71ae5a4SJacob Faibussowitsch {
125e8a65b7dSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
126cc843e7aSLisandro Dalcin PetscMPIInt rank;
127e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
128e8a65b7dSLisandro Dalcin MPI_Offset moff = 0;
129e8a65b7dSLisandro Dalcin #endif
1305c6c1daeSBarry Smith
1315c6c1daeSBarry Smith PetscFunctionBegin;
1329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
133c5853193SPierre Jolivet PetscCheck(rank == 0 || !*outviewer, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Subviewer not obtained from viewer");
134e8a65b7dSLisandro Dalcin
135e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
136e8a65b7dSLisandro Dalcin if (vbinary->usempiio && *outviewer) {
137e8a65b7dSLisandro Dalcin PetscViewer_Binary *obinary = (PetscViewer_Binary *)(*outviewer)->data;
13808401ef6SPierre Jolivet PetscCheck(obinary->mfdes == vbinary->mfsub, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Subviewer not obtained from viewer");
1399566063dSJacob Faibussowitsch if (obinary->mfsub != MPI_FILE_NULL) PetscCallMPI(MPI_File_close(&obinary->mfsub));
140e8a65b7dSLisandro Dalcin moff = obinary->moff;
141e8a65b7dSLisandro Dalcin }
142e8a65b7dSLisandro Dalcin #endif
143e8a65b7dSLisandro Dalcin
144e8a65b7dSLisandro Dalcin if (*outviewer) {
145e8a65b7dSLisandro Dalcin PetscViewer_Binary *obinary = (PetscViewer_Binary *)(*outviewer)->data;
14608401ef6SPierre Jolivet PetscCheck(obinary->fdes == vbinary->fdes, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Subviewer not obtained from viewer");
1479566063dSJacob Faibussowitsch PetscCall(PetscFree((*outviewer)->data));
1482e956fe4SStefano Zampini PetscCall(PetscViewerBinaryClearFunctionList(*outviewer));
1499566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(outviewer));
1505c6c1daeSBarry Smith }
151e8a65b7dSLisandro Dalcin
152e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
153e8a65b7dSLisandro Dalcin if (vbinary->usempiio) {
154e8a65b7dSLisandro Dalcin PetscInt64 ioff = (PetscInt64)moff; /* We could use MPI_OFFSET datatype (requires MPI 2.2) */
1559566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&ioff, 1, MPIU_INT64, 0, PetscObjectComm((PetscObject)viewer)));
156e8a65b7dSLisandro Dalcin vbinary->moff = (MPI_Offset)ioff;
157e8a65b7dSLisandro Dalcin }
158e8a65b7dSLisandro Dalcin #endif
159cc843e7aSLisandro Dalcin
160cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
1619566063dSJacob Faibussowitsch PetscCall(PetscViewerBinarySyncMPIIO(viewer));
162cc843e7aSLisandro Dalcin #endif
1633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1645c6c1daeSBarry Smith }
1655c6c1daeSBarry Smith
1665c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
1675c6c1daeSBarry Smith /*@C
16885b8072bSPatrick Sanan PetscViewerBinaryGetMPIIOOffset - Gets the current global offset that should be passed to `MPI_File_set_view()` or `MPI_File_{write|read}_at[_all]()`
1695c6c1daeSBarry Smith
170cf53795eSBarry Smith Not Collective; No Fortran Support
1715c6c1daeSBarry Smith
1725c6c1daeSBarry Smith Input Parameter:
173c410d8ccSBarry Smith . viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
1745c6c1daeSBarry Smith
1755c6c1daeSBarry Smith Output Parameter:
17622a8f86cSLisandro Dalcin . off - the current global offset
1775c6c1daeSBarry Smith
1785c6c1daeSBarry Smith Level: advanced
1795c6c1daeSBarry Smith
180811af0c4SBarry Smith Note:
181811af0c4SBarry Smith Use `PetscViewerBinaryAddMPIIOOffset()` to increase this value after you have written a view.
182811af0c4SBarry Smith
183d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryAddMPIIOOffset()`
1845c6c1daeSBarry Smith @*/
PetscViewerBinaryGetMPIIOOffset(PetscViewer viewer,MPI_Offset * off)185d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetMPIIOOffset(PetscViewer viewer, MPI_Offset *off)
186d71ae5a4SJacob Faibussowitsch {
18722a8f86cSLisandro Dalcin PetscViewer_Binary *vbinary;
1885c6c1daeSBarry Smith
1895c6c1daeSBarry Smith PetscFunctionBegin;
19022a8f86cSLisandro Dalcin PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
1914f572ea9SToby Isaac PetscAssertPointer(off, 2);
19222a8f86cSLisandro Dalcin vbinary = (PetscViewer_Binary *)viewer->data;
1935c6c1daeSBarry Smith *off = vbinary->moff;
1943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1955c6c1daeSBarry Smith }
1965c6c1daeSBarry Smith
1975c6c1daeSBarry Smith /*@C
19822a8f86cSLisandro Dalcin PetscViewerBinaryAddMPIIOOffset - Adds to the current global offset
1995c6c1daeSBarry Smith
200cf53795eSBarry Smith Logically Collective; No Fortran Support
2015c6c1daeSBarry Smith
2025c6c1daeSBarry Smith Input Parameters:
203811af0c4SBarry Smith + viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
20422a8f86cSLisandro Dalcin - off - the addition to the global offset
2055c6c1daeSBarry Smith
2065c6c1daeSBarry Smith Level: advanced
2075c6c1daeSBarry Smith
208811af0c4SBarry Smith Note:
209811af0c4SBarry Smith Use `PetscViewerBinaryGetMPIIOOffset()` to get the value that you should pass to `MPI_File_set_view()` or `MPI_File_{write|read}_at[_all]()`
210811af0c4SBarry Smith
211d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
2125c6c1daeSBarry Smith @*/
PetscViewerBinaryAddMPIIOOffset(PetscViewer viewer,MPI_Offset off)213d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryAddMPIIOOffset(PetscViewer viewer, MPI_Offset off)
214d71ae5a4SJacob Faibussowitsch {
21522a8f86cSLisandro Dalcin PetscViewer_Binary *vbinary;
2165c6c1daeSBarry Smith
2175c6c1daeSBarry Smith PetscFunctionBegin;
21822a8f86cSLisandro Dalcin PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
21922a8f86cSLisandro Dalcin PetscValidLogicalCollectiveInt(viewer, (PetscInt)off, 2);
22022a8f86cSLisandro Dalcin vbinary = (PetscViewer_Binary *)viewer->data;
2215c6c1daeSBarry Smith vbinary->moff += off;
2223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2235c6c1daeSBarry Smith }
2245c6c1daeSBarry Smith
2255c6c1daeSBarry Smith /*@C
226811af0c4SBarry Smith PetscViewerBinaryGetMPIIODescriptor - Extracts the MPI IO file descriptor from a `PetscViewer`.
2275c6c1daeSBarry Smith
228cf53795eSBarry Smith Not Collective; No Fortran Support
2295c6c1daeSBarry Smith
2305c6c1daeSBarry Smith Input Parameter:
231811af0c4SBarry Smith . viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
2325c6c1daeSBarry Smith
2335c6c1daeSBarry Smith Output Parameter:
2345c6c1daeSBarry Smith . fdes - file descriptor
2355c6c1daeSBarry Smith
2365c6c1daeSBarry Smith Level: advanced
2375c6c1daeSBarry Smith
238d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
2395c6c1daeSBarry Smith @*/
PetscViewerBinaryGetMPIIODescriptor(PetscViewer viewer,MPI_File * fdes)240d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetMPIIODescriptor(PetscViewer viewer, MPI_File *fdes)
241d71ae5a4SJacob Faibussowitsch {
24222a8f86cSLisandro Dalcin PetscViewer_Binary *vbinary;
2435c6c1daeSBarry Smith
2445c6c1daeSBarry Smith PetscFunctionBegin;
24522a8f86cSLisandro Dalcin PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
2464f572ea9SToby Isaac PetscAssertPointer(fdes, 2);
2479566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer));
24822a8f86cSLisandro Dalcin vbinary = (PetscViewer_Binary *)viewer->data;
2495c6c1daeSBarry Smith *fdes = vbinary->mfdes;
2503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2515c6c1daeSBarry Smith }
252cc843e7aSLisandro Dalcin #endif
2535c6c1daeSBarry Smith
254cc843e7aSLisandro Dalcin /*@
255cc843e7aSLisandro Dalcin PetscViewerBinarySetUseMPIIO - Sets a binary viewer to use MPI-IO for reading/writing. Must be called
256811af0c4SBarry Smith before `PetscViewerFileSetName()`
257cc843e7aSLisandro Dalcin
258c3339decSBarry Smith Logically Collective
259cc843e7aSLisandro Dalcin
260cc843e7aSLisandro Dalcin Input Parameters:
261811af0c4SBarry Smith + viewer - the `PetscViewer`; must be a `PETSCVIEWERBINARY`
262811af0c4SBarry Smith - use - `PETSC_TRUE` means MPI-IO will be used
263cc843e7aSLisandro Dalcin
264811af0c4SBarry Smith Options Database Key:
26538b83642SBarry Smith . -viewer_binary_mpiio - <true or false> flag for using MPI-IO
266cc843e7aSLisandro Dalcin
267cc843e7aSLisandro Dalcin Level: advanced
268cc843e7aSLisandro Dalcin
269d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`,
270db781477SPatrick Sanan `PetscViewerBinaryGetUseMPIIO()`
271cc843e7aSLisandro Dalcin @*/
PetscViewerBinarySetUseMPIIO(PetscViewer viewer,PetscBool use)272d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetUseMPIIO(PetscViewer viewer, PetscBool use)
273d71ae5a4SJacob Faibussowitsch {
274a8aae444SDave May PetscFunctionBegin;
275cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
276cc843e7aSLisandro Dalcin PetscValidLogicalCollectiveBool(viewer, use, 2);
277cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerBinarySetUseMPIIO_C", (PetscViewer, PetscBool), (viewer, use));
2783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
279cc843e7aSLisandro Dalcin }
280cc843e7aSLisandro Dalcin
281cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
PetscViewerBinarySetUseMPIIO_Binary(PetscViewer viewer,PetscBool use)282d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetUseMPIIO_Binary(PetscViewer viewer, PetscBool use)
283d71ae5a4SJacob Faibussowitsch {
284cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
2854d86920dSPierre Jolivet
286cc843e7aSLisandro Dalcin PetscFunctionBegin;
287cc73adaaSBarry Smith PetscCheck(!viewer->setupcalled || vbinary->usempiio == use, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Cannot change MPIIO to %s after setup", PetscBools[use]);
288cc843e7aSLisandro Dalcin vbinary->usempiio = use;
2893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
290a8aae444SDave May }
291a8aae444SDave May #endif
292a8aae444SDave May
293cc843e7aSLisandro Dalcin /*@
2943f423023SBarry Smith PetscViewerBinaryGetUseMPIIO - Returns `PETSC_TRUE` if the binary viewer uses MPI-IO.
2955c6c1daeSBarry Smith
2965c6c1daeSBarry Smith Not Collective
2975c6c1daeSBarry Smith
2985c6c1daeSBarry Smith Input Parameter:
299811af0c4SBarry Smith . viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`; must be a `PETSCVIEWERBINARY`
3005c6c1daeSBarry Smith
3015c6c1daeSBarry Smith Output Parameter:
302811af0c4SBarry Smith . use - `PETSC_TRUE` if MPI-IO is being used
3035c6c1daeSBarry Smith
3045c6c1daeSBarry Smith Level: advanced
3055c6c1daeSBarry Smith
306bc196f7cSDave May Note:
3073f423023SBarry Smith If MPI-IO is not available, this function will always return `PETSC_FALSE`
308bc196f7cSDave May
309d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
3105c6c1daeSBarry Smith @*/
PetscViewerBinaryGetUseMPIIO(PetscViewer viewer,PetscBool * use)311d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetUseMPIIO(PetscViewer viewer, PetscBool *use)
312d71ae5a4SJacob Faibussowitsch {
3135c6c1daeSBarry Smith PetscFunctionBegin;
314cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
3154f572ea9SToby Isaac PetscAssertPointer(use, 2);
316cc843e7aSLisandro Dalcin *use = PETSC_FALSE;
317cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerBinaryGetUseMPIIO_C", (PetscViewer, PetscBool *), (viewer, use));
3183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3195c6c1daeSBarry Smith }
3205c6c1daeSBarry Smith
321cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
PetscViewerBinaryGetUseMPIIO_Binary(PetscViewer viewer,PetscBool * use)322d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetUseMPIIO_Binary(PetscViewer viewer, PetscBool *use)
323d71ae5a4SJacob Faibussowitsch {
3245c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
3255c6c1daeSBarry Smith
3265c6c1daeSBarry Smith PetscFunctionBegin;
327cc843e7aSLisandro Dalcin *use = vbinary->usempiio;
3283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
329cc843e7aSLisandro Dalcin }
330cc843e7aSLisandro Dalcin #endif
331cc843e7aSLisandro Dalcin
332cc843e7aSLisandro Dalcin /*@
333c410d8ccSBarry Smith PetscViewerBinarySetFlowControl - Sets how many messages are allowed to be outstanding at the same time during parallel IO reads/writes
334cc843e7aSLisandro Dalcin
335cc843e7aSLisandro Dalcin Not Collective
336cc843e7aSLisandro Dalcin
337d8d19677SJose E. Roman Input Parameters:
338c410d8ccSBarry Smith + viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
339cc843e7aSLisandro Dalcin - fc - the number of messages, defaults to 256 if this function was not called
340cc843e7aSLisandro Dalcin
341cc843e7aSLisandro Dalcin Level: advanced
342cc843e7aSLisandro Dalcin
343d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetFlowControl()`
344cc843e7aSLisandro Dalcin @*/
PetscViewerBinarySetFlowControl(PetscViewer viewer,PetscInt fc)345d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetFlowControl(PetscViewer viewer, PetscInt fc)
346d71ae5a4SJacob Faibussowitsch {
347cc843e7aSLisandro Dalcin PetscFunctionBegin;
348cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
349cc843e7aSLisandro Dalcin PetscValidLogicalCollectiveInt(viewer, fc, 2);
350cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerBinarySetFlowControl_C", (PetscViewer, PetscInt), (viewer, fc));
3513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3525c6c1daeSBarry Smith }
3535c6c1daeSBarry Smith
PetscViewerBinarySetFlowControl_Binary(PetscViewer viewer,PetscInt fc)354d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetFlowControl_Binary(PetscViewer viewer, PetscInt fc)
355d71ae5a4SJacob Faibussowitsch {
356cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
357cc843e7aSLisandro Dalcin
358cc843e7aSLisandro Dalcin PetscFunctionBegin;
35908401ef6SPierre Jolivet PetscCheck(fc > 1, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_OUTOFRANGE, "Flow control count must be greater than 1, %" PetscInt_FMT " was set", fc);
360cc843e7aSLisandro Dalcin vbinary->flowcontrol = fc;
3613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
362cc843e7aSLisandro Dalcin }
363cc843e7aSLisandro Dalcin
364cc843e7aSLisandro Dalcin /*@
365c410d8ccSBarry Smith PetscViewerBinaryGetFlowControl - Returns how many messages are allowed to be outstanding at the same time during parallel IO reads/writes
3665c6c1daeSBarry Smith
3675c6c1daeSBarry Smith Not Collective
3685c6c1daeSBarry Smith
3695c6c1daeSBarry Smith Input Parameter:
370811af0c4SBarry Smith . viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
3715c6c1daeSBarry Smith
3725c6c1daeSBarry Smith Output Parameter:
3735c6c1daeSBarry Smith . fc - the number of messages
3745c6c1daeSBarry Smith
3755c6c1daeSBarry Smith Level: advanced
3765c6c1daeSBarry Smith
377d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinarySetFlowControl()`
3785c6c1daeSBarry Smith @*/
PetscViewerBinaryGetFlowControl(PetscViewer viewer,PetscInt * fc)379d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetFlowControl(PetscViewer viewer, PetscInt *fc)
380d71ae5a4SJacob Faibussowitsch {
3815c6c1daeSBarry Smith PetscFunctionBegin;
382cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
3834f572ea9SToby Isaac PetscAssertPointer(fc, 2);
384cac4c232SBarry Smith PetscUseMethod(viewer, "PetscViewerBinaryGetFlowControl_C", (PetscViewer, PetscInt *), (viewer, fc));
3853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3865c6c1daeSBarry Smith }
3875c6c1daeSBarry Smith
PetscViewerBinaryGetFlowControl_Binary(PetscViewer viewer,PetscInt * fc)38876667918SBarry Smith PETSC_INTERN PetscErrorCode PetscViewerBinaryGetFlowControl_Binary(PetscViewer viewer, PetscInt *fc)
389d71ae5a4SJacob Faibussowitsch {
3905c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
3915c6c1daeSBarry Smith
3925c6c1daeSBarry Smith PetscFunctionBegin;
393cc843e7aSLisandro Dalcin *fc = vbinary->flowcontrol;
3943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3955c6c1daeSBarry Smith }
3965c6c1daeSBarry Smith
3975c6c1daeSBarry Smith /*@C
398811af0c4SBarry Smith PetscViewerBinaryGetDescriptor - Extracts the file descriptor from a `PetscViewer` of `PetscViewerType` `PETSCVIEWERBINARY`.
3995c6c1daeSBarry Smith
400cd05f99aSJacob Faibussowitsch Collective because it may trigger a `PetscViewerSetUp()` call; No Fortran Support
4015c6c1daeSBarry Smith
4025c6c1daeSBarry Smith Input Parameter:
403811af0c4SBarry Smith . viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
4045c6c1daeSBarry Smith
4055c6c1daeSBarry Smith Output Parameter:
4065c6c1daeSBarry Smith . fdes - file descriptor
4075c6c1daeSBarry Smith
4085c6c1daeSBarry Smith Level: advanced
4095c6c1daeSBarry Smith
410811af0c4SBarry Smith Note:
411811af0c4SBarry Smith For writable binary `PetscViewer`s, the descriptor will only be valid for the
412811af0c4SBarry Smith first processor in the communicator that shares the `PetscViewer`. For readable
413c410d8ccSBarry Smith files it will only be valid on processes that have the file. If MPI rank 0 does not
414c410d8ccSBarry Smith have the file it generates an error even if another MPI process does have the file.
4155c6c1daeSBarry Smith
416d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`
4175c6c1daeSBarry Smith @*/
PetscViewerBinaryGetDescriptor(PetscViewer viewer,int * fdes)418d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetDescriptor(PetscViewer viewer, int *fdes)
419d71ae5a4SJacob Faibussowitsch {
42022a8f86cSLisandro Dalcin PetscViewer_Binary *vbinary;
4215c6c1daeSBarry Smith
4225c6c1daeSBarry Smith PetscFunctionBegin;
42322a8f86cSLisandro Dalcin PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
4244f572ea9SToby Isaac PetscAssertPointer(fdes, 2);
4259566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer));
42622a8f86cSLisandro Dalcin vbinary = (PetscViewer_Binary *)viewer->data;
4275c6c1daeSBarry Smith *fdes = vbinary->fdes;
4283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4295c6c1daeSBarry Smith }
4305c6c1daeSBarry Smith
4315c6c1daeSBarry Smith /*@
432c410d8ccSBarry Smith PetscViewerBinarySkipInfo - Binary file will not have `.info` file created with it
4335c6c1daeSBarry Smith
4345c6c1daeSBarry Smith Not Collective
4355c6c1daeSBarry Smith
436fd292e60Sprj- Input Parameter:
437811af0c4SBarry Smith . viewer - `PetscViewer` context, obtained from `PetscViewerCreate()`
4385c6c1daeSBarry Smith
4395c6c1daeSBarry Smith Options Database Key:
440c410d8ccSBarry Smith . -viewer_binary_skip_info - true indicates do not generate `.info` file
4415c6c1daeSBarry Smith
4425c6c1daeSBarry Smith Level: advanced
4435c6c1daeSBarry Smith
44495452b02SPatrick Sanan Notes:
445811af0c4SBarry Smith This must be called after `PetscViewerSetType()`. If you use `PetscViewerBinaryOpen()` then
4463f423023SBarry Smith you can only skip the info file with the `-viewer_binary_skip_info` flag. To use the function you must open the
447811af0c4SBarry Smith viewer with `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinarySkipInfo()`.
4485c6c1daeSBarry Smith
449c410d8ccSBarry Smith The `.info` files contains meta information about the data in the binary file, for example the block size if it was
4505c6c1daeSBarry Smith set for a vector or matrix.
4515c6c1daeSBarry Smith
452811af0c4SBarry Smith This routine is deprecated, use `PetscViewerBinarySetSkipInfo()`
453811af0c4SBarry Smith
454d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySetSkipOptions()`,
455db781477SPatrick Sanan `PetscViewerBinaryGetSkipOptions()`, `PetscViewerBinaryGetSkipInfo()`
4565c6c1daeSBarry Smith @*/
PetscViewerBinarySkipInfo(PetscViewer viewer)457d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySkipInfo(PetscViewer viewer)
458d71ae5a4SJacob Faibussowitsch {
4595c6c1daeSBarry Smith PetscFunctionBegin;
4609566063dSJacob Faibussowitsch PetscCall(PetscViewerBinarySetSkipInfo(viewer, PETSC_TRUE));
4613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
462807ea322SDave May }
463807ea322SDave May
464807ea322SDave May /*@
465c410d8ccSBarry Smith PetscViewerBinarySetSkipInfo - Binary file will not have `.info` file created with it
466807ea322SDave May
467807ea322SDave May Not Collective
468807ea322SDave May
469d8d19677SJose E. Roman Input Parameters:
4703f423023SBarry Smith + viewer - PetscViewer context, obtained from `PetscViewerCreate()`
471c410d8ccSBarry Smith - skip - `PETSC_TRUE` implies the `.info` file will not be generated
472807ea322SDave May
473807ea322SDave May Options Database Key:
474c410d8ccSBarry Smith . -viewer_binary_skip_info - true indicates do not generate `.info` file
475807ea322SDave May
476807ea322SDave May Level: advanced
477807ea322SDave May
478c410d8ccSBarry Smith Notes:
479c410d8ccSBarry Smith This must be called after `PetscViewerSetType()`. If you use `PetscViewerBinaryOpen()` then
480c410d8ccSBarry Smith you can only skip the info file with the `-viewer_binary_skip_info` flag. To use the function you must open the
481c410d8ccSBarry Smith viewer with `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinarySkipInfo()`.
482c410d8ccSBarry Smith
483c410d8ccSBarry Smith The `.info` file contains meta information about the data in the binary file, for example the block size if it was
484c410d8ccSBarry Smith set for a vector or matrix.
485c410d8ccSBarry Smith
486d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySetSkipOptions()`,
487c410d8ccSBarry Smith `PetscViewerBinaryGetSkipOptions()`, `PetscViewerBinaryGetSkipInfo()`, `PetscViewerBinaryGetInfoPointer()`
488807ea322SDave May @*/
PetscViewerBinarySetSkipInfo(PetscViewer viewer,PetscBool skip)489d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetSkipInfo(PetscViewer viewer, PetscBool skip)
490d71ae5a4SJacob Faibussowitsch {
491807ea322SDave May PetscFunctionBegin;
492cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
493cc843e7aSLisandro Dalcin PetscValidLogicalCollectiveBool(viewer, skip, 2);
494cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerBinarySetSkipInfo_C", (PetscViewer, PetscBool), (viewer, skip));
4953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
496807ea322SDave May }
497807ea322SDave May
PetscViewerBinarySetSkipInfo_Binary(PetscViewer viewer,PetscBool skip)498d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetSkipInfo_Binary(PetscViewer viewer, PetscBool skip)
499d71ae5a4SJacob Faibussowitsch {
500807ea322SDave May PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
501807ea322SDave May
502807ea322SDave May PetscFunctionBegin;
503cc843e7aSLisandro Dalcin vbinary->skipinfo = skip;
5043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
505807ea322SDave May }
506807ea322SDave May
507807ea322SDave May /*@
508c410d8ccSBarry Smith PetscViewerBinaryGetSkipInfo - check if viewer wrote a `.info` file
509807ea322SDave May
510807ea322SDave May Not Collective
511807ea322SDave May
512807ea322SDave May Input Parameter:
513811af0c4SBarry Smith . viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
514807ea322SDave May
515807ea322SDave May Output Parameter:
516c410d8ccSBarry Smith . skip - `PETSC_TRUE` implies the `.info` file was not generated
517807ea322SDave May
518807ea322SDave May Level: advanced
519807ea322SDave May
520811af0c4SBarry Smith Note:
521811af0c4SBarry Smith This must be called after `PetscViewerSetType()`
522807ea322SDave May
523d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`,
524c410d8ccSBarry Smith `PetscViewerBinarySetSkipOptions()`, `PetscViewerBinarySetSkipInfo()`, `PetscViewerBinaryGetInfoPointer()`
525807ea322SDave May @*/
PetscViewerBinaryGetSkipInfo(PetscViewer viewer,PetscBool * skip)526d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetSkipInfo(PetscViewer viewer, PetscBool *skip)
527d71ae5a4SJacob Faibussowitsch {
528807ea322SDave May PetscFunctionBegin;
529cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
5304f572ea9SToby Isaac PetscAssertPointer(skip, 2);
531cac4c232SBarry Smith PetscUseMethod(viewer, "PetscViewerBinaryGetSkipInfo_C", (PetscViewer, PetscBool *), (viewer, skip));
5323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
533807ea322SDave May }
534807ea322SDave May
PetscViewerBinaryGetSkipInfo_Binary(PetscViewer viewer,PetscBool * skip)535d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetSkipInfo_Binary(PetscViewer viewer, PetscBool *skip)
536d71ae5a4SJacob Faibussowitsch {
537807ea322SDave May PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
538807ea322SDave May
539807ea322SDave May PetscFunctionBegin;
540cc843e7aSLisandro Dalcin *skip = vbinary->skipinfo;
5413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
542807ea322SDave May }
543807ea322SDave May
5445c6c1daeSBarry Smith /*@
545c410d8ccSBarry Smith PetscViewerBinarySetSkipOptions - do not use values in the PETSc options database when loading objects
5465c6c1daeSBarry Smith
5475c6c1daeSBarry Smith Not Collective
5485c6c1daeSBarry Smith
5495c6c1daeSBarry Smith Input Parameters:
550811af0c4SBarry Smith + viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
551811af0c4SBarry Smith - skip - `PETSC_TRUE` means do not use the options from the options database
5525c6c1daeSBarry Smith
5535c6c1daeSBarry Smith Options Database Key:
554811af0c4SBarry Smith . -viewer_binary_skip_options <true or false> - true means do not use the options from the options database
5555c6c1daeSBarry Smith
5565c6c1daeSBarry Smith Level: advanced
5575c6c1daeSBarry Smith
558811af0c4SBarry Smith Note:
559811af0c4SBarry Smith This must be called after `PetscViewerSetType()`
5605c6c1daeSBarry Smith
561d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
562db781477SPatrick Sanan `PetscViewerBinaryGetSkipOptions()`
5635c6c1daeSBarry Smith @*/
PetscViewerBinarySetSkipOptions(PetscViewer viewer,PetscBool skip)564d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetSkipOptions(PetscViewer viewer, PetscBool skip)
565d71ae5a4SJacob Faibussowitsch {
5665c6c1daeSBarry Smith PetscFunctionBegin;
567cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
568cc843e7aSLisandro Dalcin PetscValidLogicalCollectiveBool(viewer, skip, 2);
569cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerBinarySetSkipOptions_C", (PetscViewer, PetscBool), (viewer, skip));
5703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
571807ea322SDave May }
572807ea322SDave May
PetscViewerBinarySetSkipOptions_Binary(PetscViewer viewer,PetscBool skip)573d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetSkipOptions_Binary(PetscViewer viewer, PetscBool skip)
574d71ae5a4SJacob Faibussowitsch {
575cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
576807ea322SDave May
577807ea322SDave May PetscFunctionBegin;
578cc843e7aSLisandro Dalcin vbinary->skipoptions = skip;
5793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5805c6c1daeSBarry Smith }
5815c6c1daeSBarry Smith
5825c6c1daeSBarry Smith /*@
5835c6c1daeSBarry Smith PetscViewerBinaryGetSkipOptions - checks if viewer uses the PETSc options database when loading objects
5845c6c1daeSBarry Smith
5855c6c1daeSBarry Smith Not Collective
5865c6c1daeSBarry Smith
5875c6c1daeSBarry Smith Input Parameter:
588811af0c4SBarry Smith . viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
5895c6c1daeSBarry Smith
5905c6c1daeSBarry Smith Output Parameter:
591811af0c4SBarry Smith . skip - `PETSC_TRUE` means do not use
5925c6c1daeSBarry Smith
5935c6c1daeSBarry Smith Level: advanced
5945c6c1daeSBarry Smith
595811af0c4SBarry Smith Note:
596811af0c4SBarry Smith This must be called after `PetscViewerSetType()`
5975c6c1daeSBarry Smith
598d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
599db781477SPatrick Sanan `PetscViewerBinarySetSkipOptions()`
6005c6c1daeSBarry Smith @*/
PetscViewerBinaryGetSkipOptions(PetscViewer viewer,PetscBool * skip)601d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetSkipOptions(PetscViewer viewer, PetscBool *skip)
602d71ae5a4SJacob Faibussowitsch {
6035c6c1daeSBarry Smith PetscFunctionBegin;
604cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
6054f572ea9SToby Isaac PetscAssertPointer(skip, 2);
606cac4c232SBarry Smith PetscUseMethod(viewer, "PetscViewerBinaryGetSkipOptions_C", (PetscViewer, PetscBool *), (viewer, skip));
6073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6085c6c1daeSBarry Smith }
6095c6c1daeSBarry Smith
PetscViewerBinaryGetSkipOptions_Binary(PetscViewer viewer,PetscBool * skip)610d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetSkipOptions_Binary(PetscViewer viewer, PetscBool *skip)
611d71ae5a4SJacob Faibussowitsch {
612d21b9a37SPierre Jolivet PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
6135c6c1daeSBarry Smith
6145c6c1daeSBarry Smith PetscFunctionBegin;
615cc843e7aSLisandro Dalcin *skip = vbinary->skipoptions;
6163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6175c6c1daeSBarry Smith }
6185c6c1daeSBarry Smith
6195c6c1daeSBarry Smith /*@
6205c6c1daeSBarry Smith PetscViewerBinarySetSkipHeader - do not write a header with size information on output, just raw data
6215c6c1daeSBarry Smith
6225c6c1daeSBarry Smith Not Collective
6235c6c1daeSBarry Smith
6245c6c1daeSBarry Smith Input Parameters:
625811af0c4SBarry Smith + viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
626811af0c4SBarry Smith - skip - `PETSC_TRUE` means do not write header
6275c6c1daeSBarry Smith
6285c6c1daeSBarry Smith Options Database Key:
629811af0c4SBarry Smith . -viewer_binary_skip_header <true or false> - true means do not write header
6305c6c1daeSBarry Smith
6315c6c1daeSBarry Smith Level: advanced
6325c6c1daeSBarry Smith
633c410d8ccSBarry Smith Notes:
634811af0c4SBarry Smith This must be called after `PetscViewerSetType()`
6355c6c1daeSBarry Smith
636c410d8ccSBarry Smith If this option is selected, the output file cannot be read with the `XXXLoad()` such as `VecLoad()`
637c410d8ccSBarry Smith
638d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
639db781477SPatrick Sanan `PetscViewerBinaryGetSkipHeader()`
6405c6c1daeSBarry Smith @*/
PetscViewerBinarySetSkipHeader(PetscViewer viewer,PetscBool skip)641d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinarySetSkipHeader(PetscViewer viewer, PetscBool skip)
642d71ae5a4SJacob Faibussowitsch {
6435c6c1daeSBarry Smith PetscFunctionBegin;
644cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
645cc843e7aSLisandro Dalcin PetscValidLogicalCollectiveBool(viewer, skip, 2);
646cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerBinarySetSkipHeader_C", (PetscViewer, PetscBool), (viewer, skip));
6473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6485c6c1daeSBarry Smith }
6495c6c1daeSBarry Smith
PetscViewerBinarySetSkipHeader_Binary(PetscViewer viewer,PetscBool skip)650d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinarySetSkipHeader_Binary(PetscViewer viewer, PetscBool skip)
651d71ae5a4SJacob Faibussowitsch {
6525c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
6535c6c1daeSBarry Smith
6545c6c1daeSBarry Smith PetscFunctionBegin;
655cc843e7aSLisandro Dalcin vbinary->skipheader = skip;
6563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6575c6c1daeSBarry Smith }
6585c6c1daeSBarry Smith
6595c6c1daeSBarry Smith /*@
6605c6c1daeSBarry Smith PetscViewerBinaryGetSkipHeader - checks whether to write a header with size information on output, or just raw data
6615c6c1daeSBarry Smith
6625c6c1daeSBarry Smith Not Collective
6635c6c1daeSBarry Smith
6645c6c1daeSBarry Smith Input Parameter:
665811af0c4SBarry Smith . viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
6665c6c1daeSBarry Smith
6675c6c1daeSBarry Smith Output Parameter:
668811af0c4SBarry Smith . skip - `PETSC_TRUE` means do not write header
6695c6c1daeSBarry Smith
6705c6c1daeSBarry Smith Level: advanced
6715c6c1daeSBarry Smith
67295452b02SPatrick Sanan Notes:
67395452b02SPatrick Sanan This must be called after PetscViewerSetType()
6745c6c1daeSBarry Smith
675811af0c4SBarry Smith Returns `PETSC_FALSE` for `PETSCSOCKETVIEWER`, you cannot skip the header for it.
6765c6c1daeSBarry Smith
677d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`,
678db781477SPatrick Sanan `PetscViewerBinarySetSkipHeader()`
6795c6c1daeSBarry Smith @*/
PetscViewerBinaryGetSkipHeader(PetscViewer viewer,PetscBool * skip)680d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetSkipHeader(PetscViewer viewer, PetscBool *skip)
681d71ae5a4SJacob Faibussowitsch {
6825c6c1daeSBarry Smith PetscFunctionBegin;
683cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
6844f572ea9SToby Isaac PetscAssertPointer(skip, 2);
685cac4c232SBarry Smith PetscUseMethod(viewer, "PetscViewerBinaryGetSkipHeader_C", (PetscViewer, PetscBool *), (viewer, skip));
6863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6875c6c1daeSBarry Smith }
6885c6c1daeSBarry Smith
PetscViewerBinaryGetSkipHeader_Binary(PetscViewer viewer,PetscBool * skip)689d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetSkipHeader_Binary(PetscViewer viewer, PetscBool *skip)
690d71ae5a4SJacob Faibussowitsch {
691f3b211e4SSatish Balay PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
6925c6c1daeSBarry Smith
6935c6c1daeSBarry Smith PetscFunctionBegin;
694cc843e7aSLisandro Dalcin *skip = vbinary->skipheader;
6953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6965c6c1daeSBarry Smith }
6975c6c1daeSBarry Smith
6985c6c1daeSBarry Smith /*@C
6995c6c1daeSBarry Smith PetscViewerBinaryGetInfoPointer - Extracts the file pointer for the ASCII
700c410d8ccSBarry Smith `.info` file associated with a binary file.
7015c6c1daeSBarry Smith
702cf53795eSBarry Smith Not Collective; No Fortran Support
7035c6c1daeSBarry Smith
7045c6c1daeSBarry Smith Input Parameter:
705811af0c4SBarry Smith . viewer - `PetscViewer` context, obtained from `PetscViewerBinaryOpen()`
7065c6c1daeSBarry Smith
7075c6c1daeSBarry Smith Output Parameter:
708c410d8ccSBarry Smith . file - file pointer Always returns `NULL` if not a binary viewer
7095c6c1daeSBarry Smith
7105c6c1daeSBarry Smith Level: advanced
7115c6c1daeSBarry Smith
712811af0c4SBarry Smith Note:
713811af0c4SBarry Smith For writable binary `PetscViewer`s, the file pointer will only be valid for the
714c410d8ccSBarry Smith first processor in the MPI communicator that shares the `PetscViewer`.
7155c6c1daeSBarry Smith
716c410d8ccSBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinaryGetSkipInfo()`,
717c410d8ccSBarry Smith `PetscViewerBinarySetSkipInfo()`
7185c6c1daeSBarry Smith @*/
PetscViewerBinaryGetInfoPointer(PetscViewer viewer,FILE ** file)719d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryGetInfoPointer(PetscViewer viewer, FILE **file)
720d71ae5a4SJacob Faibussowitsch {
7215c6c1daeSBarry Smith PetscFunctionBegin;
722cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
7234f572ea9SToby Isaac PetscAssertPointer(file, 2);
7240298fd71SBarry Smith *file = NULL;
725cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerBinaryGetInfoPointer_C", (PetscViewer, FILE **), (viewer, file));
7263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7275c6c1daeSBarry Smith }
7285c6c1daeSBarry Smith
PetscViewerBinaryGetInfoPointer_Binary(PetscViewer viewer,FILE ** file)729d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryGetInfoPointer_Binary(PetscViewer viewer, FILE **file)
730d71ae5a4SJacob Faibussowitsch {
731cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
7325c6c1daeSBarry Smith
7335c6c1daeSBarry Smith PetscFunctionBegin;
7349566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer));
735cc843e7aSLisandro Dalcin *file = vbinary->fdes_info;
736cc843e7aSLisandro Dalcin if (viewer->format == PETSC_VIEWER_BINARY_MATLAB && !vbinary->matlabheaderwritten) {
7375c6c1daeSBarry Smith if (vbinary->fdes_info) {
738cc843e7aSLisandro Dalcin FILE *info = vbinary->fdes_info;
7399566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- begin code written by PetscViewerBinary for MATLAB format ---#\n"));
7409566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#$$ Set.filename = '%s';\n", vbinary->filename));
7419566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#$$ fd = PetscOpenFile(Set.filename);\n"));
7429566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n"));
743cc843e7aSLisandro Dalcin }
744cc843e7aSLisandro Dalcin vbinary->matlabheaderwritten = PETSC_TRUE;
7455c6c1daeSBarry Smith }
7463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7475c6c1daeSBarry Smith }
7485c6c1daeSBarry Smith
749e0385b85SDave May #if defined(PETSC_HAVE_MPIIO)
PetscViewerFileClose_BinaryMPIIO(PetscViewer v)750d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_BinaryMPIIO(PetscViewer v)
751d71ae5a4SJacob Faibussowitsch {
752e0385b85SDave May PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
753e0385b85SDave May
754e0385b85SDave May PetscFunctionBegin;
75548a46eb9SPierre Jolivet if (vbinary->mfdes != MPI_FILE_NULL) PetscCallMPI(MPI_File_close(&vbinary->mfdes));
75648a46eb9SPierre Jolivet if (vbinary->mfsub != MPI_FILE_NULL) PetscCallMPI(MPI_File_close(&vbinary->mfsub));
757cc843e7aSLisandro Dalcin vbinary->moff = 0;
7583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
759e0385b85SDave May }
760e0385b85SDave May #endif
761e0385b85SDave May
PetscViewerFileClose_BinarySTDIO(PetscViewer v)762d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_BinarySTDIO(PetscViewer v)
763d71ae5a4SJacob Faibussowitsch {
764cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
765cc843e7aSLisandro Dalcin
766cc843e7aSLisandro Dalcin PetscFunctionBegin;
767cc843e7aSLisandro Dalcin if (vbinary->fdes != -1) {
7689566063dSJacob Faibussowitsch PetscCall(PetscBinaryClose(vbinary->fdes));
769cc843e7aSLisandro Dalcin vbinary->fdes = -1;
770cc843e7aSLisandro Dalcin if (vbinary->storecompressed) {
771cc843e7aSLisandro Dalcin char cmd[8 + PETSC_MAX_PATH_LEN], out[64 + PETSC_MAX_PATH_LEN] = "";
772cc843e7aSLisandro Dalcin const char *gzfilename = vbinary->ogzfilename ? vbinary->ogzfilename : vbinary->filename;
773cc843e7aSLisandro Dalcin /* compress the file */
7749566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(cmd, "gzip -f ", sizeof(cmd)));
7759566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(cmd, gzfilename, sizeof(cmd)));
776cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_POPEN)
777cc843e7aSLisandro Dalcin {
778cc843e7aSLisandro Dalcin FILE *fp;
7799566063dSJacob Faibussowitsch PetscCall(PetscPOpen(PETSC_COMM_SELF, NULL, cmd, "r", &fp));
78000045ab3SPierre Jolivet PetscCheck(!fgets(out, (int)(sizeof(out) - 1), fp), PETSC_COMM_SELF, PETSC_ERR_LIB, "Error from command %s %s", cmd, out);
7819566063dSJacob Faibussowitsch PetscCall(PetscPClose(PETSC_COMM_SELF, fp));
782cc843e7aSLisandro Dalcin }
783cc843e7aSLisandro Dalcin #endif
784cc843e7aSLisandro Dalcin }
785cc843e7aSLisandro Dalcin }
7869566063dSJacob Faibussowitsch PetscCall(PetscFree(vbinary->ogzfilename));
7873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
788cc843e7aSLisandro Dalcin }
789cc843e7aSLisandro Dalcin
PetscViewerFileClose_BinaryInfo(PetscViewer v)790d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_BinaryInfo(PetscViewer v)
791d71ae5a4SJacob Faibussowitsch {
792cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
793cc843e7aSLisandro Dalcin
794cc843e7aSLisandro Dalcin PetscFunctionBegin;
795cc843e7aSLisandro Dalcin if (v->format == PETSC_VIEWER_BINARY_MATLAB && vbinary->matlabheaderwritten) {
796cc843e7aSLisandro Dalcin if (vbinary->fdes_info) {
797cc843e7aSLisandro Dalcin FILE *info = vbinary->fdes_info;
7989566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- begin code written by PetscViewerBinary for MATLAB format ---#\n"));
7999566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#$$ close(fd);\n"));
8009566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n"));
801cc843e7aSLisandro Dalcin }
802cc843e7aSLisandro Dalcin }
803cc843e7aSLisandro Dalcin if (vbinary->fdes_info) {
804cc843e7aSLisandro Dalcin FILE *info = vbinary->fdes_info;
805cc843e7aSLisandro Dalcin vbinary->fdes_info = NULL;
806cc73adaaSBarry Smith PetscCheck(!fclose(info), PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file");
807cc843e7aSLisandro Dalcin }
8083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
809cc843e7aSLisandro Dalcin }
810cc843e7aSLisandro Dalcin
PetscViewerFileClose_Binary(PetscViewer v)811d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_Binary(PetscViewer v)
812d71ae5a4SJacob Faibussowitsch {
813cc843e7aSLisandro Dalcin PetscFunctionBegin;
814cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
8159566063dSJacob Faibussowitsch PetscCall(PetscViewerFileClose_BinaryMPIIO(v));
816cc843e7aSLisandro Dalcin #endif
8179566063dSJacob Faibussowitsch PetscCall(PetscViewerFileClose_BinarySTDIO(v));
8189566063dSJacob Faibussowitsch PetscCall(PetscViewerFileClose_BinaryInfo(v));
8193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
820cc843e7aSLisandro Dalcin }
821cc843e7aSLisandro Dalcin
PetscViewerDestroy_Binary(PetscViewer v)822d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerDestroy_Binary(PetscViewer v)
823d71ae5a4SJacob Faibussowitsch {
8245c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
8255c6c1daeSBarry Smith
8265c6c1daeSBarry Smith PetscFunctionBegin;
8279566063dSJacob Faibussowitsch PetscCall(PetscViewerFileClose_Binary(v));
8289566063dSJacob Faibussowitsch PetscCall(PetscFree(vbinary->filename));
8299566063dSJacob Faibussowitsch PetscCall(PetscFree(vbinary));
8302e956fe4SStefano Zampini PetscCall(PetscViewerBinaryClearFunctionList(v));
8313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
832e0385b85SDave May }
8335c6c1daeSBarry Smith
834cc4c1da9SBarry Smith /*@
8355c6c1daeSBarry Smith PetscViewerBinaryOpen - Opens a file for binary input/output.
8365c6c1daeSBarry Smith
837d083f849SBarry Smith Collective
8385c6c1daeSBarry Smith
8395c6c1daeSBarry Smith Input Parameters:
8405c6c1daeSBarry Smith + comm - MPI communicator
8415c6c1daeSBarry Smith . name - name of file
842cc843e7aSLisandro Dalcin - mode - open mode of file
8433f423023SBarry Smith .vb
8443f423023SBarry Smith FILE_MODE_WRITE - create new file for binary output
8453f423023SBarry Smith FILE_MODE_READ - open existing file for binary input
8463f423023SBarry Smith FILE_MODE_APPEND - open existing file for binary output
8473f423023SBarry Smith .ve
8485c6c1daeSBarry Smith
8495c6c1daeSBarry Smith Output Parameter:
850cc843e7aSLisandro Dalcin . viewer - PetscViewer for binary input/output to use with the specified file
8515c6c1daeSBarry Smith
8525c6c1daeSBarry Smith Options Database Keys:
853811af0c4SBarry Smith + -viewer_binary_filename <name> - name of file to use
854811af0c4SBarry Smith . -viewer_binary_skip_info - true to skip opening an info file
855811af0c4SBarry Smith . -viewer_binary_skip_options - true to not use options database while creating viewer
856811af0c4SBarry Smith . -viewer_binary_skip_header - true to skip output object headers to the file
857811af0c4SBarry Smith - -viewer_binary_mpiio - true to use MPI-IO for input and output to the file (more scalable for large problems)
8585c6c1daeSBarry Smith
8595c6c1daeSBarry Smith Level: beginner
8605c6c1daeSBarry Smith
8615c6c1daeSBarry Smith Note:
862811af0c4SBarry Smith This `PetscViewer` should be destroyed with `PetscViewerDestroy()`.
8635c6c1daeSBarry Smith
8645c6c1daeSBarry Smith For reading files, the filename may begin with ftp:// or http:// and/or
8655c6c1daeSBarry Smith end with .gz; in this case file is brought over and uncompressed.
8665c6c1daeSBarry Smith
8675c6c1daeSBarry Smith For creating files, if the file name ends with .gz it is automatically
8685c6c1daeSBarry Smith compressed when closed.
8695c6c1daeSBarry Smith
870d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
871db781477SPatrick Sanan `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
872db781477SPatrick Sanan `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`, `PetscViewerBinarySetUseMPIIO()`,
873db781477SPatrick Sanan `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()`
8745c6c1daeSBarry Smith @*/
PetscViewerBinaryOpen(MPI_Comm comm,const char name[],PetscFileMode mode,PetscViewer * viewer)875d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryOpen(MPI_Comm comm, const char name[], PetscFileMode mode, PetscViewer *viewer)
876d71ae5a4SJacob Faibussowitsch {
8775c6c1daeSBarry Smith PetscFunctionBegin;
8789566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(comm, viewer));
8799566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(*viewer, PETSCVIEWERBINARY));
8809566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(*viewer, mode));
8819566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(*viewer, name));
8829566063dSJacob Faibussowitsch PetscCall(PetscViewerSetFromOptions(*viewer));
8833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
8845c6c1daeSBarry Smith }
8855c6c1daeSBarry Smith
8865c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
PetscViewerBinaryWriteReadMPIIO(PetscViewer viewer,void * data,PetscInt num,PetscInt * count,PetscDataType dtype,PetscBool write)887d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerBinaryWriteReadMPIIO(PetscViewer viewer, void *data, PetscInt num, PetscInt *count, PetscDataType dtype, PetscBool write)
888d71ae5a4SJacob Faibussowitsch {
88922a8f86cSLisandro Dalcin MPI_Comm comm = PetscObjectComm((PetscObject)viewer);
8905c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
89122a8f86cSLisandro Dalcin MPI_File mfdes = vbinary->mfdes;
8925c6c1daeSBarry Smith MPI_Datatype mdtype;
89369e10bbaSLisandro Dalcin PetscMPIInt rank, cnt;
8945c6c1daeSBarry Smith MPI_Status status;
8955c6c1daeSBarry Smith MPI_Aint ul, dsize;
8965c6c1daeSBarry Smith
8975c6c1daeSBarry Smith PetscFunctionBegin;
8989566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank));
8999566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(num, &cnt));
9009566063dSJacob Faibussowitsch PetscCall(PetscDataTypeToMPIDataType(dtype, &mdtype));
9015c6c1daeSBarry Smith if (write) {
90248a46eb9SPierre Jolivet if (rank == 0) PetscCall(MPIU_File_write_at(mfdes, vbinary->moff, data, cnt, mdtype, &status));
9035c6c1daeSBarry Smith } else {
904dd400576SPatrick Sanan if (rank == 0) {
9059566063dSJacob Faibussowitsch PetscCall(MPIU_File_read_at(mfdes, vbinary->moff, data, cnt, mdtype, &status));
9069566063dSJacob Faibussowitsch if (cnt > 0) PetscCallMPI(MPI_Get_count(&status, mdtype, &cnt));
9075c6c1daeSBarry Smith }
9089566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&cnt, 1, MPI_INT, 0, comm));
9099566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(data, cnt, mdtype, 0, comm));
91069e10bbaSLisandro Dalcin }
9119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_get_extent(mdtype, &ul, &dsize));
9125c6c1daeSBarry Smith vbinary->moff += dsize * cnt;
9139860990eSLisandro Dalcin if (count) *count = cnt;
9143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9155c6c1daeSBarry Smith }
9165c6c1daeSBarry Smith #endif
9175c6c1daeSBarry Smith
9185c6c1daeSBarry Smith /*@C
9195c6c1daeSBarry Smith PetscViewerBinaryRead - Reads from a binary file, all processors get the same result
9205c6c1daeSBarry Smith
92138b83642SBarry Smith Collective; No Fortran Support
9225c6c1daeSBarry Smith
9235c6c1daeSBarry Smith Input Parameters:
924c410d8ccSBarry Smith + viewer - the `PETSCVIEWERBINARY` viewer
925060da220SMatthew G. Knepley . num - number of items of data to read
926907376e6SBarry Smith - dtype - type of data to read
9275c6c1daeSBarry Smith
928c410d8ccSBarry Smith Output Parameters:
929c410d8ccSBarry Smith + data - location of the read data, treated as an array of the type indicated by `dtype`
930c410d8ccSBarry Smith - count - number of items of data actually read, or `NULL`.
931f8e4bde8SMatthew G. Knepley
9325c6c1daeSBarry Smith Level: beginner
9335c6c1daeSBarry Smith
934d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
935db781477SPatrick Sanan `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
93642747ad1SJacob Faibussowitsch `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`
9375c6c1daeSBarry Smith @*/
PetscViewerBinaryRead(PetscViewer viewer,void * data,PetscInt num,PetscInt * count,PetscDataType dtype)938d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryRead(PetscViewer viewer, void *data, PetscInt num, PetscInt *count, PetscDataType dtype)
939d71ae5a4SJacob Faibussowitsch {
94022a8f86cSLisandro Dalcin PetscViewer_Binary *vbinary;
9415c6c1daeSBarry Smith
94222a8f86cSLisandro Dalcin PetscFunctionBegin;
94322a8f86cSLisandro Dalcin PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
94422a8f86cSLisandro Dalcin PetscValidLogicalCollectiveInt(viewer, num, 3);
9459566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer));
94622a8f86cSLisandro Dalcin vbinary = (PetscViewer_Binary *)viewer->data;
9475c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
948bc196f7cSDave May if (vbinary->usempiio) {
9499566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWriteReadMPIIO(viewer, data, num, count, dtype, PETSC_FALSE));
9505c6c1daeSBarry Smith } else {
9515c6c1daeSBarry Smith #endif
9529566063dSJacob Faibussowitsch PetscCall(PetscBinarySynchronizedRead(PetscObjectComm((PetscObject)viewer), vbinary->fdes, data, num, count, dtype));
9535c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
9545c6c1daeSBarry Smith }
9555c6c1daeSBarry Smith #endif
9563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9575c6c1daeSBarry Smith }
9585c6c1daeSBarry Smith
9595c6c1daeSBarry Smith /*@C
960811af0c4SBarry Smith PetscViewerBinaryWrite - writes to a binary file, only from the first MPI rank
9615c6c1daeSBarry Smith
96238b83642SBarry Smith Collective; No Fortran Support
9635c6c1daeSBarry Smith
9645c6c1daeSBarry Smith Input Parameters:
965c410d8ccSBarry Smith + viewer - the `PETSCVIEWERBINARY` viewer
966c410d8ccSBarry Smith . data - location of data, treated as an array of the type indicated by `dtype`
9675824c9d0SPatrick Sanan . count - number of items of data to write
968f253e43cSLisandro Dalcin - dtype - type of data to write
9695c6c1daeSBarry Smith
9705c6c1daeSBarry Smith Level: beginner
9715c6c1daeSBarry Smith
972d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
973db781477SPatrick Sanan `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`, `PetscDataType`
974db781477SPatrick Sanan `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
9755c6c1daeSBarry Smith @*/
PetscViewerBinaryWrite(PetscViewer viewer,const void * data,PetscInt count,PetscDataType dtype)976d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryWrite(PetscViewer viewer, const void *data, PetscInt count, PetscDataType dtype)
977d71ae5a4SJacob Faibussowitsch {
97822a8f86cSLisandro Dalcin PetscViewer_Binary *vbinary;
9795c6c1daeSBarry Smith
9805c6c1daeSBarry Smith PetscFunctionBegin;
98122a8f86cSLisandro Dalcin PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
98222a8f86cSLisandro Dalcin PetscValidLogicalCollectiveInt(viewer, count, 3);
9839566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer));
98422a8f86cSLisandro Dalcin vbinary = (PetscViewer_Binary *)viewer->data;
9855c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
986bc196f7cSDave May if (vbinary->usempiio) {
9879566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWriteReadMPIIO(viewer, (void *)data, count, NULL, dtype, PETSC_TRUE));
9885c6c1daeSBarry Smith } else {
9895c6c1daeSBarry Smith #endif
9909566063dSJacob Faibussowitsch PetscCall(PetscBinarySynchronizedWrite(PetscObjectComm((PetscObject)viewer), vbinary->fdes, data, count, dtype));
9915c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
9925c6c1daeSBarry Smith }
9935c6c1daeSBarry Smith #endif
9943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9955c6c1daeSBarry Smith }
9965c6c1daeSBarry Smith
PetscViewerBinaryWriteReadAll(PetscViewer viewer,PetscBool write,void * data,PetscCount count,PetscCount start,PetscCount total,PetscDataType dtype)9976497c311SBarry Smith static PetscErrorCode PetscViewerBinaryWriteReadAll(PetscViewer viewer, PetscBool write, void *data, PetscCount count, PetscCount start, PetscCount total, PetscDataType dtype)
998d71ae5a4SJacob Faibussowitsch {
9995972f5f3SLisandro Dalcin MPI_Comm comm = PetscObjectComm((PetscObject)viewer);
10005972f5f3SLisandro Dalcin PetscMPIInt size, rank;
10015972f5f3SLisandro Dalcin MPI_Datatype mdtype;
1002ec4bef21SJose E. Roman PETSC_UNUSED MPI_Aint lb;
1003ec4bef21SJose E. Roman MPI_Aint dsize;
10045972f5f3SLisandro Dalcin PetscBool useMPIIO;
10055972f5f3SLisandro Dalcin
10065972f5f3SLisandro Dalcin PetscFunctionBegin;
10075972f5f3SLisandro Dalcin PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 1, PETSCVIEWERBINARY);
100857508eceSPierre Jolivet PetscValidLogicalCollectiveBool(viewer, (start >= 0) || (start == PETSC_DETERMINE), 5);
100957508eceSPierre Jolivet PetscValidLogicalCollectiveBool(viewer, (total >= 0) || (total == PETSC_DETERMINE), 6);
10106497c311SBarry Smith PetscValidLogicalCollectiveCount(viewer, total, 6);
10119566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer));
10125972f5f3SLisandro Dalcin
10139566063dSJacob Faibussowitsch PetscCall(PetscDataTypeToMPIDataType(dtype, &mdtype));
10149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_get_extent(mdtype, &lb, &dsize));
10159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank));
10169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size));
10175972f5f3SLisandro Dalcin
10189566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetUseMPIIO(viewer, &useMPIIO));
10195972f5f3SLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
10205972f5f3SLisandro Dalcin if (useMPIIO) {
10215972f5f3SLisandro Dalcin MPI_File mfdes;
10225972f5f3SLisandro Dalcin MPI_Offset off;
10235972f5f3SLisandro Dalcin PetscMPIInt cnt;
10245972f5f3SLisandro Dalcin
10255972f5f3SLisandro Dalcin if (start == PETSC_DETERMINE) {
10266497c311SBarry Smith PetscCallMPI(MPI_Scan(&count, &start, 1, MPIU_COUNT, MPI_SUM, comm));
10275972f5f3SLisandro Dalcin start -= count;
10285972f5f3SLisandro Dalcin }
10295972f5f3SLisandro Dalcin if (total == PETSC_DETERMINE) {
10305972f5f3SLisandro Dalcin total = start + count;
10316497c311SBarry Smith PetscCallMPI(MPI_Bcast(&total, 1, MPIU_COUNT, size - 1, comm));
10325972f5f3SLisandro Dalcin }
10339566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(count, &cnt));
10349566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetMPIIODescriptor(viewer, &mfdes));
10359566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetMPIIOOffset(viewer, &off));
10365972f5f3SLisandro Dalcin off += (MPI_Offset)(start * dsize);
10375972f5f3SLisandro Dalcin if (write) {
10389566063dSJacob Faibussowitsch PetscCall(MPIU_File_write_at_all(mfdes, off, data, cnt, mdtype, MPI_STATUS_IGNORE));
10395972f5f3SLisandro Dalcin } else {
10409566063dSJacob Faibussowitsch PetscCall(MPIU_File_read_at_all(mfdes, off, data, cnt, mdtype, MPI_STATUS_IGNORE));
10415972f5f3SLisandro Dalcin }
10425972f5f3SLisandro Dalcin off = (MPI_Offset)(total * dsize);
10439566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryAddMPIIOOffset(viewer, off));
10443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
10455972f5f3SLisandro Dalcin }
10465972f5f3SLisandro Dalcin #endif
10475972f5f3SLisandro Dalcin {
10485972f5f3SLisandro Dalcin int fdes;
10495972f5f3SLisandro Dalcin char *workbuf = NULL;
10506497c311SBarry Smith PetscCount tcount = rank == 0 ? 0 : count, maxcount = 0;
10516497c311SBarry Smith PetscInt message_count, flowcontrolcount;
10525972f5f3SLisandro Dalcin PetscMPIInt tag, cnt, maxcnt, scnt = 0, rcnt = 0, j;
10535972f5f3SLisandro Dalcin MPI_Status status;
10545972f5f3SLisandro Dalcin
10559566063dSJacob Faibussowitsch PetscCall(PetscCommGetNewTag(comm, &tag));
10566497c311SBarry Smith PetscCallMPI(MPI_Reduce(&tcount, &maxcount, 1, MPIU_COUNT, MPI_MAX, 0, comm));
10579566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(maxcount, &maxcnt));
10585972f5f3SLisandro Dalcin
10599566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetDescriptor(viewer, &fdes));
10602fb3221eSBarry Smith if (rank == 0) {
10612fb3221eSBarry Smith if (write) PetscCall(PetscBinaryWrite(fdes, data, count, dtype));
10622fb3221eSBarry Smith else PetscCall(PetscBinaryRead(fdes, data, count, NULL, dtype));
10632fb3221eSBarry Smith }
10642fb3221eSBarry Smith
10652fb3221eSBarry Smith if (size > 1) {
10669566063dSJacob Faibussowitsch PetscCall(PetscViewerFlowControlStart(viewer, &message_count, &flowcontrolcount));
1067dd400576SPatrick Sanan if (rank == 0) {
10689566063dSJacob Faibussowitsch PetscCall(PetscMalloc(maxcnt * dsize, &workbuf));
10695972f5f3SLisandro Dalcin for (j = 1; j < size; j++) {
10709566063dSJacob Faibussowitsch PetscCall(PetscViewerFlowControlStepMain(viewer, j, &message_count, flowcontrolcount));
10715972f5f3SLisandro Dalcin if (write) {
10729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(workbuf, maxcnt, mdtype, j, tag, comm, &status));
10739566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status, mdtype, &rcnt));
10749566063dSJacob Faibussowitsch PetscCall(PetscBinaryWrite(fdes, workbuf, rcnt, dtype));
10755972f5f3SLisandro Dalcin } else {
10769566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(&scnt, 1, MPI_INT, j, tag, comm, MPI_STATUS_IGNORE));
10779566063dSJacob Faibussowitsch PetscCall(PetscBinaryRead(fdes, workbuf, scnt, NULL, dtype));
10789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(workbuf, scnt, mdtype, j, tag, comm));
10795972f5f3SLisandro Dalcin }
10805972f5f3SLisandro Dalcin }
10819566063dSJacob Faibussowitsch PetscCall(PetscFree(workbuf));
10829566063dSJacob Faibussowitsch PetscCall(PetscViewerFlowControlEndMain(viewer, &message_count));
10835972f5f3SLisandro Dalcin } else {
10842fb3221eSBarry Smith PetscCall(PetscMPIIntCast(count, &cnt));
10859566063dSJacob Faibussowitsch PetscCall(PetscViewerFlowControlStepWorker(viewer, rank, &message_count));
10865972f5f3SLisandro Dalcin if (write) {
10879566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(data, cnt, mdtype, 0, tag, comm));
10885972f5f3SLisandro Dalcin } else {
10899566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(&cnt, 1, MPI_INT, 0, tag, comm));
10909566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(data, cnt, mdtype, 0, tag, comm, MPI_STATUS_IGNORE));
10915972f5f3SLisandro Dalcin }
10929566063dSJacob Faibussowitsch PetscCall(PetscViewerFlowControlEndWorker(viewer, &message_count));
10935972f5f3SLisandro Dalcin }
10945972f5f3SLisandro Dalcin }
10952fb3221eSBarry Smith }
10963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
10975972f5f3SLisandro Dalcin }
10985972f5f3SLisandro Dalcin
10995972f5f3SLisandro Dalcin /*@C
1100c410d8ccSBarry Smith PetscViewerBinaryReadAll - reads from a binary file from all MPI processes, each rank receives its own portion of the data
11015972f5f3SLisandro Dalcin
110238b83642SBarry Smith Collective; No Fortran Support
11035972f5f3SLisandro Dalcin
11045972f5f3SLisandro Dalcin Input Parameters:
1105c410d8ccSBarry Smith + viewer - the `PETSCVIEWERBINARY` viewer
11065972f5f3SLisandro Dalcin . count - local number of items of data to read
1107811af0c4SBarry Smith . start - local start, can be `PETSC_DETERMINE`
1108811af0c4SBarry Smith . total - global number of items of data to read, can be `PETSC_DETERMINE`
11095972f5f3SLisandro Dalcin - dtype - type of data to read
11105972f5f3SLisandro Dalcin
1111c410d8ccSBarry Smith Output Parameter:
1112c410d8ccSBarry Smith . data - location of data, treated as an array of type indicated by `dtype`
1113c410d8ccSBarry Smith
11145972f5f3SLisandro Dalcin Level: advanced
11155972f5f3SLisandro Dalcin
1116d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryRead()`, `PetscViewerBinaryWriteAll()`
11175972f5f3SLisandro Dalcin @*/
PetscViewerBinaryReadAll(PetscViewer viewer,void * data,PetscCount count,PetscCount start,PetscCount total,PetscDataType dtype)11183e1d7bceSPierre Jolivet PetscErrorCode PetscViewerBinaryReadAll(PetscViewer viewer, void *data, PetscCount count, PetscCount start, PetscCount total, PetscDataType dtype)
1119d71ae5a4SJacob Faibussowitsch {
11205972f5f3SLisandro Dalcin PetscFunctionBegin;
11219566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWriteReadAll(viewer, PETSC_FALSE, data, count, start, total, dtype));
11223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11235972f5f3SLisandro Dalcin }
11245972f5f3SLisandro Dalcin
11255972f5f3SLisandro Dalcin /*@C
1126c410d8ccSBarry Smith PetscViewerBinaryWriteAll - writes to a binary file from all MPI processes, each rank writes its own portion of the data
11275972f5f3SLisandro Dalcin
112838b83642SBarry Smith Collective; No Fortran Support
11295972f5f3SLisandro Dalcin
11305972f5f3SLisandro Dalcin Input Parameters:
1131c410d8ccSBarry Smith + viewer - the `PETSCVIEWERBINARY` viewer
11325972f5f3SLisandro Dalcin . data - location of data
1133c410d8ccSBarry Smith . count - local number of items of data to write, treated as an array of type indicated by `dtype`
1134811af0c4SBarry Smith . start - local start, can be `PETSC_DETERMINE`
1135811af0c4SBarry Smith . total - global number of items of data to write, can be `PETSC_DETERMINE`
11365972f5f3SLisandro Dalcin - dtype - type of data to write
11375972f5f3SLisandro Dalcin
11385972f5f3SLisandro Dalcin Level: advanced
11395972f5f3SLisandro Dalcin
114042747ad1SJacob Faibussowitsch .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerBinaryOpen()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryReadAll()`
11415972f5f3SLisandro Dalcin @*/
PetscViewerBinaryWriteAll(PetscViewer viewer,const void * data,PetscCount count,PetscCount start,PetscCount total,PetscDataType dtype)11426497c311SBarry Smith PetscErrorCode PetscViewerBinaryWriteAll(PetscViewer viewer, const void *data, PetscCount count, PetscCount start, PetscCount total, PetscDataType dtype)
1143d71ae5a4SJacob Faibussowitsch {
11445972f5f3SLisandro Dalcin PetscFunctionBegin;
11459566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWriteReadAll(viewer, PETSC_TRUE, (void *)data, count, start, total, dtype));
11463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11475972f5f3SLisandro Dalcin }
11485972f5f3SLisandro Dalcin
11495c6c1daeSBarry Smith /*@C
1150811af0c4SBarry Smith PetscViewerBinaryWriteStringArray - writes to a binary file, only from the first MPI rank, an array of strings
11515c6c1daeSBarry Smith
115238b83642SBarry Smith Collective; No Fortran Support
11535c6c1daeSBarry Smith
11545c6c1daeSBarry Smith Input Parameters:
1155c410d8ccSBarry Smith + viewer - the `PETSCVIEWERBINARY` viewer
11565c6c1daeSBarry Smith - data - location of the array of strings
11575c6c1daeSBarry Smith
11585c6c1daeSBarry Smith Level: intermediate
11595c6c1daeSBarry Smith
1160811af0c4SBarry Smith Note:
11613f423023SBarry Smith The array of strings must be `NULL` terminated
11625c6c1daeSBarry Smith
1163d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
1164db781477SPatrick Sanan `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
1165db781477SPatrick Sanan `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
11665c6c1daeSBarry Smith @*/
PetscViewerBinaryWriteStringArray(PetscViewer viewer,const char * const data[])1167cc4c1da9SBarry Smith PetscErrorCode PetscViewerBinaryWriteStringArray(PetscViewer viewer, const char *const data[])
1168d71ae5a4SJacob Faibussowitsch {
11695c6c1daeSBarry Smith PetscInt i, n = 0, *sizes;
1170cc843e7aSLisandro Dalcin size_t len;
11715c6c1daeSBarry Smith
117222a8f86cSLisandro Dalcin PetscFunctionBegin;
11739566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer));
11745c6c1daeSBarry Smith /* count number of strings */
1175fbccb6d4SPierre Jolivet while (data[n++]);
11765c6c1daeSBarry Smith n--;
11779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n + 1, &sizes));
11785c6c1daeSBarry Smith sizes[0] = n;
11795c6c1daeSBarry Smith for (i = 0; i < n; i++) {
11809566063dSJacob Faibussowitsch PetscCall(PetscStrlen(data[i], &len));
1181cc843e7aSLisandro Dalcin sizes[i + 1] = (PetscInt)len + 1; /* size includes space for the null terminator */
11825c6c1daeSBarry Smith }
11839566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWrite(viewer, sizes, n + 1, PETSC_INT));
118448a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall(PetscViewerBinaryWrite(viewer, (void *)data[i], sizes[i + 1], PETSC_CHAR));
11859566063dSJacob Faibussowitsch PetscCall(PetscFree(sizes));
11863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11875c6c1daeSBarry Smith }
11885c6c1daeSBarry Smith
11895c6c1daeSBarry Smith /*@C
1190c410d8ccSBarry Smith PetscViewerBinaryReadStringArray - reads a binary file an array of strings to all MPI processes
11915c6c1daeSBarry Smith
119238b83642SBarry Smith Collective; No Fortran Support
11935c6c1daeSBarry Smith
11945c6c1daeSBarry Smith Input Parameter:
1195c410d8ccSBarry Smith . viewer - the `PETSCVIEWERBINARY` viewer
11965c6c1daeSBarry Smith
11975c6c1daeSBarry Smith Output Parameter:
11985c6c1daeSBarry Smith . data - location of the array of strings
11995c6c1daeSBarry Smith
12005c6c1daeSBarry Smith Level: intermediate
12015c6c1daeSBarry Smith
1202811af0c4SBarry Smith Note:
12033f423023SBarry Smith The array of strings must `NULL` terminated
12045c6c1daeSBarry Smith
1205d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
1206db781477SPatrick Sanan `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`,
1207db781477SPatrick Sanan `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`
12085c6c1daeSBarry Smith @*/
PetscViewerBinaryReadStringArray(PetscViewer viewer,char *** data)1209d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerBinaryReadStringArray(PetscViewer viewer, char ***data)
1210d71ae5a4SJacob Faibussowitsch {
1211060da220SMatthew G. Knepley PetscInt i, n, *sizes, N = 0;
12125c6c1daeSBarry Smith
121322a8f86cSLisandro Dalcin PetscFunctionBegin;
12149566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer));
12155c6c1daeSBarry Smith /* count number of strings */
12169566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &n, 1, NULL, PETSC_INT));
12179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &sizes));
12189566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, sizes, n, NULL, PETSC_INT));
1219a297a907SKarl Rupp for (i = 0; i < n; i++) N += sizes[i];
12209566063dSJacob Faibussowitsch PetscCall(PetscMalloc((n + 1) * sizeof(char *) + N * sizeof(char), data));
12215c6c1daeSBarry Smith (*data)[0] = (char *)((*data) + n + 1);
1222a297a907SKarl Rupp for (i = 1; i < n; i++) (*data)[i] = (*data)[i - 1] + sizes[i - 1];
12239566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, (*data)[0], N, NULL, PETSC_CHAR));
122402c9f0b5SLisandro Dalcin (*data)[n] = NULL;
12259566063dSJacob Faibussowitsch PetscCall(PetscFree(sizes));
12263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
12275c6c1daeSBarry Smith }
12285c6c1daeSBarry Smith
12295d83a8b1SBarry Smith /*@
1230cc843e7aSLisandro Dalcin PetscViewerFileSetMode - Sets the open mode of file
1231cc843e7aSLisandro Dalcin
1232c3339decSBarry Smith Logically Collective
1233cc843e7aSLisandro Dalcin
1234cc843e7aSLisandro Dalcin Input Parameters:
123515229ffcSPierre Jolivet + viewer - the `PetscViewer`; must be a `PETSCVIEWERBINARY`, `PETSCVIEWERMATLAB`, `PETSCVIEWERHDF5`, or `PETSCVIEWERASCII` `PetscViewer`
1236cc843e7aSLisandro Dalcin - mode - open mode of file
1237cf53795eSBarry Smith .vb
1238cf53795eSBarry Smith FILE_MODE_WRITE - create new file for output
1239cf53795eSBarry Smith FILE_MODE_READ - open existing file for input
1240cf53795eSBarry Smith FILE_MODE_APPEND - open existing file for output
1241cf53795eSBarry Smith .ve
1242cc843e7aSLisandro Dalcin
1243cc843e7aSLisandro Dalcin Level: advanced
1244cc843e7aSLisandro Dalcin
124542747ad1SJacob Faibussowitsch .seealso: [](sec_viewers), `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
1246cc843e7aSLisandro Dalcin @*/
PetscViewerFileSetMode(PetscViewer viewer,PetscFileMode mode)1247d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerFileSetMode(PetscViewer viewer, PetscFileMode mode)
1248d71ae5a4SJacob Faibussowitsch {
1249cc843e7aSLisandro Dalcin PetscFunctionBegin;
1250cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1251cc843e7aSLisandro Dalcin PetscValidLogicalCollectiveEnum(viewer, mode, 2);
125208401ef6SPierre Jolivet PetscCheck(mode != FILE_MODE_UNDEFINED, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Cannot set FILE_MODE_UNDEFINED");
1253f7d195e4SLawrence Mitchell PetscCheck(mode >= FILE_MODE_UNDEFINED && mode <= FILE_MODE_APPEND_UPDATE, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_OUTOFRANGE, "Invalid file mode %d", (int)mode);
1254cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerFileSetMode_C", (PetscViewer, PetscFileMode), (viewer, mode));
12553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1256cc843e7aSLisandro Dalcin }
1257cc843e7aSLisandro Dalcin
PetscViewerFileSetMode_Binary(PetscViewer viewer,PetscFileMode mode)1258d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetMode_Binary(PetscViewer viewer, PetscFileMode mode)
1259d71ae5a4SJacob Faibussowitsch {
1260cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1261cc843e7aSLisandro Dalcin
1262cc843e7aSLisandro Dalcin PetscFunctionBegin;
1263cc73adaaSBarry Smith PetscCheck(!viewer->setupcalled || vbinary->filemode == mode, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Cannot change mode to %s after setup", PetscFileModes[mode]);
1264cc843e7aSLisandro Dalcin vbinary->filemode = mode;
12653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1266cc843e7aSLisandro Dalcin }
1267cc843e7aSLisandro Dalcin
12685d83a8b1SBarry Smith /*@
1269c410d8ccSBarry Smith PetscViewerFileGetMode - Gets the open mode of a file associated with a `PetscViewer`
1270cc843e7aSLisandro Dalcin
1271cc843e7aSLisandro Dalcin Not Collective
1272cc843e7aSLisandro Dalcin
1273cc843e7aSLisandro Dalcin Input Parameter:
1274811af0c4SBarry Smith . viewer - the `PetscViewer`; must be a `PETSCVIEWERBINARY`, `PETSCVIEWERMATLAB`, `PETSCVIEWERHDF5`, or `PETSCVIEWERASCII` `PetscViewer`
1275cc843e7aSLisandro Dalcin
1276cc843e7aSLisandro Dalcin Output Parameter:
1277cc843e7aSLisandro Dalcin . mode - open mode of file
1278cf53795eSBarry Smith .vb
1279cf53795eSBarry Smith FILE_MODE_WRITE - create new file for binary output
1280cf53795eSBarry Smith FILE_MODE_READ - open existing file for binary input
1281cf53795eSBarry Smith FILE_MODE_APPEND - open existing file for binary output
1282cf53795eSBarry Smith .ve
1283cc843e7aSLisandro Dalcin
1284cc843e7aSLisandro Dalcin Level: advanced
1285cc843e7aSLisandro Dalcin
1286d1f92df0SBarry Smith .seealso: [](sec_viewers), `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
1287cc843e7aSLisandro Dalcin @*/
PetscViewerFileGetMode(PetscViewer viewer,PetscFileMode * mode)1288d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerFileGetMode(PetscViewer viewer, PetscFileMode *mode)
1289d71ae5a4SJacob Faibussowitsch {
1290cc843e7aSLisandro Dalcin PetscFunctionBegin;
1291cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
12924f572ea9SToby Isaac PetscAssertPointer(mode, 2);
1293cac4c232SBarry Smith PetscUseMethod(viewer, "PetscViewerFileGetMode_C", (PetscViewer, PetscFileMode *), (viewer, mode));
12943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1295cc843e7aSLisandro Dalcin }
1296cc843e7aSLisandro Dalcin
PetscViewerFileGetMode_Binary(PetscViewer viewer,PetscFileMode * mode)1297d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetMode_Binary(PetscViewer viewer, PetscFileMode *mode)
1298d71ae5a4SJacob Faibussowitsch {
1299cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1300cc843e7aSLisandro Dalcin
1301cc843e7aSLisandro Dalcin PetscFunctionBegin;
1302cc843e7aSLisandro Dalcin *mode = vbinary->filemode;
13033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1304cc843e7aSLisandro Dalcin }
1305cc843e7aSLisandro Dalcin
PetscViewerFileSetName_Binary(PetscViewer viewer,const char name[])1306d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetName_Binary(PetscViewer viewer, const char name[])
1307d71ae5a4SJacob Faibussowitsch {
1308cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1309cc843e7aSLisandro Dalcin
1310cc843e7aSLisandro Dalcin PetscFunctionBegin;
1311cc843e7aSLisandro Dalcin if (viewer->setupcalled && vbinary->filename) {
1312cc843e7aSLisandro Dalcin /* gzip can be run after the file with the previous filename has been closed */
13139566063dSJacob Faibussowitsch PetscCall(PetscFree(vbinary->ogzfilename));
13149566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(vbinary->filename, &vbinary->ogzfilename));
1315cc843e7aSLisandro Dalcin }
13169566063dSJacob Faibussowitsch PetscCall(PetscFree(vbinary->filename));
13179566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &vbinary->filename));
1318cc843e7aSLisandro Dalcin viewer->setupcalled = PETSC_FALSE;
13193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1320cc843e7aSLisandro Dalcin }
1321cc843e7aSLisandro Dalcin
PetscViewerFileGetName_Binary(PetscViewer viewer,const char ** name)1322d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetName_Binary(PetscViewer viewer, const char **name)
1323d71ae5a4SJacob Faibussowitsch {
13245c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
13255c6c1daeSBarry Smith
13265c6c1daeSBarry Smith PetscFunctionBegin;
13275c6c1daeSBarry Smith *name = vbinary->filename;
13283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
13295c6c1daeSBarry Smith }
13305c6c1daeSBarry Smith
13315c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
PetscViewerFileSetUp_BinaryMPIIO(PetscViewer viewer)1332d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetUp_BinaryMPIIO(PetscViewer viewer)
1333d71ae5a4SJacob Faibussowitsch {
13345c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1335e8a65b7dSLisandro Dalcin int amode;
13365c6c1daeSBarry Smith
13375c6c1daeSBarry Smith PetscFunctionBegin;
1338a297a907SKarl Rupp vbinary->storecompressed = PETSC_FALSE;
13395c6c1daeSBarry Smith
1340cc843e7aSLisandro Dalcin vbinary->moff = 0;
1341cc843e7aSLisandro Dalcin switch (vbinary->filemode) {
1342d71ae5a4SJacob Faibussowitsch case FILE_MODE_READ:
1343d71ae5a4SJacob Faibussowitsch amode = MPI_MODE_RDONLY;
1344d71ae5a4SJacob Faibussowitsch break;
1345d71ae5a4SJacob Faibussowitsch case FILE_MODE_WRITE:
1346d71ae5a4SJacob Faibussowitsch amode = MPI_MODE_WRONLY | MPI_MODE_CREATE;
1347d71ae5a4SJacob Faibussowitsch break;
1348d71ae5a4SJacob Faibussowitsch case FILE_MODE_APPEND:
1349d71ae5a4SJacob Faibussowitsch amode = MPI_MODE_WRONLY | MPI_MODE_CREATE | MPI_MODE_APPEND;
1350d71ae5a4SJacob Faibussowitsch break;
1351d71ae5a4SJacob Faibussowitsch case FILE_MODE_UNDEFINED:
1352d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerSetUp()");
1353d71ae5a4SJacob Faibussowitsch default:
1354d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[vbinary->filemode]);
13555c6c1daeSBarry Smith }
13569566063dSJacob Faibussowitsch PetscCallMPI(MPI_File_open(PetscObjectComm((PetscObject)viewer), vbinary->filename, amode, MPI_INFO_NULL, &vbinary->mfdes));
135722a8f86cSLisandro Dalcin /*
135822a8f86cSLisandro Dalcin The MPI standard does not have MPI_MODE_TRUNCATE. We emulate this behavior by setting the file size to zero.
135922a8f86cSLisandro Dalcin */
13609566063dSJacob Faibussowitsch if (vbinary->filemode == FILE_MODE_WRITE) PetscCallMPI(MPI_File_set_size(vbinary->mfdes, 0));
136122a8f86cSLisandro Dalcin /*
136222a8f86cSLisandro Dalcin Initially, all processes view the file as a linear byte stream. Therefore, for files opened with MPI_MODE_APPEND,
136322a8f86cSLisandro Dalcin MPI_File_get_position[_shared](fh, &offset) returns the absolute byte position at the end of file.
136422a8f86cSLisandro Dalcin Otherwise, we would need to call MPI_File_get_byte_offset(fh, offset, &byte_offset) to convert
136522a8f86cSLisandro Dalcin the offset in etype units to an absolute byte position.
136622a8f86cSLisandro Dalcin */
13679566063dSJacob Faibussowitsch if (vbinary->filemode == FILE_MODE_APPEND) PetscCallMPI(MPI_File_get_position(vbinary->mfdes, &vbinary->moff));
13683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1369cc843e7aSLisandro Dalcin }
1370cc843e7aSLisandro Dalcin #endif
13715c6c1daeSBarry Smith
PetscViewerFileSetUp_BinarySTDIO(PetscViewer viewer)1372d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetUp_BinarySTDIO(PetscViewer viewer)
1373d71ae5a4SJacob Faibussowitsch {
1374cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1375cc843e7aSLisandro Dalcin const char *fname;
1376bbcf679cSJacob Faibussowitsch char bname[PETSC_MAX_PATH_LEN], *gz = NULL;
1377cc843e7aSLisandro Dalcin PetscBool found;
1378cc843e7aSLisandro Dalcin PetscMPIInt rank;
13795c6c1daeSBarry Smith
1380cc843e7aSLisandro Dalcin PetscFunctionBegin;
13819566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
13825c6c1daeSBarry Smith
1383cc843e7aSLisandro Dalcin /* if file name ends in .gz strip that off and note user wants file compressed */
1384cc843e7aSLisandro Dalcin vbinary->storecompressed = PETSC_FALSE;
1385cc843e7aSLisandro Dalcin if (vbinary->filemode == FILE_MODE_WRITE) {
13869566063dSJacob Faibussowitsch PetscCall(PetscStrstr(vbinary->filename, ".gz", &gz));
13879371c9d4SSatish Balay if (gz && gz[3] == 0) {
13889371c9d4SSatish Balay *gz = 0;
13899371c9d4SSatish Balay vbinary->storecompressed = PETSC_TRUE;
13909371c9d4SSatish Balay }
1391cc843e7aSLisandro Dalcin }
1392cc843e7aSLisandro Dalcin #if !defined(PETSC_HAVE_POPEN)
139328b400f6SJacob Faibussowitsch PetscCheck(!vbinary->storecompressed, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP_SYS, "Cannot run gzip on this machine");
1394cc843e7aSLisandro Dalcin #endif
1395cc843e7aSLisandro Dalcin
1396cc843e7aSLisandro Dalcin fname = vbinary->filename;
1397cc843e7aSLisandro Dalcin if (vbinary->filemode == FILE_MODE_READ) { /* possibly get the file from remote site or compressed file */
13989566063dSJacob Faibussowitsch PetscCall(PetscFileRetrieve(PetscObjectComm((PetscObject)viewer), fname, bname, PETSC_MAX_PATH_LEN, &found));
139928b400f6SJacob Faibussowitsch PetscCheck(found, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_OPEN, "Cannot locate file: %s", fname);
1400cc843e7aSLisandro Dalcin fname = bname;
14015c6c1daeSBarry Smith }
14025c6c1daeSBarry Smith
1403cc843e7aSLisandro Dalcin vbinary->fdes = -1;
1404dd400576SPatrick Sanan if (rank == 0) { /* only first processor opens file*/
1405cc843e7aSLisandro Dalcin PetscFileMode mode = vbinary->filemode;
1406cc843e7aSLisandro Dalcin if (mode == FILE_MODE_APPEND) {
1407cc843e7aSLisandro Dalcin /* check if asked to append to a non-existing file */
14089566063dSJacob Faibussowitsch PetscCall(PetscTestFile(fname, '\0', &found));
1409cc843e7aSLisandro Dalcin if (!found) mode = FILE_MODE_WRITE;
1410cc843e7aSLisandro Dalcin }
14119566063dSJacob Faibussowitsch PetscCall(PetscBinaryOpen(fname, mode, &vbinary->fdes));
1412cc843e7aSLisandro Dalcin }
14133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1414cc843e7aSLisandro Dalcin }
1415cc843e7aSLisandro Dalcin
1416bf31d2d3SBarry Smith #include <errno.h>
PetscViewerFileSetUp_BinaryInfo(PetscViewer viewer)1417d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetUp_BinaryInfo(PetscViewer viewer)
1418d71ae5a4SJacob Faibussowitsch {
1419cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1420cc843e7aSLisandro Dalcin PetscMPIInt rank;
1421cc843e7aSLisandro Dalcin PetscBool found;
1422cc843e7aSLisandro Dalcin
1423cc843e7aSLisandro Dalcin PetscFunctionBegin;
1424cc843e7aSLisandro Dalcin vbinary->fdes_info = NULL;
14259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
1426dd400576SPatrick Sanan if (!vbinary->skipinfo && (vbinary->filemode == FILE_MODE_READ || rank == 0)) {
1427cc843e7aSLisandro Dalcin char infoname[PETSC_MAX_PATH_LEN], iname[PETSC_MAX_PATH_LEN], *gz;
1428cc843e7aSLisandro Dalcin
14299566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(infoname, vbinary->filename, sizeof(infoname)));
1430cc843e7aSLisandro Dalcin /* remove .gz if it ends file name */
14319566063dSJacob Faibussowitsch PetscCall(PetscStrstr(infoname, ".gz", &gz));
1432cc843e7aSLisandro Dalcin if (gz && gz[3] == 0) *gz = 0;
1433cc843e7aSLisandro Dalcin
14349566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(infoname, ".info", sizeof(infoname)));
1435cc843e7aSLisandro Dalcin if (vbinary->filemode == FILE_MODE_READ) {
14369566063dSJacob Faibussowitsch PetscCall(PetscFixFilename(infoname, iname));
14379566063dSJacob Faibussowitsch PetscCall(PetscFileRetrieve(PetscObjectComm((PetscObject)viewer), iname, infoname, PETSC_MAX_PATH_LEN, &found));
14389566063dSJacob Faibussowitsch if (found) PetscCall(PetscOptionsInsertFile(PetscObjectComm((PetscObject)viewer), ((PetscObject)viewer)->options, infoname, PETSC_FALSE));
1439dd400576SPatrick Sanan } else if (rank == 0) { /* write or append */
1440cc843e7aSLisandro Dalcin const char *omode = (vbinary->filemode == FILE_MODE_APPEND) ? "a" : "w";
1441cc843e7aSLisandro Dalcin vbinary->fdes_info = fopen(infoname, omode);
1442bf31d2d3SBarry Smith PetscCheck(vbinary->fdes_info, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open .info file %s for writing due to \"%s\"", infoname, strerror(errno));
14435c6c1daeSBarry Smith }
14445c6c1daeSBarry Smith }
14453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
14465c6c1daeSBarry Smith }
14475c6c1daeSBarry Smith
PetscViewerSetUp_Binary(PetscViewer viewer)1448d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerSetUp_Binary(PetscViewer viewer)
1449d71ae5a4SJacob Faibussowitsch {
14505c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary *)viewer->data;
1451cc843e7aSLisandro Dalcin PetscBool usempiio;
1452cc843e7aSLisandro Dalcin
14535c6c1daeSBarry Smith PetscFunctionBegin;
14549566063dSJacob Faibussowitsch if (!vbinary->setfromoptionscalled) PetscCall(PetscViewerSetFromOptions(viewer));
145528b400f6SJacob Faibussowitsch PetscCheck(vbinary->filename, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetName()");
145608401ef6SPierre Jolivet PetscCheck(vbinary->filemode != (PetscFileMode)-1, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode()");
14579566063dSJacob Faibussowitsch PetscCall(PetscViewerFileClose_Binary(viewer));
1458cc843e7aSLisandro Dalcin
14599566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetUseMPIIO(viewer, &usempiio));
1460cc843e7aSLisandro Dalcin if (usempiio) {
1461cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
14629566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetUp_BinaryMPIIO(viewer));
1463cc843e7aSLisandro Dalcin #endif
1464cc843e7aSLisandro Dalcin } else {
14659566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetUp_BinarySTDIO(viewer));
1466cc843e7aSLisandro Dalcin }
14679566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetUp_BinaryInfo(viewer));
1468cc843e7aSLisandro Dalcin
14699566063dSJacob Faibussowitsch PetscCall(PetscLogObjectState((PetscObject)viewer, "File: %s", vbinary->filename));
14703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
14715c6c1daeSBarry Smith }
14725c6c1daeSBarry Smith
PetscViewerView_Binary(PetscViewer v,PetscViewer viewer)1473d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerView_Binary(PetscViewer v, PetscViewer viewer)
1474d71ae5a4SJacob Faibussowitsch {
1475cb6ad94fSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary *)v->data;
1476cb6ad94fSLisandro Dalcin const char *fname = vbinary->filename ? vbinary->filename : "not yet set";
1477cc843e7aSLisandro Dalcin const char *fmode = vbinary->filemode != (PetscFileMode)-1 ? PetscFileModes[vbinary->filemode] : "not yet set";
1478cb6ad94fSLisandro Dalcin PetscBool usempiio;
14792bf49c77SBarry Smith
14802bf49c77SBarry Smith PetscFunctionBegin;
14819566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetUseMPIIO(v, &usempiio));
14829566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Filename: %s\n", fname));
14839566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Mode: %s (%s)\n", fmode, usempiio ? "mpiio" : "stdio"));
14843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
14852bf49c77SBarry Smith }
14862bf49c77SBarry Smith
PetscViewerSetFromOptions_Binary(PetscViewer viewer,PetscOptionItems PetscOptionsObject)1487ce78bad3SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_Binary(PetscViewer viewer, PetscOptionItems PetscOptionsObject)
1488d71ae5a4SJacob Faibussowitsch {
148922a8f86cSLisandro Dalcin PetscViewer_Binary *binary = (PetscViewer_Binary *)viewer->data;
1490d22fd6bcSDave May char defaultname[PETSC_MAX_PATH_LEN];
149103a1d59fSDave May PetscBool flg;
1492e0385b85SDave May
149303a1d59fSDave May PetscFunctionBegin;
14943ba16761SJacob Faibussowitsch if (viewer->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
1495d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Binary PetscViewer Options");
14969566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(defaultname, PETSC_MAX_PATH_LEN - 1, "binaryoutput"));
14979566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-viewer_binary_filename", "Specify filename", "PetscViewerFileSetName", defaultname, defaultname, sizeof(defaultname), &flg));
14989566063dSJacob Faibussowitsch if (flg) PetscCall(PetscViewerFileSetName_Binary(viewer, defaultname));
14999566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_binary_skip_info", "Skip writing/reading .info file", "PetscViewerBinarySetSkipInfo", binary->skipinfo, &binary->skipinfo, NULL));
15009566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_binary_skip_options", "Skip parsing Vec/Mat load options", "PetscViewerBinarySetSkipOptions", binary->skipoptions, &binary->skipoptions, NULL));
15019566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_binary_skip_header", "Skip writing/reading header information", "PetscViewerBinarySetSkipHeader", binary->skipheader, &binary->skipheader, NULL));
150203a1d59fSDave May #if defined(PETSC_HAVE_MPIIO)
15039566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_binary_mpiio", "Use MPI-IO functionality to write/read binary file", "PetscViewerBinarySetUseMPIIO", binary->usempiio, &binary->usempiio, NULL));
15045972f5f3SLisandro Dalcin #else
15057983b6f9SJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_binary_mpiio", "Use MPI-IO functionality to write/read binary file (NOT AVAILABLE)", "PetscViewerBinarySetUseMPIIO", PETSC_FALSE, &flg, NULL));
150603a1d59fSDave May #endif
1507d0609cedSBarry Smith PetscOptionsHeadEnd();
1508bc196f7cSDave May binary->setfromoptionscalled = PETSC_TRUE;
15093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
151003a1d59fSDave May }
151103a1d59fSDave May
15128556b5ebSBarry Smith /*MC
15138556b5ebSBarry Smith PETSCVIEWERBINARY - A viewer that saves to binary files
15148556b5ebSBarry Smith
15151b266c99SBarry Smith Level: beginner
15161b266c99SBarry Smith
1517d1f92df0SBarry Smith .seealso: [](sec_viewers), `PetscViewerBinaryOpen()`, `PETSC_VIEWER_STDOUT_()`, `PETSC_VIEWER_STDOUT_SELF`, `PETSC_VIEWER_STDOUT_WORLD`, `PetscViewerCreate()`, `PetscViewerASCIIOpen()`,
1518811af0c4SBarry Smith `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`, `PETSCVIEWERDRAW`, `PETSCVIEWERSOCKET`
1519811af0c4SBarry Smith `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()`,
1520811af0c4SBarry Smith `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`
15218556b5ebSBarry Smith M*/
15228556b5ebSBarry Smith
PetscViewerCreate_Binary(PetscViewer v)1523d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscViewerCreate_Binary(PetscViewer v)
1524d71ae5a4SJacob Faibussowitsch {
15255c6c1daeSBarry Smith PetscViewer_Binary *vbinary;
15265c6c1daeSBarry Smith
15275c6c1daeSBarry Smith PetscFunctionBegin;
15284dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&vbinary));
15295c6c1daeSBarry Smith v->data = (void *)vbinary;
1530cc843e7aSLisandro Dalcin
153103a1d59fSDave May v->ops->setfromoptions = PetscViewerSetFromOptions_Binary;
15325c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_Binary;
15332bf49c77SBarry Smith v->ops->view = PetscViewerView_Binary;
1534c98fd787SBarry Smith v->ops->setup = PetscViewerSetUp_Binary;
1535cc843e7aSLisandro Dalcin v->ops->flush = NULL; /* Should we support Flush() ? */
1536cc843e7aSLisandro Dalcin v->ops->getsubviewer = PetscViewerGetSubViewer_Binary;
1537cc843e7aSLisandro Dalcin v->ops->restoresubviewer = PetscViewerRestoreSubViewer_Binary;
1538cc843e7aSLisandro Dalcin v->ops->read = PetscViewerBinaryRead;
1539cc843e7aSLisandro Dalcin
1540cc843e7aSLisandro Dalcin vbinary->fdes = -1;
1541e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
1542cc843e7aSLisandro Dalcin vbinary->usempiio = PETSC_FALSE;
1543e8a65b7dSLisandro Dalcin vbinary->mfdes = MPI_FILE_NULL;
1544e8a65b7dSLisandro Dalcin vbinary->mfsub = MPI_FILE_NULL;
1545e8a65b7dSLisandro Dalcin #endif
1546cc843e7aSLisandro Dalcin vbinary->filename = NULL;
15477e4fd573SVaclav Hapla vbinary->filemode = FILE_MODE_UNDEFINED;
154802c9f0b5SLisandro Dalcin vbinary->fdes_info = NULL;
15495c6c1daeSBarry Smith vbinary->skipinfo = PETSC_FALSE;
15505c6c1daeSBarry Smith vbinary->skipoptions = PETSC_TRUE;
15515c6c1daeSBarry Smith vbinary->skipheader = PETSC_FALSE;
15525c6c1daeSBarry Smith vbinary->storecompressed = PETSC_FALSE;
1553f90597f1SStefano Zampini vbinary->ogzfilename = NULL;
15545c6c1daeSBarry Smith vbinary->flowcontrol = 256; /* seems a good number for Cray XT-5 */
15555c6c1daeSBarry Smith
1556cc843e7aSLisandro Dalcin vbinary->setfromoptionscalled = PETSC_FALSE;
1557cc843e7aSLisandro Dalcin
15589566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetFlowControl_C", PetscViewerBinaryGetFlowControl_Binary));
15599566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetFlowControl_C", PetscViewerBinarySetFlowControl_Binary));
15609566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipHeader_C", PetscViewerBinaryGetSkipHeader_Binary));
15619566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipHeader_C", PetscViewerBinarySetSkipHeader_Binary));
15629566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipOptions_C", PetscViewerBinaryGetSkipOptions_Binary));
15639566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipOptions_C", PetscViewerBinarySetSkipOptions_Binary));
15649566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetSkipInfo_C", PetscViewerBinaryGetSkipInfo_Binary));
15659566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetSkipInfo_C", PetscViewerBinarySetSkipInfo_Binary));
15669566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetInfoPointer_C", PetscViewerBinaryGetInfoPointer_Binary));
15679566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_Binary));
15689566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_Binary));
15699566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_Binary));
15709566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_Binary));
15715c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
15729566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinaryGetUseMPIIO_C", PetscViewerBinaryGetUseMPIIO_Binary));
15739566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerBinarySetUseMPIIO_C", PetscViewerBinarySetUseMPIIO_Binary));
15745c6c1daeSBarry Smith #endif
15753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
15765c6c1daeSBarry Smith }
15775c6c1daeSBarry Smith
15785c6c1daeSBarry Smith /*
15795c6c1daeSBarry Smith The variable Petsc_Viewer_Binary_keyval is used to indicate an MPI attribute that
15805c6c1daeSBarry Smith is attached to a communicator, in this case the attribute is a PetscViewer.
15815c6c1daeSBarry Smith */
1582d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_Binary_keyval = MPI_KEYVAL_INVALID;
15835c6c1daeSBarry Smith
15845c6c1daeSBarry Smith /*@C
1585b44f4de4SBarry Smith PETSC_VIEWER_BINARY_ - Creates a `PETSCVIEWERBINARY` `PetscViewer` shared by all processes in a communicator.
15865c6c1daeSBarry Smith
1587d083f849SBarry Smith Collective
15885c6c1daeSBarry Smith
15895c6c1daeSBarry Smith Input Parameter:
1590811af0c4SBarry Smith . comm - the MPI communicator to share the `PETSCVIEWERBINARY`
15915c6c1daeSBarry Smith
15925c6c1daeSBarry Smith Level: intermediate
15935c6c1daeSBarry Smith
15945c6c1daeSBarry Smith Options Database Keys:
159510699b91SBarry Smith + -viewer_binary_filename <name> - filename in which to store the binary data, defaults to binaryoutput
159610699b91SBarry Smith . -viewer_binary_skip_info - true means do not create .info file for this viewer
159710699b91SBarry Smith . -viewer_binary_skip_options - true means do not use the options database for this viewer
159810699b91SBarry Smith . -viewer_binary_skip_header - true means do not store the usual header information in the binary file
159910699b91SBarry Smith - -viewer_binary_mpiio - true means use the file via MPI-IO, maybe faster for large files and many MPI ranks
16005c6c1daeSBarry Smith
1601811af0c4SBarry Smith Environmental variable:
160210699b91SBarry Smith - PETSC_VIEWER_BINARY_FILENAME - filename in which to store the binary data, defaults to binaryoutput
16035c6c1daeSBarry Smith
160434fa283eSBarry Smith Notes:
160534fa283eSBarry Smith This object is destroyed in `PetscFinalize()`, `PetscViewerDestroy()` should never be called on it
160634fa283eSBarry Smith
1607811af0c4SBarry Smith Unlike almost all other PETSc routines, `PETSC_VIEWER_BINARY_` does not return
16085c6c1daeSBarry Smith an error code. The binary PetscViewer is usually used in the form
1609b44f4de4SBarry Smith .vb
1610b44f4de4SBarry Smith XXXView(XXX object, PETSC_VIEWER_BINARY_(comm));
1611b44f4de4SBarry Smith .ve
16125c6c1daeSBarry Smith
1613d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERBINARY`, `PETSC_VIEWER_BINARY_WORLD`, `PETSC_VIEWER_BINARY_SELF`, `PetscViewerBinaryOpen()`, `PetscViewerCreate()`,
1614db781477SPatrick Sanan `PetscViewerDestroy()`
16155c6c1daeSBarry Smith @*/
PETSC_VIEWER_BINARY_(MPI_Comm comm)1616d71ae5a4SJacob Faibussowitsch PetscViewer PETSC_VIEWER_BINARY_(MPI_Comm comm)
1617d71ae5a4SJacob Faibussowitsch {
16185c6c1daeSBarry Smith PetscBool flg;
1619*b8b5be36SMartin Diehl PetscMPIInt iflg;
16205c6c1daeSBarry Smith PetscViewer viewer;
16215c6c1daeSBarry Smith char fname[PETSC_MAX_PATH_LEN];
16225c6c1daeSBarry Smith MPI_Comm ncomm;
16235c6c1daeSBarry Smith
16245c6c1daeSBarry Smith PetscFunctionBegin;
1625648c30bcSBarry Smith PetscCallNull(PetscCommDuplicate(comm, &ncomm, NULL));
1626dd460d27SBarry Smith if (Petsc_Viewer_Binary_keyval == MPI_KEYVAL_INVALID) PetscCallMPINull(MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Viewer_Binary_keyval, NULL));
1627*b8b5be36SMartin Diehl PetscCallMPINull(MPI_Comm_get_attr(ncomm, Petsc_Viewer_Binary_keyval, (void **)&viewer, &iflg));
1628*b8b5be36SMartin Diehl if (!iflg) { /* PetscViewer not yet created */
1629648c30bcSBarry Smith PetscCallNull(PetscOptionsGetenv(ncomm, "PETSC_VIEWER_BINARY_FILENAME", fname, PETSC_MAX_PATH_LEN, &flg));
1630dd460d27SBarry Smith if (!flg) PetscCallNull(PetscStrncpy(fname, "binaryoutput", sizeof(fname)));
1631648c30bcSBarry Smith PetscCallNull(PetscViewerBinaryOpen(ncomm, fname, FILE_MODE_WRITE, &viewer));
1632648c30bcSBarry Smith PetscCallNull(PetscObjectRegisterDestroy((PetscObject)viewer));
1633648c30bcSBarry Smith PetscCallMPINull(MPI_Comm_set_attr(ncomm, Petsc_Viewer_Binary_keyval, (void *)viewer));
16349371c9d4SSatish Balay }
1635648c30bcSBarry Smith PetscCallNull(PetscCommDestroy(&ncomm));
16365c6c1daeSBarry Smith PetscFunctionReturn(viewer);
16375c6c1daeSBarry Smith }
1638