xref: /petsc/src/sys/classes/viewer/impls/binary/binv.c (revision 9dddd24924da2034e9ad37bd0330bf8579e05078)
1af0996ceSBarry Smith #include <petsc/private/viewerimpl.h>    /*I   "petscviewer.h"   I*/
25c6c1daeSBarry Smith 
35c6c1daeSBarry Smith typedef struct  {
45c6c1daeSBarry Smith   int           fdes;                 /* file descriptor, ignored if using MPI IO */
55c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
6bc196f7cSDave May   PetscBool     usempiio;
75c6c1daeSBarry Smith   MPI_File      mfdes;                /* ignored unless using MPI IO */
8e8a65b7dSLisandro Dalcin   MPI_File      mfsub;                /* subviewer support */
95c6c1daeSBarry Smith   MPI_Offset    moff;
105c6c1daeSBarry Smith #endif
11cc843e7aSLisandro Dalcin   char          *filename;            /* file name */
12cc843e7aSLisandro Dalcin   PetscFileMode filemode;             /* read/write/append mode */
135c6c1daeSBarry Smith   FILE          *fdes_info;           /* optional file containing info on binary file*/
145c6c1daeSBarry Smith   PetscBool     storecompressed;      /* gzip the write binary file when closing it*/
15f90597f1SStefano Zampini   char          *ogzfilename;         /* gzip can be run after the filename has been updated */
165c6c1daeSBarry Smith   PetscBool     skipinfo;             /* Don't create info file for writing; don't use for reading */
175c6c1daeSBarry Smith   PetscBool     skipoptions;          /* don't use PETSc options database when loading */
185c6c1daeSBarry Smith   PetscInt      flowcontrol;          /* allow only <flowcontrol> messages outstanding at a time while doing IO */
195c6c1daeSBarry Smith   PetscBool     skipheader;           /* don't write header, only raw data */
20a261c58fSBarry Smith   PetscBool     matlabheaderwritten;  /* if format is PETSC_VIEWER_BINARY_MATLAB has the MATLAB .info header been written yet */
21c98fd787SBarry Smith   PetscBool     setfromoptionscalled;
225c6c1daeSBarry Smith } PetscViewer_Binary;
235c6c1daeSBarry Smith 
24cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
25cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinarySyncMPIIO(PetscViewer viewer)
26cc843e7aSLisandro Dalcin {
27cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data;
28cc843e7aSLisandro Dalcin   PetscErrorCode     ierr;
29cc843e7aSLisandro Dalcin 
30cc843e7aSLisandro Dalcin   PetscFunctionBegin;
31cc843e7aSLisandro Dalcin   if (vbinary->filemode == FILE_MODE_READ) PetscFunctionReturn(0);
32cc843e7aSLisandro Dalcin   if (vbinary->mfsub != MPI_FILE_NULL) {
33ffc4695bSBarry Smith     ierr = MPI_File_sync(vbinary->mfsub);CHKERRMPI(ierr);
34cc843e7aSLisandro Dalcin   }
35cc843e7aSLisandro Dalcin   if (vbinary->mfdes != MPI_FILE_NULL) {
36ffc4695bSBarry Smith     ierr = MPI_Barrier(PetscObjectComm((PetscObject)viewer));CHKERRMPI(ierr);
37ffc4695bSBarry Smith     ierr = MPI_File_sync(vbinary->mfdes);CHKERRMPI(ierr);
38cc843e7aSLisandro Dalcin   }
39cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
40cc843e7aSLisandro Dalcin }
41cc843e7aSLisandro Dalcin #endif
42cc843e7aSLisandro Dalcin 
4381f0254dSBarry Smith static PetscErrorCode PetscViewerGetSubViewer_Binary(PetscViewer viewer,MPI_Comm comm,PetscViewer *outviewer)
445c6c1daeSBarry Smith {
45e8a65b7dSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data;
46cc843e7aSLisandro Dalcin   PetscMPIInt        rank;
47cc843e7aSLisandro Dalcin   PetscErrorCode     ierr;
485c6c1daeSBarry Smith 
495c6c1daeSBarry Smith   PetscFunctionBegin;
50c98fd787SBarry Smith   ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr);
51e8a65b7dSLisandro Dalcin 
52e8a65b7dSLisandro Dalcin   /* Return subviewer in process zero */
53ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank);CHKERRMPI(ierr);
545c6c1daeSBarry Smith   if (!rank) {
55e5afcf28SBarry Smith     PetscMPIInt flg;
56e5afcf28SBarry Smith 
57ffc4695bSBarry Smith     ierr = MPI_Comm_compare(PETSC_COMM_SELF,comm,&flg);CHKERRMPI(ierr);
58e5afcf28SBarry Smith     if (flg != MPI_IDENT && flg != MPI_CONGRUENT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscViewerGetSubViewer() for PETSCVIEWERBINARY requires a singleton MPI_Comm");
59e5afcf28SBarry Smith     ierr = PetscViewerCreate(comm,outviewer);CHKERRQ(ierr);
605c6c1daeSBarry Smith     ierr = PetscViewerSetType(*outviewer,PETSCVIEWERBINARY);CHKERRQ(ierr);
61e8a65b7dSLisandro Dalcin     ierr = PetscMemcpy((*outviewer)->data,vbinary,sizeof(PetscViewer_Binary));CHKERRQ(ierr);
6212f4c3a9SLisandro Dalcin     (*outviewer)->setupcalled = PETSC_TRUE;
6396509d17SLisandro Dalcin   } else {
6496509d17SLisandro Dalcin     *outviewer = NULL;
6596509d17SLisandro Dalcin   }
66e8a65b7dSLisandro Dalcin 
67e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
68e8a65b7dSLisandro Dalcin   if (vbinary->usempiio && *outviewer) {
69e8a65b7dSLisandro Dalcin     PetscViewer_Binary *obinary = (PetscViewer_Binary*)(*outviewer)->data;
70e8a65b7dSLisandro Dalcin     /* Parent viewer opens a new MPI file handle on PETSC_COMM_SELF and keeps track of it for future reuse */
71e8a65b7dSLisandro Dalcin     if (vbinary->mfsub == MPI_FILE_NULL) {
72e8a65b7dSLisandro Dalcin       int amode;
73cc843e7aSLisandro Dalcin       switch (vbinary->filemode) {
74e8a65b7dSLisandro Dalcin       case FILE_MODE_READ:   amode = MPI_MODE_RDONLY; break;
75e8a65b7dSLisandro Dalcin       case FILE_MODE_WRITE:  amode = MPI_MODE_WRONLY; break;
7622a8f86cSLisandro Dalcin       case FILE_MODE_APPEND: amode = MPI_MODE_WRONLY; break;
77cc843e7aSLisandro Dalcin       default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported file mode %s",PetscFileModes[vbinary->filemode]);
78e8a65b7dSLisandro Dalcin       }
79ffc4695bSBarry Smith       ierr = MPI_File_open(PETSC_COMM_SELF,vbinary->filename,amode,MPI_INFO_NULL,&vbinary->mfsub);CHKERRMPI(ierr);
80e8a65b7dSLisandro Dalcin     }
81e8a65b7dSLisandro Dalcin     /* Subviewer gets the MPI file handle on PETSC_COMM_SELF */
82e8a65b7dSLisandro Dalcin     obinary->mfdes = vbinary->mfsub;
83e8a65b7dSLisandro Dalcin     obinary->mfsub = MPI_FILE_NULL;
84e8a65b7dSLisandro Dalcin     obinary->moff  = vbinary->moff;
85e8a65b7dSLisandro Dalcin   }
86e8a65b7dSLisandro Dalcin #endif
87cc843e7aSLisandro Dalcin 
88cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
89cc843e7aSLisandro Dalcin   ierr = PetscViewerBinarySyncMPIIO(viewer);CHKERRQ(ierr);
90cc843e7aSLisandro Dalcin #endif
915c6c1daeSBarry Smith   PetscFunctionReturn(0);
925c6c1daeSBarry Smith }
935c6c1daeSBarry Smith 
9481f0254dSBarry Smith static PetscErrorCode PetscViewerRestoreSubViewer_Binary(PetscViewer viewer,MPI_Comm comm,PetscViewer *outviewer)
955c6c1daeSBarry Smith {
96e8a65b7dSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data;
97cc843e7aSLisandro Dalcin   PetscMPIInt        rank;
98e8a65b7dSLisandro Dalcin   PetscErrorCode     ierr;
99e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
100e8a65b7dSLisandro Dalcin   MPI_Offset         moff = 0;
101e8a65b7dSLisandro Dalcin #endif
1025c6c1daeSBarry Smith 
1035c6c1daeSBarry Smith   PetscFunctionBegin;
104ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank);CHKERRMPI(ierr);
105e8a65b7dSLisandro Dalcin   if (rank && *outviewer) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Subviewer not obtained from viewer");
106e8a65b7dSLisandro Dalcin 
107e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
108e8a65b7dSLisandro Dalcin   if (vbinary->usempiio && *outviewer) {
109e8a65b7dSLisandro Dalcin     PetscViewer_Binary *obinary = (PetscViewer_Binary*)(*outviewer)->data;
110e8a65b7dSLisandro Dalcin     if (obinary->mfdes != vbinary->mfsub) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Subviewer not obtained from viewer");
111ffc4695bSBarry Smith     if (obinary->mfsub != MPI_FILE_NULL) {ierr = MPI_File_close(&obinary->mfsub);CHKERRMPI(ierr);}
112e8a65b7dSLisandro Dalcin     moff = obinary->moff;
113e8a65b7dSLisandro Dalcin   }
114e8a65b7dSLisandro Dalcin #endif
115e8a65b7dSLisandro Dalcin 
116e8a65b7dSLisandro Dalcin   if (*outviewer) {
117e8a65b7dSLisandro Dalcin     PetscViewer_Binary *obinary = (PetscViewer_Binary*)(*outviewer)->data;
118e8a65b7dSLisandro Dalcin     if (obinary->fdes != vbinary->fdes) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Subviewer not obtained from viewer");
1195c6c1daeSBarry Smith     ierr = PetscFree((*outviewer)->data);CHKERRQ(ierr);
1205c6c1daeSBarry Smith     ierr = PetscHeaderDestroy(outviewer);CHKERRQ(ierr);
1215c6c1daeSBarry Smith   }
122e8a65b7dSLisandro Dalcin 
123e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
124e8a65b7dSLisandro Dalcin   if (vbinary->usempiio) {
125e8a65b7dSLisandro Dalcin     PetscInt64 ioff = (PetscInt64)moff; /* We could use MPI_OFFSET datatype (requires MPI 2.2) */
126ffc4695bSBarry Smith     ierr = MPI_Bcast(&ioff,1,MPIU_INT64,0,PetscObjectComm((PetscObject)viewer));CHKERRMPI(ierr);
127e8a65b7dSLisandro Dalcin     vbinary->moff = (MPI_Offset)ioff;
128e8a65b7dSLisandro Dalcin   }
129e8a65b7dSLisandro Dalcin #endif
130cc843e7aSLisandro Dalcin 
131cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
132cc843e7aSLisandro Dalcin   ierr = PetscViewerBinarySyncMPIIO(viewer);CHKERRQ(ierr);
133cc843e7aSLisandro Dalcin #endif
1345c6c1daeSBarry Smith   PetscFunctionReturn(0);
1355c6c1daeSBarry Smith }
1365c6c1daeSBarry Smith 
1375c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
1385c6c1daeSBarry Smith /*@C
13922a8f86cSLisandro Dalcin     PetscViewerBinaryGetMPIIOOffset - Gets the current global offset that should be passed to MPI_File_set_view() or MPI_File_{write|read}_at[_all]()
1405c6c1daeSBarry Smith 
1415c6c1daeSBarry Smith     Not Collective
1425c6c1daeSBarry Smith 
1435c6c1daeSBarry Smith     Input Parameter:
1445c6c1daeSBarry Smith .   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
1455c6c1daeSBarry Smith 
1465c6c1daeSBarry Smith     Output Parameter:
14722a8f86cSLisandro Dalcin .   off - the current global offset
1485c6c1daeSBarry Smith 
1495c6c1daeSBarry Smith     Level: advanced
1505c6c1daeSBarry Smith 
1515c6c1daeSBarry Smith     Fortran Note:
1525c6c1daeSBarry Smith     This routine is not supported in Fortran.
1535c6c1daeSBarry Smith 
1545c6c1daeSBarry Smith     Use PetscViewerBinaryAddMPIIOOffset() to increase this value after you have written a view.
1555c6c1daeSBarry Smith 
1565c6c1daeSBarry Smith 
15767918a83SBarry Smith .seealso: PetscViewerBinaryOpen(), PetscViewerBinaryGetInfoPointer(), PetscViewerBinaryGetUseMPIIO(), PetscViewerBinarySetUseMPIIO(), PetscViewerBinaryAddMPIIOOffset()
1585c6c1daeSBarry Smith @*/
1595c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryGetMPIIOOffset(PetscViewer viewer,MPI_Offset *off)
1605c6c1daeSBarry Smith {
16122a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
1625c6c1daeSBarry Smith 
1635c6c1daeSBarry Smith   PetscFunctionBegin;
16422a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer,PETSC_VIEWER_CLASSID,1,PETSCVIEWERBINARY);
16522a8f86cSLisandro Dalcin   PetscValidPointer(off,2);
16622a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary*)viewer->data;
1675c6c1daeSBarry Smith   *off = vbinary->moff;
1685c6c1daeSBarry Smith   PetscFunctionReturn(0);
1695c6c1daeSBarry Smith }
1705c6c1daeSBarry Smith 
1715c6c1daeSBarry Smith /*@C
17222a8f86cSLisandro Dalcin     PetscViewerBinaryAddMPIIOOffset - Adds to the current global offset
1735c6c1daeSBarry Smith 
17422a8f86cSLisandro Dalcin     Logically Collective
1755c6c1daeSBarry Smith 
1765c6c1daeSBarry Smith     Input Parameters:
1775c6c1daeSBarry Smith +   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
17822a8f86cSLisandro Dalcin -   off - the addition to the global offset
1795c6c1daeSBarry Smith 
1805c6c1daeSBarry Smith     Level: advanced
1815c6c1daeSBarry Smith 
1825c6c1daeSBarry Smith     Fortran Note:
1835c6c1daeSBarry Smith     This routine is not supported in Fortran.
1845c6c1daeSBarry Smith 
18522a8f86cSLisandro Dalcin     Use PetscViewerBinaryGetMPIIOOffset() to get the value that you should pass to MPI_File_set_view() or MPI_File_{write|read}_at[_all]()
1865c6c1daeSBarry Smith 
1875c6c1daeSBarry Smith 
18867918a83SBarry Smith .seealso: PetscViewerBinaryOpen(),PetscViewerBinaryGetInfoPointer(), PetscViewerBinaryGetUseMPIIO(), PetscViewerBinarySetUseMPIIO(), PetscViewerBinaryGetMPIIOOffset()
1895c6c1daeSBarry Smith @*/
1905c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryAddMPIIOOffset(PetscViewer viewer,MPI_Offset off)
1915c6c1daeSBarry Smith {
19222a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
1935c6c1daeSBarry Smith 
1945c6c1daeSBarry Smith   PetscFunctionBegin;
19522a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer,PETSC_VIEWER_CLASSID,1,PETSCVIEWERBINARY);
19622a8f86cSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer,(PetscInt)off,2);
19722a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary*)viewer->data;
1985c6c1daeSBarry Smith   vbinary->moff += off;
1995c6c1daeSBarry Smith   PetscFunctionReturn(0);
2005c6c1daeSBarry Smith }
2015c6c1daeSBarry Smith 
2025c6c1daeSBarry Smith /*@C
2035c6c1daeSBarry Smith     PetscViewerBinaryGetMPIIODescriptor - Extracts the MPI IO file descriptor from a PetscViewer.
2045c6c1daeSBarry Smith 
2055c6c1daeSBarry Smith     Not Collective
2065c6c1daeSBarry Smith 
2075c6c1daeSBarry Smith     Input Parameter:
2085c6c1daeSBarry Smith .   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
2095c6c1daeSBarry Smith 
2105c6c1daeSBarry Smith     Output Parameter:
2115c6c1daeSBarry Smith .   fdes - file descriptor
2125c6c1daeSBarry Smith 
2135c6c1daeSBarry Smith     Level: advanced
2145c6c1daeSBarry Smith 
2155c6c1daeSBarry Smith     Fortran Note:
2165c6c1daeSBarry Smith     This routine is not supported in Fortran.
2175c6c1daeSBarry Smith 
2185c6c1daeSBarry Smith 
21967918a83SBarry Smith .seealso: PetscViewerBinaryOpen(),PetscViewerBinaryGetInfoPointer(), PetscViewerBinaryGetUseMPIIO(), PetscViewerBinarySetUseMPIIO(), PetscViewerBinaryGetMPIIOOffset()
2205c6c1daeSBarry Smith @*/
2215c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryGetMPIIODescriptor(PetscViewer viewer,MPI_File *fdes)
2225c6c1daeSBarry Smith {
22303a1d59fSDave May   PetscErrorCode     ierr;
22422a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
2255c6c1daeSBarry Smith 
2265c6c1daeSBarry Smith   PetscFunctionBegin;
22722a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer,PETSC_VIEWER_CLASSID,1,PETSCVIEWERBINARY);
22822a8f86cSLisandro Dalcin   PetscValidPointer(fdes,2);
229c98fd787SBarry Smith   ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr);
23022a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary*)viewer->data;
2315c6c1daeSBarry Smith   *fdes = vbinary->mfdes;
2325c6c1daeSBarry Smith   PetscFunctionReturn(0);
2335c6c1daeSBarry Smith }
234cc843e7aSLisandro Dalcin #endif
2355c6c1daeSBarry Smith 
236cc843e7aSLisandro Dalcin /*@
237cc843e7aSLisandro Dalcin     PetscViewerBinarySetUseMPIIO - Sets a binary viewer to use MPI-IO for reading/writing. Must be called
238cc843e7aSLisandro Dalcin         before PetscViewerFileSetName()
239cc843e7aSLisandro Dalcin 
240cc843e7aSLisandro Dalcin     Logically Collective on PetscViewer
241cc843e7aSLisandro Dalcin 
242cc843e7aSLisandro Dalcin     Input Parameters:
243cc843e7aSLisandro Dalcin +   viewer - the PetscViewer; must be a binary
244cc843e7aSLisandro Dalcin -   use - PETSC_TRUE means MPI-IO will be used
245cc843e7aSLisandro Dalcin 
246cc843e7aSLisandro Dalcin     Options Database:
247cc843e7aSLisandro Dalcin     -viewer_binary_mpiio : Flag for using MPI-IO
248cc843e7aSLisandro Dalcin 
249cc843e7aSLisandro Dalcin     Level: advanced
250cc843e7aSLisandro Dalcin 
251cc843e7aSLisandro Dalcin .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
252cc843e7aSLisandro Dalcin           PetscViewerBinaryGetUseMPIIO()
253cc843e7aSLisandro Dalcin 
254cc843e7aSLisandro Dalcin @*/
255cc843e7aSLisandro Dalcin PetscErrorCode PetscViewerBinarySetUseMPIIO(PetscViewer viewer,PetscBool use)
256a8aae444SDave May {
257cc843e7aSLisandro Dalcin   PetscErrorCode ierr;
258a8aae444SDave May 
259a8aae444SDave May   PetscFunctionBegin;
260cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
261cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer,use,2);
262cc843e7aSLisandro Dalcin   ierr = PetscTryMethod(viewer,"PetscViewerBinarySetUseMPIIO_C",(PetscViewer,PetscBool),(viewer,use));CHKERRQ(ierr);
263cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
264cc843e7aSLisandro Dalcin }
265cc843e7aSLisandro Dalcin 
266cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
267cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinarySetUseMPIIO_Binary(PetscViewer viewer,PetscBool use)
268cc843e7aSLisandro Dalcin {
269cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data;
270cc843e7aSLisandro Dalcin   PetscFunctionBegin;
271cc843e7aSLisandro Dalcin   if (viewer->setupcalled && vbinary->usempiio != use) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER,"Cannot change MPIIO to %s after setup",PetscBools[use]);
272cc843e7aSLisandro Dalcin   vbinary->usempiio = use;
273a8aae444SDave May   PetscFunctionReturn(0);
274a8aae444SDave May }
275a8aae444SDave May #endif
276a8aae444SDave May 
277cc843e7aSLisandro Dalcin /*@
278bc196f7cSDave May     PetscViewerBinaryGetUseMPIIO - Returns PETSC_TRUE if the binary viewer uses MPI-IO.
2795c6c1daeSBarry Smith 
2805c6c1daeSBarry Smith     Not Collective
2815c6c1daeSBarry Smith 
2825c6c1daeSBarry Smith     Input Parameter:
2835c6c1daeSBarry Smith .   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
2845c6c1daeSBarry Smith 
2855c6c1daeSBarry Smith     Output Parameter:
286cc843e7aSLisandro Dalcin -   use - PETSC_TRUE if MPI-IO is being used
2875c6c1daeSBarry Smith 
2885c6c1daeSBarry Smith     Options Database:
2895c6c1daeSBarry Smith     -viewer_binary_mpiio : Flag for using MPI-IO
2905c6c1daeSBarry Smith 
2915c6c1daeSBarry Smith     Level: advanced
2925c6c1daeSBarry Smith 
293bc196f7cSDave May     Note:
294bc196f7cSDave May     If MPI-IO is not available, this function will always return PETSC_FALSE
295bc196f7cSDave May 
2965c6c1daeSBarry Smith     Fortran Note:
2975c6c1daeSBarry Smith     This routine is not supported in Fortran.
2985c6c1daeSBarry Smith 
2995c6c1daeSBarry Smith 
30067918a83SBarry Smith .seealso: PetscViewerBinaryOpen(), PetscViewerBinaryGetInfoPointer(), PetscViewerBinarySetUseMPIIO(), PetscViewerBinaryGetMPIIOOffset()
3015c6c1daeSBarry Smith @*/
302cc843e7aSLisandro Dalcin PetscErrorCode PetscViewerBinaryGetUseMPIIO(PetscViewer viewer,PetscBool *use)
3035c6c1daeSBarry Smith {
304a8aae444SDave May   PetscErrorCode ierr;
3055c6c1daeSBarry Smith 
3065c6c1daeSBarry Smith   PetscFunctionBegin;
307cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
308cc843e7aSLisandro Dalcin   PetscValidBoolPointer(use,2);
309cc843e7aSLisandro Dalcin   *use = PETSC_FALSE;
310cc843e7aSLisandro Dalcin   ierr = PetscTryMethod(viewer,"PetscViewerBinaryGetUseMPIIO_C",(PetscViewer,PetscBool*),(viewer,use));CHKERRQ(ierr);
3115c6c1daeSBarry Smith   PetscFunctionReturn(0);
3125c6c1daeSBarry Smith }
3135c6c1daeSBarry Smith 
314cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
315cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinaryGetUseMPIIO_Binary(PetscViewer viewer,PetscBool  *use)
3165c6c1daeSBarry Smith {
3175c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data;
3185c6c1daeSBarry Smith 
3195c6c1daeSBarry Smith   PetscFunctionBegin;
320cc843e7aSLisandro Dalcin   *use = vbinary->usempiio;
321cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
322cc843e7aSLisandro Dalcin }
323cc843e7aSLisandro Dalcin #endif
324cc843e7aSLisandro Dalcin 
325cc843e7aSLisandro Dalcin /*@
326cc843e7aSLisandro Dalcin     PetscViewerBinarySetFlowControl - Sets how many messages are allowed to outstanding at the same time during parallel IO reads/writes
327cc843e7aSLisandro Dalcin 
328cc843e7aSLisandro Dalcin     Not Collective
329cc843e7aSLisandro Dalcin 
330cc843e7aSLisandro Dalcin     Input Parameter:
331cc843e7aSLisandro Dalcin +   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
332cc843e7aSLisandro Dalcin -   fc - the number of messages, defaults to 256 if this function was not called
333cc843e7aSLisandro Dalcin 
334cc843e7aSLisandro Dalcin     Level: advanced
335cc843e7aSLisandro Dalcin 
336cc843e7aSLisandro Dalcin .seealso: PetscViewerBinaryOpen(),PetscViewerBinaryGetInfoPointer(), PetscViewerBinaryGetFlowControl()
337cc843e7aSLisandro Dalcin 
338cc843e7aSLisandro Dalcin @*/
339cc843e7aSLisandro Dalcin PetscErrorCode  PetscViewerBinarySetFlowControl(PetscViewer viewer,PetscInt fc)
340cc843e7aSLisandro Dalcin {
341cc843e7aSLisandro Dalcin   PetscErrorCode ierr;
342cc843e7aSLisandro Dalcin 
343cc843e7aSLisandro Dalcin   PetscFunctionBegin;
344cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
345cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer,fc,2);
346cc843e7aSLisandro Dalcin   ierr = PetscTryMethod(viewer,"PetscViewerBinarySetFlowControl_C",(PetscViewer,PetscInt),(viewer,fc));CHKERRQ(ierr);
3475c6c1daeSBarry Smith   PetscFunctionReturn(0);
3485c6c1daeSBarry Smith }
3495c6c1daeSBarry Smith 
350cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinarySetFlowControl_Binary(PetscViewer viewer,PetscInt fc)
351cc843e7aSLisandro Dalcin {
352cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data;
353cc843e7aSLisandro Dalcin 
354cc843e7aSLisandro Dalcin   PetscFunctionBegin;
355cc843e7aSLisandro Dalcin   if (fc <= 1) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_OUTOFRANGE,"Flow control count must be greater than 1, %D was set",fc);
356cc843e7aSLisandro Dalcin   vbinary->flowcontrol = fc;
357cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
358cc843e7aSLisandro Dalcin }
359cc843e7aSLisandro Dalcin 
360cc843e7aSLisandro Dalcin /*@
3615c6c1daeSBarry Smith     PetscViewerBinaryGetFlowControl - Returns how many messages are allowed to outstanding at the same time during parallel IO reads/writes
3625c6c1daeSBarry Smith 
3635c6c1daeSBarry Smith     Not Collective
3645c6c1daeSBarry Smith 
3655c6c1daeSBarry Smith     Input Parameter:
3665c6c1daeSBarry Smith .   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
3675c6c1daeSBarry Smith 
3685c6c1daeSBarry Smith     Output Parameter:
3695c6c1daeSBarry Smith .   fc - the number of messages
3705c6c1daeSBarry Smith 
3715c6c1daeSBarry Smith     Level: advanced
3725c6c1daeSBarry Smith 
3735c6c1daeSBarry Smith .seealso: PetscViewerBinaryOpen(),PetscViewerBinaryGetInfoPointer(), PetscViewerBinarySetFlowControl()
3745c6c1daeSBarry Smith 
3755c6c1daeSBarry Smith @*/
3765c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryGetFlowControl(PetscViewer viewer,PetscInt *fc)
3775c6c1daeSBarry Smith {
3785c6c1daeSBarry Smith   PetscErrorCode ierr;
3795c6c1daeSBarry Smith 
3805c6c1daeSBarry Smith   PetscFunctionBegin;
381cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
382cc843e7aSLisandro Dalcin   PetscValidIntPointer(fc,2);
383163d334eSBarry Smith   ierr = PetscUseMethod(viewer,"PetscViewerBinaryGetFlowControl_C",(PetscViewer,PetscInt*),(viewer,fc));CHKERRQ(ierr);
3845c6c1daeSBarry Smith   PetscFunctionReturn(0);
3855c6c1daeSBarry Smith }
3865c6c1daeSBarry Smith 
387cc843e7aSLisandro Dalcin static PetscErrorCode  PetscViewerBinaryGetFlowControl_Binary(PetscViewer viewer,PetscInt *fc)
3885c6c1daeSBarry Smith {
3895c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data;
3905c6c1daeSBarry Smith 
3915c6c1daeSBarry Smith   PetscFunctionBegin;
392cc843e7aSLisandro Dalcin   *fc = vbinary->flowcontrol;
3935c6c1daeSBarry Smith   PetscFunctionReturn(0);
3945c6c1daeSBarry Smith }
3955c6c1daeSBarry Smith 
3965c6c1daeSBarry Smith /*@C
3975c6c1daeSBarry Smith     PetscViewerBinaryGetDescriptor - Extracts the file descriptor from a PetscViewer.
3985c6c1daeSBarry Smith 
3995872f025SBarry Smith     Collective On PetscViewer
4005c6c1daeSBarry Smith 
4015c6c1daeSBarry Smith     Input Parameter:
4025c6c1daeSBarry Smith .   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
4035c6c1daeSBarry Smith 
4045c6c1daeSBarry Smith     Output Parameter:
4055c6c1daeSBarry Smith .   fdes - file descriptor
4065c6c1daeSBarry Smith 
4075c6c1daeSBarry Smith     Level: advanced
4085c6c1daeSBarry Smith 
4095c6c1daeSBarry Smith     Notes:
4105c6c1daeSBarry Smith       For writable binary PetscViewers, the descriptor will only be valid for the
4115c6c1daeSBarry Smith     first processor in the communicator that shares the PetscViewer. For readable
4125c6c1daeSBarry Smith     files it will only be valid on nodes that have the file. If node 0 does not
4135c6c1daeSBarry Smith     have the file it generates an error even if another node does have the file.
4145c6c1daeSBarry Smith 
4155c6c1daeSBarry Smith     Fortran Note:
4165c6c1daeSBarry Smith     This routine is not supported in Fortran.
4175c6c1daeSBarry Smith 
41895452b02SPatrick Sanan     Developer Notes:
41995452b02SPatrick Sanan     This must be called on all processes because Dave May changed
4205872f025SBarry Smith     the source code that this may be trigger a PetscViewerSetUp() call if it was not previously triggered.
4215872f025SBarry Smith 
4225872f025SBarry Smith 
4235c6c1daeSBarry Smith 
4245c6c1daeSBarry Smith .seealso: PetscViewerBinaryOpen(),PetscViewerBinaryGetInfoPointer()
4255c6c1daeSBarry Smith @*/
4265c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryGetDescriptor(PetscViewer viewer,int *fdes)
4275c6c1daeSBarry Smith {
42803a1d59fSDave May   PetscErrorCode     ierr;
42922a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
4305c6c1daeSBarry Smith 
4315c6c1daeSBarry Smith   PetscFunctionBegin;
43222a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer,PETSC_VIEWER_CLASSID,1,PETSCVIEWERBINARY);
43322a8f86cSLisandro Dalcin   PetscValidPointer(fdes,2);
434c98fd787SBarry Smith   ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr);
43522a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary*)viewer->data;
4365c6c1daeSBarry Smith   *fdes = vbinary->fdes;
4375c6c1daeSBarry Smith   PetscFunctionReturn(0);
4385c6c1daeSBarry Smith }
4395c6c1daeSBarry Smith 
4405c6c1daeSBarry Smith /*@
4415c6c1daeSBarry Smith     PetscViewerBinarySkipInfo - Binary file will not have .info file created with it
4425c6c1daeSBarry Smith 
4435c6c1daeSBarry Smith     Not Collective
4445c6c1daeSBarry Smith 
445fd292e60Sprj-     Input Parameter:
4465c6c1daeSBarry Smith .   viewer - PetscViewer context, obtained from PetscViewerCreate()
4475c6c1daeSBarry Smith 
4485c6c1daeSBarry Smith     Options Database Key:
4495c6c1daeSBarry Smith .   -viewer_binary_skip_info
4505c6c1daeSBarry Smith 
4515c6c1daeSBarry Smith     Level: advanced
4525c6c1daeSBarry Smith 
45395452b02SPatrick Sanan     Notes:
45495452b02SPatrick Sanan     This must be called after PetscViewerSetType(). If you use PetscViewerBinaryOpen() then
4555c6c1daeSBarry Smith     you can only skip the info file with the -viewer_binary_skip_info flag. To use the function you must open the
456a2d7db39SDave May     viewer with PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinarySkipInfo().
4575c6c1daeSBarry Smith 
4585c6c1daeSBarry Smith     The .info contains meta information about the data in the binary file, for example the block size if it was
4595c6c1daeSBarry Smith     set for a vector or matrix.
4605c6c1daeSBarry Smith 
4615c6c1daeSBarry Smith .seealso: PetscViewerBinaryOpen(), PetscViewerBinaryGetDescriptor(), PetscViewerBinarySetSkipOptions(),
462807ea322SDave May           PetscViewerBinaryGetSkipOptions(), PetscViewerBinaryGetSkipInfo()
4635c6c1daeSBarry Smith @*/
4645c6c1daeSBarry Smith PetscErrorCode PetscViewerBinarySkipInfo(PetscViewer viewer)
4655c6c1daeSBarry Smith {
466cc843e7aSLisandro Dalcin   PetscErrorCode ierr;
4675c6c1daeSBarry Smith 
4685c6c1daeSBarry Smith   PetscFunctionBegin;
469cc843e7aSLisandro Dalcin   ierr = PetscViewerBinarySetSkipInfo(viewer,PETSC_TRUE);CHKERRQ(ierr);
470807ea322SDave May   PetscFunctionReturn(0);
471807ea322SDave May }
472807ea322SDave May 
473807ea322SDave May /*@
474807ea322SDave May     PetscViewerBinarySetSkipInfo - Binary file will not have .info file created with it
475807ea322SDave May 
476807ea322SDave May     Not Collective
477807ea322SDave May 
478fd292e60Sprj-     Input Parameter:
479cc843e7aSLisandro Dalcin +   viewer - PetscViewer context, obtained from PetscViewerCreate()
480cc843e7aSLisandro Dalcin -   skip - PETSC_TRUE implies the .info file will not be generated
481807ea322SDave May 
482807ea322SDave May     Options Database Key:
483807ea322SDave May .   -viewer_binary_skip_info
484807ea322SDave May 
485807ea322SDave May     Level: advanced
486807ea322SDave May 
487807ea322SDave May .seealso: PetscViewerBinaryOpen(), PetscViewerBinaryGetDescriptor(), PetscViewerBinarySetSkipOptions(),
488807ea322SDave May           PetscViewerBinaryGetSkipOptions(), PetscViewerBinaryGetSkipInfo()
489807ea322SDave May @*/
490807ea322SDave May PetscErrorCode PetscViewerBinarySetSkipInfo(PetscViewer viewer,PetscBool skip)
491807ea322SDave May {
492807ea322SDave May   PetscErrorCode ierr;
493807ea322SDave May 
494807ea322SDave May   PetscFunctionBegin;
495cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
496cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer,skip,2);
497cc843e7aSLisandro Dalcin   ierr = PetscTryMethod(viewer,"PetscViewerBinarySetSkipInfo_C",(PetscViewer,PetscBool),(viewer,skip));CHKERRQ(ierr);
498807ea322SDave May   PetscFunctionReturn(0);
499807ea322SDave May }
500807ea322SDave May 
501cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinarySetSkipInfo_Binary(PetscViewer viewer,PetscBool skip)
502807ea322SDave May {
503807ea322SDave May   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data;
504807ea322SDave May 
505807ea322SDave May   PetscFunctionBegin;
506cc843e7aSLisandro Dalcin   vbinary->skipinfo = skip;
507807ea322SDave May   PetscFunctionReturn(0);
508807ea322SDave May }
509807ea322SDave May 
510807ea322SDave May /*@
511807ea322SDave May     PetscViewerBinaryGetSkipInfo - check if viewer wrote a .info file
512807ea322SDave May 
513807ea322SDave May     Not Collective
514807ea322SDave May 
515807ea322SDave May     Input Parameter:
516807ea322SDave May .   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
517807ea322SDave May 
518807ea322SDave May     Output Parameter:
519807ea322SDave May .   skip - PETSC_TRUE implies the .info file was not generated
520807ea322SDave May 
521807ea322SDave May     Level: advanced
522807ea322SDave May 
52395452b02SPatrick Sanan     Notes:
52495452b02SPatrick Sanan     This must be called after PetscViewerSetType()
525807ea322SDave May 
526807ea322SDave May .seealso: PetscViewerBinaryOpen(), PetscViewerBinaryGetDescriptor(), PetscViewerBinarySkipInfo(),
527807ea322SDave May           PetscViewerBinarySetSkipOptions(), PetscViewerBinarySetSkipInfo()
528807ea322SDave May @*/
529807ea322SDave May PetscErrorCode PetscViewerBinaryGetSkipInfo(PetscViewer viewer,PetscBool *skip)
530807ea322SDave May {
531807ea322SDave May   PetscErrorCode ierr;
532807ea322SDave May 
533807ea322SDave May   PetscFunctionBegin;
534cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
535cc843e7aSLisandro Dalcin   PetscValidBoolPointer(skip,2);
536807ea322SDave May   ierr = PetscUseMethod(viewer,"PetscViewerBinaryGetSkipInfo_C",(PetscViewer,PetscBool*),(viewer,skip));CHKERRQ(ierr);
537807ea322SDave May   PetscFunctionReturn(0);
538807ea322SDave May }
539807ea322SDave May 
540cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinaryGetSkipInfo_Binary(PetscViewer viewer,PetscBool *skip)
541807ea322SDave May {
542807ea322SDave May   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data;
543807ea322SDave May 
544807ea322SDave May   PetscFunctionBegin;
545cc843e7aSLisandro Dalcin   *skip  = vbinary->skipinfo;
546807ea322SDave May   PetscFunctionReturn(0);
547807ea322SDave May }
548807ea322SDave May 
5495c6c1daeSBarry Smith /*@
5505c6c1daeSBarry Smith     PetscViewerBinarySetSkipOptions - do not use the PETSc options database when loading objects
5515c6c1daeSBarry Smith 
5525c6c1daeSBarry Smith     Not Collective
5535c6c1daeSBarry Smith 
5545c6c1daeSBarry Smith     Input Parameters:
5555c6c1daeSBarry Smith +   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
5565c6c1daeSBarry Smith -   skip - PETSC_TRUE means do not use
5575c6c1daeSBarry Smith 
5585c6c1daeSBarry Smith     Options Database Key:
5595c6c1daeSBarry Smith .   -viewer_binary_skip_options
5605c6c1daeSBarry Smith 
5615c6c1daeSBarry Smith     Level: advanced
5625c6c1daeSBarry Smith 
56395452b02SPatrick Sanan     Notes:
56495452b02SPatrick Sanan     This must be called after PetscViewerSetType()
5655c6c1daeSBarry Smith 
5665c6c1daeSBarry Smith .seealso: PetscViewerBinaryOpen(), PetscViewerBinaryGetDescriptor(), PetscViewerBinarySkipInfo(),
5675c6c1daeSBarry Smith           PetscViewerBinaryGetSkipOptions()
5685c6c1daeSBarry Smith @*/
5695c6c1daeSBarry Smith PetscErrorCode PetscViewerBinarySetSkipOptions(PetscViewer viewer,PetscBool skip)
5705c6c1daeSBarry Smith {
571807ea322SDave May   PetscErrorCode ierr;
5725c6c1daeSBarry Smith 
5735c6c1daeSBarry Smith   PetscFunctionBegin;
574cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
575cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer,skip,2);
576cc843e7aSLisandro Dalcin   ierr = PetscTryMethod(viewer,"PetscViewerBinarySetSkipOptions_C",(PetscViewer,PetscBool),(viewer,skip));CHKERRQ(ierr);
577807ea322SDave May   PetscFunctionReturn(0);
578807ea322SDave May }
579807ea322SDave May 
580cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinarySetSkipOptions_Binary(PetscViewer viewer,PetscBool skip)
581807ea322SDave May {
582cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data;
583807ea322SDave May 
584807ea322SDave May   PetscFunctionBegin;
585cc843e7aSLisandro Dalcin   vbinary->skipoptions = skip;
5865c6c1daeSBarry Smith   PetscFunctionReturn(0);
5875c6c1daeSBarry Smith }
5885c6c1daeSBarry Smith 
5895c6c1daeSBarry Smith /*@
5905c6c1daeSBarry Smith     PetscViewerBinaryGetSkipOptions - checks if viewer uses the PETSc options database when loading objects
5915c6c1daeSBarry Smith 
5925c6c1daeSBarry Smith     Not Collective
5935c6c1daeSBarry Smith 
5945c6c1daeSBarry Smith     Input Parameter:
5955c6c1daeSBarry Smith .   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
5965c6c1daeSBarry Smith 
5975c6c1daeSBarry Smith     Output Parameter:
5985c6c1daeSBarry Smith .   skip - PETSC_TRUE means do not use
5995c6c1daeSBarry Smith 
6005c6c1daeSBarry Smith     Level: advanced
6015c6c1daeSBarry Smith 
60295452b02SPatrick Sanan     Notes:
60395452b02SPatrick Sanan     This must be called after PetscViewerSetType()
6045c6c1daeSBarry Smith 
6055c6c1daeSBarry Smith .seealso: PetscViewerBinaryOpen(), PetscViewerBinaryGetDescriptor(), PetscViewerBinarySkipInfo(),
6065c6c1daeSBarry Smith           PetscViewerBinarySetSkipOptions()
6075c6c1daeSBarry Smith @*/
6085c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryGetSkipOptions(PetscViewer viewer,PetscBool *skip)
6095c6c1daeSBarry Smith {
610807ea322SDave May   PetscErrorCode ierr;
6115c6c1daeSBarry Smith 
6125c6c1daeSBarry Smith   PetscFunctionBegin;
613cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
614cc843e7aSLisandro Dalcin   PetscValidBoolPointer(skip,2);
615807ea322SDave May   ierr = PetscUseMethod(viewer,"PetscViewerBinaryGetSkipOptions_C",(PetscViewer,PetscBool*),(viewer,skip));CHKERRQ(ierr);
6165c6c1daeSBarry Smith   PetscFunctionReturn(0);
6175c6c1daeSBarry Smith }
6185c6c1daeSBarry Smith 
619cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinaryGetSkipOptions_Binary(PetscViewer viewer,PetscBool *skip)
6205c6c1daeSBarry Smith {
621d21b9a37SPierre Jolivet   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data;
6225c6c1daeSBarry Smith 
6235c6c1daeSBarry Smith   PetscFunctionBegin;
624cc843e7aSLisandro Dalcin   *skip = vbinary->skipoptions;
6255c6c1daeSBarry Smith   PetscFunctionReturn(0);
6265c6c1daeSBarry Smith }
6275c6c1daeSBarry Smith 
6285c6c1daeSBarry Smith /*@
6295c6c1daeSBarry Smith     PetscViewerBinarySetSkipHeader - do not write a header with size information on output, just raw data
6305c6c1daeSBarry Smith 
6315c6c1daeSBarry Smith     Not Collective
6325c6c1daeSBarry Smith 
6335c6c1daeSBarry Smith     Input Parameters:
6345c6c1daeSBarry Smith +   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
6355c6c1daeSBarry Smith -   skip - PETSC_TRUE means do not write header
6365c6c1daeSBarry Smith 
6375c6c1daeSBarry Smith     Options Database Key:
6385c6c1daeSBarry Smith .   -viewer_binary_skip_header
6395c6c1daeSBarry Smith 
6405c6c1daeSBarry Smith     Level: advanced
6415c6c1daeSBarry Smith 
64295452b02SPatrick Sanan     Notes:
64395452b02SPatrick Sanan     This must be called after PetscViewerSetType()
6445c6c1daeSBarry Smith 
6455c6c1daeSBarry Smith            Can ONLY be called on a binary viewer
6465c6c1daeSBarry Smith 
6475c6c1daeSBarry Smith .seealso: PetscViewerBinaryOpen(), PetscViewerBinaryGetDescriptor(), PetscViewerBinarySkipInfo(),
6485c6c1daeSBarry Smith           PetscViewerBinaryGetSkipHeader()
6495c6c1daeSBarry Smith @*/
6505c6c1daeSBarry Smith PetscErrorCode PetscViewerBinarySetSkipHeader(PetscViewer viewer,PetscBool skip)
6515c6c1daeSBarry Smith {
6525c6c1daeSBarry Smith   PetscErrorCode ierr;
6535c6c1daeSBarry Smith 
6545c6c1daeSBarry Smith   PetscFunctionBegin;
655cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
656cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer,skip,2);
657cc843e7aSLisandro Dalcin   ierr = PetscTryMethod(viewer,"PetscViewerBinarySetSkipHeader_C",(PetscViewer,PetscBool),(viewer,skip));CHKERRQ(ierr);
6585c6c1daeSBarry Smith   PetscFunctionReturn(0);
6595c6c1daeSBarry Smith }
6605c6c1daeSBarry Smith 
661cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinarySetSkipHeader_Binary(PetscViewer viewer,PetscBool skip)
6625c6c1daeSBarry Smith {
6635c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data;
6645c6c1daeSBarry Smith 
6655c6c1daeSBarry Smith   PetscFunctionBegin;
666cc843e7aSLisandro Dalcin   vbinary->skipheader = skip;
6675c6c1daeSBarry Smith   PetscFunctionReturn(0);
6685c6c1daeSBarry Smith }
6695c6c1daeSBarry Smith 
6705c6c1daeSBarry Smith /*@
6715c6c1daeSBarry Smith     PetscViewerBinaryGetSkipHeader - checks whether to write a header with size information on output, or just raw data
6725c6c1daeSBarry Smith 
6735c6c1daeSBarry Smith     Not Collective
6745c6c1daeSBarry Smith 
6755c6c1daeSBarry Smith     Input Parameter:
6765c6c1daeSBarry Smith .   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
6775c6c1daeSBarry Smith 
6785c6c1daeSBarry Smith     Output Parameter:
6795c6c1daeSBarry Smith .   skip - PETSC_TRUE means do not write header
6805c6c1daeSBarry Smith 
6815c6c1daeSBarry Smith     Level: advanced
6825c6c1daeSBarry Smith 
68395452b02SPatrick Sanan     Notes:
68495452b02SPatrick Sanan     This must be called after PetscViewerSetType()
6855c6c1daeSBarry Smith 
6865c6c1daeSBarry Smith             Returns false for PETSCSOCKETVIEWER, you cannot skip the header for it.
6875c6c1daeSBarry Smith 
6885c6c1daeSBarry Smith .seealso: PetscViewerBinaryOpen(), PetscViewerBinaryGetDescriptor(), PetscViewerBinarySkipInfo(),
6895c6c1daeSBarry Smith           PetscViewerBinarySetSkipHeader()
6905c6c1daeSBarry Smith @*/
6915c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryGetSkipHeader(PetscViewer viewer,PetscBool  *skip)
6925c6c1daeSBarry Smith {
6935c6c1daeSBarry Smith   PetscErrorCode ierr;
6945c6c1daeSBarry Smith 
6955c6c1daeSBarry Smith   PetscFunctionBegin;
696cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
697cc843e7aSLisandro Dalcin   PetscValidBoolPointer(skip,2);
698163d334eSBarry Smith   ierr = PetscUseMethod(viewer,"PetscViewerBinaryGetSkipHeader_C",(PetscViewer,PetscBool*),(viewer,skip));CHKERRQ(ierr);
6995c6c1daeSBarry Smith   PetscFunctionReturn(0);
7005c6c1daeSBarry Smith }
7015c6c1daeSBarry Smith 
702cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinaryGetSkipHeader_Binary(PetscViewer viewer,PetscBool  *skip)
7035c6c1daeSBarry Smith {
704f3b211e4SSatish Balay   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data;
7055c6c1daeSBarry Smith 
7065c6c1daeSBarry Smith   PetscFunctionBegin;
707cc843e7aSLisandro Dalcin   *skip = vbinary->skipheader;
7085c6c1daeSBarry Smith   PetscFunctionReturn(0);
7095c6c1daeSBarry Smith }
7105c6c1daeSBarry Smith 
7115c6c1daeSBarry Smith /*@C
7125c6c1daeSBarry Smith     PetscViewerBinaryGetInfoPointer - Extracts the file pointer for the ASCII
7135c6c1daeSBarry Smith           info file associated with a binary file.
7145c6c1daeSBarry Smith 
7155c6c1daeSBarry Smith     Not Collective
7165c6c1daeSBarry Smith 
7175c6c1daeSBarry Smith     Input Parameter:
7185c6c1daeSBarry Smith .   viewer - PetscViewer context, obtained from PetscViewerBinaryOpen()
7195c6c1daeSBarry Smith 
7205c6c1daeSBarry Smith     Output Parameter:
7210298fd71SBarry Smith .   file - file pointer  Always returns NULL if not a binary viewer
7225c6c1daeSBarry Smith 
7235c6c1daeSBarry Smith     Level: advanced
7245c6c1daeSBarry Smith 
7255c6c1daeSBarry Smith     Notes:
7265c6c1daeSBarry Smith       For writable binary PetscViewers, the descriptor will only be valid for the
7275c6c1daeSBarry Smith     first processor in the communicator that shares the PetscViewer.
7285c6c1daeSBarry Smith 
7295c6c1daeSBarry Smith     Fortran Note:
7305c6c1daeSBarry Smith     This routine is not supported in Fortran.
7315c6c1daeSBarry Smith 
7325c6c1daeSBarry Smith .seealso: PetscViewerBinaryOpen(),PetscViewerBinaryGetDescriptor()
7335c6c1daeSBarry Smith @*/
7345c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryGetInfoPointer(PetscViewer viewer,FILE **file)
7355c6c1daeSBarry Smith {
7365c6c1daeSBarry Smith   PetscErrorCode ierr;
7375c6c1daeSBarry Smith 
7385c6c1daeSBarry Smith   PetscFunctionBegin;
739cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
740cc843e7aSLisandro Dalcin   PetscValidPointer(file,2);
7410298fd71SBarry Smith   *file = NULL;
7425c6c1daeSBarry Smith   ierr = PetscTryMethod(viewer,"PetscViewerBinaryGetInfoPointer_C",(PetscViewer,FILE **),(viewer,file));CHKERRQ(ierr);
7435c6c1daeSBarry Smith   PetscFunctionReturn(0);
7445c6c1daeSBarry Smith }
7455c6c1daeSBarry Smith 
746cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinaryGetInfoPointer_Binary(PetscViewer viewer,FILE **file)
7475c6c1daeSBarry Smith {
748cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data;
7495c6c1daeSBarry Smith   PetscErrorCode     ierr;
7505c6c1daeSBarry Smith 
7515c6c1daeSBarry Smith   PetscFunctionBegin;
752cc843e7aSLisandro Dalcin   ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr);
753cc843e7aSLisandro Dalcin   *file = vbinary->fdes_info;
754cc843e7aSLisandro Dalcin   if (viewer->format == PETSC_VIEWER_BINARY_MATLAB && !vbinary->matlabheaderwritten) {
7555c6c1daeSBarry Smith     if (vbinary->fdes_info) {
756cc843e7aSLisandro Dalcin       FILE *info = vbinary->fdes_info;
757cc843e7aSLisandro Dalcin       ierr = PetscFPrintf(PETSC_COMM_SELF,info,"#--- begin code written by PetscViewerBinary for MATLAB format ---#\n");CHKERRQ(ierr);
758cc843e7aSLisandro Dalcin       ierr = PetscFPrintf(PETSC_COMM_SELF,info,"#$$ Set.filename = '%s';\n",vbinary->filename);CHKERRQ(ierr);
759cc843e7aSLisandro Dalcin       ierr = PetscFPrintf(PETSC_COMM_SELF,info,"#$$ fd = PetscOpenFile(Set.filename);\n");CHKERRQ(ierr);
760cc843e7aSLisandro Dalcin       ierr = PetscFPrintf(PETSC_COMM_SELF,info,"#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n");CHKERRQ(ierr);
761cc843e7aSLisandro Dalcin     }
762cc843e7aSLisandro Dalcin     vbinary->matlabheaderwritten = PETSC_TRUE;
7635c6c1daeSBarry Smith   }
7645c6c1daeSBarry Smith   PetscFunctionReturn(0);
7655c6c1daeSBarry Smith }
7665c6c1daeSBarry Smith 
767cc843e7aSLisandro Dalcin 
768e0385b85SDave May #if defined(PETSC_HAVE_MPIIO)
769e0385b85SDave May static PetscErrorCode PetscViewerFileClose_BinaryMPIIO(PetscViewer v)
770e0385b85SDave May {
771e0385b85SDave May   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)v->data;
772e0385b85SDave May   PetscErrorCode     ierr;
773e0385b85SDave May 
774e0385b85SDave May   PetscFunctionBegin;
775e8a65b7dSLisandro Dalcin   if (vbinary->mfdes != MPI_FILE_NULL) {
776ffc4695bSBarry Smith     ierr = MPI_File_close(&vbinary->mfdes);CHKERRMPI(ierr);
777e0385b85SDave May   }
778e8a65b7dSLisandro Dalcin   if (vbinary->mfsub != MPI_FILE_NULL) {
779ffc4695bSBarry Smith     ierr = MPI_File_close(&vbinary->mfsub);CHKERRMPI(ierr);
780e8a65b7dSLisandro Dalcin   }
781cc843e7aSLisandro Dalcin   vbinary->moff = 0;
782e0385b85SDave May   PetscFunctionReturn(0);
783e0385b85SDave May }
784e0385b85SDave May #endif
785e0385b85SDave May 
786cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileClose_BinarySTDIO(PetscViewer v)
787cc843e7aSLisandro Dalcin {
788cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)v->data;
789cc843e7aSLisandro Dalcin   PetscErrorCode     ierr;
790cc843e7aSLisandro Dalcin 
791cc843e7aSLisandro Dalcin   PetscFunctionBegin;
792cc843e7aSLisandro Dalcin   if (vbinary->fdes != -1) {
793cc843e7aSLisandro Dalcin     ierr = PetscBinaryClose(vbinary->fdes);CHKERRQ(ierr);
794cc843e7aSLisandro Dalcin     vbinary->fdes = -1;
795cc843e7aSLisandro Dalcin     if (vbinary->storecompressed) {
796cc843e7aSLisandro Dalcin       char cmd[8+PETSC_MAX_PATH_LEN],out[64+PETSC_MAX_PATH_LEN] = "";
797cc843e7aSLisandro Dalcin       const char *gzfilename = vbinary->ogzfilename ? vbinary->ogzfilename : vbinary->filename;
798cc843e7aSLisandro Dalcin       /* compress the file */
799cc843e7aSLisandro Dalcin       ierr = PetscStrncpy(cmd,"gzip -f ",sizeof(cmd));CHKERRQ(ierr);
800cc843e7aSLisandro Dalcin       ierr = PetscStrlcat(cmd,gzfilename,sizeof(cmd));CHKERRQ(ierr);
801cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_POPEN)
802cc843e7aSLisandro Dalcin       {
803cc843e7aSLisandro Dalcin         FILE *fp;
804cc843e7aSLisandro Dalcin         ierr = PetscPOpen(PETSC_COMM_SELF,NULL,cmd,"r",&fp);CHKERRQ(ierr);
805cc843e7aSLisandro Dalcin         if (fgets(out,(int)(sizeof(out)-1),fp)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error from command %s\n%s",cmd,out);
806cc843e7aSLisandro Dalcin         ierr = PetscPClose(PETSC_COMM_SELF,fp);CHKERRQ(ierr);
807cc843e7aSLisandro Dalcin       }
808cc843e7aSLisandro Dalcin #endif
809cc843e7aSLisandro Dalcin     }
810cc843e7aSLisandro Dalcin   }
811cc843e7aSLisandro Dalcin   ierr = PetscFree(vbinary->ogzfilename);CHKERRQ(ierr);
812cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
813cc843e7aSLisandro Dalcin }
814cc843e7aSLisandro Dalcin 
815cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileClose_BinaryInfo(PetscViewer v)
816cc843e7aSLisandro Dalcin {
817cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)v->data;
818cc843e7aSLisandro Dalcin   PetscErrorCode     ierr;
819cc843e7aSLisandro Dalcin 
820cc843e7aSLisandro Dalcin   PetscFunctionBegin;
821cc843e7aSLisandro Dalcin   if (v->format == PETSC_VIEWER_BINARY_MATLAB && vbinary->matlabheaderwritten) {
822cc843e7aSLisandro Dalcin     if (vbinary->fdes_info) {
823cc843e7aSLisandro Dalcin       FILE *info = vbinary->fdes_info;
824cc843e7aSLisandro Dalcin       ierr = PetscFPrintf(PETSC_COMM_SELF,info,"#--- begin code written by PetscViewerBinary for MATLAB format ---#\n");CHKERRQ(ierr);
825cc843e7aSLisandro Dalcin       ierr = PetscFPrintf(PETSC_COMM_SELF,info,"#$$ close(fd);\n");CHKERRQ(ierr);
826cc843e7aSLisandro Dalcin       ierr = PetscFPrintf(PETSC_COMM_SELF,info,"#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n");CHKERRQ(ierr);
827cc843e7aSLisandro Dalcin     }
828cc843e7aSLisandro Dalcin   }
829cc843e7aSLisandro Dalcin   if (vbinary->fdes_info) {
830cc843e7aSLisandro Dalcin     FILE *info = vbinary->fdes_info;
831cc843e7aSLisandro Dalcin     vbinary->fdes_info = NULL;
832cc843e7aSLisandro Dalcin     if (fclose(info)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
833cc843e7aSLisandro Dalcin   }
834cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
835cc843e7aSLisandro Dalcin }
836cc843e7aSLisandro Dalcin 
837cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileClose_Binary(PetscViewer v)
838cc843e7aSLisandro Dalcin {
839cc843e7aSLisandro Dalcin   PetscErrorCode     ierr;
840cc843e7aSLisandro Dalcin 
841cc843e7aSLisandro Dalcin   PetscFunctionBegin;
842cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
843cc843e7aSLisandro Dalcin   ierr = PetscViewerFileClose_BinaryMPIIO(v);CHKERRQ(ierr);
844cc843e7aSLisandro Dalcin #endif
845cc843e7aSLisandro Dalcin   ierr = PetscViewerFileClose_BinarySTDIO(v);CHKERRQ(ierr);
846cc843e7aSLisandro Dalcin   ierr = PetscViewerFileClose_BinaryInfo(v);CHKERRQ(ierr);
847cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
848cc843e7aSLisandro Dalcin }
849cc843e7aSLisandro Dalcin 
85081f0254dSBarry Smith static PetscErrorCode PetscViewerDestroy_Binary(PetscViewer v)
8515c6c1daeSBarry Smith {
8525c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)v->data;
8535c6c1daeSBarry Smith   PetscErrorCode     ierr;
8545c6c1daeSBarry Smith 
8555c6c1daeSBarry Smith   PetscFunctionBegin;
856e0385b85SDave May   ierr = PetscViewerFileClose_Binary(v);CHKERRQ(ierr);
857f90597f1SStefano Zampini   ierr = PetscFree(vbinary->filename);CHKERRQ(ierr);
858e0385b85SDave May   ierr = PetscFree(vbinary);CHKERRQ(ierr);
85922a8f86cSLisandro Dalcin 
86022a8f86cSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetFlowControl_C",NULL);CHKERRQ(ierr);
86122a8f86cSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetFlowControl_C",NULL);CHKERRQ(ierr);
86222a8f86cSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetSkipHeader_C",NULL);CHKERRQ(ierr);
863cc843e7aSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetSkipHeader_C",NULL);CHKERRQ(ierr);
86422a8f86cSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetSkipOptions_C",NULL);CHKERRQ(ierr);
86522a8f86cSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetSkipOptions_C",NULL);CHKERRQ(ierr);
86622a8f86cSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetSkipInfo_C",NULL);CHKERRQ(ierr);
86722a8f86cSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetSkipInfo_C",NULL);CHKERRQ(ierr);
86822a8f86cSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetInfoPointer_C",NULL);CHKERRQ(ierr);
86922a8f86cSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr);
870cc843e7aSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr);
871cc843e7aSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",NULL);CHKERRQ(ierr);
872cc843e7aSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr);
87322a8f86cSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
87422a8f86cSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetUseMPIIO_C",NULL);CHKERRQ(ierr);
87522a8f86cSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetUseMPIIO_C",NULL);CHKERRQ(ierr);
87622a8f86cSLisandro Dalcin #endif
877e0385b85SDave May   PetscFunctionReturn(0);
878e0385b85SDave May }
8795c6c1daeSBarry Smith 
8805c6c1daeSBarry Smith /*@C
8815c6c1daeSBarry Smith    PetscViewerBinaryOpen - Opens a file for binary input/output.
8825c6c1daeSBarry Smith 
883d083f849SBarry Smith    Collective
8845c6c1daeSBarry Smith 
8855c6c1daeSBarry Smith    Input Parameters:
8865c6c1daeSBarry Smith +  comm - MPI communicator
8875c6c1daeSBarry Smith .  name - name of file
888cc843e7aSLisandro Dalcin -  mode - open mode of file
8895c6c1daeSBarry Smith $    FILE_MODE_WRITE - create new file for binary output
8905c6c1daeSBarry Smith $    FILE_MODE_READ - open existing file for binary input
8915c6c1daeSBarry Smith $    FILE_MODE_APPEND - open existing file for binary output
8925c6c1daeSBarry Smith 
8935c6c1daeSBarry Smith    Output Parameter:
894cc843e7aSLisandro Dalcin .  viewer - PetscViewer for binary input/output to use with the specified file
8955c6c1daeSBarry Smith 
8965c6c1daeSBarry Smith     Options Database Keys:
89763c55180SPatrick Sanan +    -viewer_binary_filename <name> -
89863c55180SPatrick Sanan .    -viewer_binary_skip_info -
89963c55180SPatrick Sanan .    -viewer_binary_skip_options -
90063c55180SPatrick Sanan .    -viewer_binary_skip_header -
90163c55180SPatrick Sanan -    -viewer_binary_mpiio -
9025c6c1daeSBarry Smith 
9035c6c1daeSBarry Smith    Level: beginner
9045c6c1daeSBarry Smith 
9055c6c1daeSBarry Smith    Note:
9065c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
9075c6c1daeSBarry Smith 
9085c6c1daeSBarry Smith     For reading files, the filename may begin with ftp:// or http:// and/or
9095c6c1daeSBarry Smith     end with .gz; in this case file is brought over and uncompressed.
9105c6c1daeSBarry Smith 
9115c6c1daeSBarry Smith     For creating files, if the file name ends with .gz it is automatically
9125c6c1daeSBarry Smith     compressed when closed.
9135c6c1daeSBarry Smith 
9146a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(),
9155c6c1daeSBarry Smith           VecView(), MatView(), VecLoad(), MatLoad(), PetscViewerBinaryGetDescriptor(),
91667918a83SBarry Smith           PetscViewerBinaryGetInfoPointer(), PetscFileMode, PetscViewer, PetscViewerBinaryRead(), PetscViewerBinarySetUseMPIIO(),
91767918a83SBarry Smith           PetscViewerBinaryGetUseMPIIO(), PetscViewerBinaryGetMPIIOOffset()
9185c6c1daeSBarry Smith @*/
919cc843e7aSLisandro Dalcin PetscErrorCode PetscViewerBinaryOpen(MPI_Comm comm,const char name[],PetscFileMode mode,PetscViewer *viewer)
9205c6c1daeSBarry Smith {
9215c6c1daeSBarry Smith   PetscErrorCode ierr;
9225c6c1daeSBarry Smith 
9235c6c1daeSBarry Smith   PetscFunctionBegin;
924cc843e7aSLisandro Dalcin   ierr = PetscViewerCreate(comm,viewer);CHKERRQ(ierr);
925cc843e7aSLisandro Dalcin   ierr = PetscViewerSetType(*viewer,PETSCVIEWERBINARY);CHKERRQ(ierr);
926cc843e7aSLisandro Dalcin   ierr = PetscViewerFileSetMode(*viewer,mode);CHKERRQ(ierr);
927cc843e7aSLisandro Dalcin   ierr = PetscViewerFileSetName(*viewer,name);CHKERRQ(ierr);
928cc843e7aSLisandro Dalcin   ierr = PetscViewerSetFromOptions(*viewer);CHKERRQ(ierr);
9295c6c1daeSBarry Smith   PetscFunctionReturn(0);
9305c6c1daeSBarry Smith }
9315c6c1daeSBarry Smith 
9325c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
933060da220SMatthew G. Knepley static PetscErrorCode PetscViewerBinaryWriteReadMPIIO(PetscViewer viewer,void *data,PetscInt num,PetscInt *count,PetscDataType dtype,PetscBool write)
9345c6c1daeSBarry Smith {
93522a8f86cSLisandro Dalcin   MPI_Comm           comm = PetscObjectComm((PetscObject)viewer);
9365c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data;
93722a8f86cSLisandro Dalcin   MPI_File           mfdes = vbinary->mfdes;
9385c6c1daeSBarry Smith   PetscErrorCode     ierr;
9395c6c1daeSBarry Smith   MPI_Datatype       mdtype;
94069e10bbaSLisandro Dalcin   PetscMPIInt        rank,cnt;
9415c6c1daeSBarry Smith   MPI_Status         status;
9425c6c1daeSBarry Smith   MPI_Aint           ul,dsize;
9435c6c1daeSBarry Smith 
9445c6c1daeSBarry Smith   PetscFunctionBegin;
945ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
946060da220SMatthew G. Knepley   ierr = PetscMPIIntCast(num,&cnt);CHKERRQ(ierr);
9475c6c1daeSBarry Smith   ierr = PetscDataTypeToMPIDataType(dtype,&mdtype);CHKERRQ(ierr);
9485c6c1daeSBarry Smith   if (write) {
94969e10bbaSLisandro Dalcin     if (!rank) {
95069e10bbaSLisandro Dalcin       ierr = MPIU_File_write_at(mfdes,vbinary->moff,data,cnt,mdtype,&status);CHKERRQ(ierr);
95169e10bbaSLisandro Dalcin     }
9525c6c1daeSBarry Smith   } else {
95369e10bbaSLisandro Dalcin     if (!rank) {
95469e10bbaSLisandro Dalcin       ierr = MPIU_File_read_at(mfdes,vbinary->moff,data,cnt,mdtype,&status);CHKERRQ(ierr);
955ffc4695bSBarry Smith       if (cnt > 0) {ierr = MPI_Get_count(&status,mdtype,&cnt);CHKERRMPI(ierr);}
9565c6c1daeSBarry Smith     }
957ffc4695bSBarry Smith     ierr = MPI_Bcast(&cnt,1,MPI_INT,0,comm);CHKERRMPI(ierr);
958ffc4695bSBarry Smith     ierr = MPI_Bcast(data,cnt,mdtype,0,comm);CHKERRMPI(ierr);
95969e10bbaSLisandro Dalcin   }
960ffc4695bSBarry Smith   ierr = MPI_Type_get_extent(mdtype,&ul,&dsize);CHKERRMPI(ierr);
9615c6c1daeSBarry Smith   vbinary->moff += dsize*cnt;
9629860990eSLisandro Dalcin   if (count) *count = cnt;
9635c6c1daeSBarry Smith   PetscFunctionReturn(0);
9645c6c1daeSBarry Smith }
9655c6c1daeSBarry Smith #endif
9665c6c1daeSBarry Smith 
9675c6c1daeSBarry Smith /*@C
9685c6c1daeSBarry Smith    PetscViewerBinaryRead - Reads from a binary file, all processors get the same result
9695c6c1daeSBarry Smith 
970d083f849SBarry Smith    Collective
9715c6c1daeSBarry Smith 
9725c6c1daeSBarry Smith    Input Parameters:
9735c6c1daeSBarry Smith +  viewer - the binary viewer
974907376e6SBarry Smith .  data - location of the data to be written
975060da220SMatthew G. Knepley .  num - number of items of data to read
976907376e6SBarry Smith -  dtype - type of data to read
9775c6c1daeSBarry Smith 
978f8e4bde8SMatthew G. Knepley    Output Parameters:
9795972f5f3SLisandro Dalcin .  count - number of items of data actually read, or NULL.
980f8e4bde8SMatthew G. Knepley 
9815c6c1daeSBarry Smith    Level: beginner
9825c6c1daeSBarry Smith 
9836a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(),
9845c6c1daeSBarry Smith           VecView(), MatView(), VecLoad(), MatLoad(), PetscViewerBinaryGetDescriptor(),
985d23dbae5SPatrick Sanan           PetscViewerBinaryGetInfoPointer(), PetscFileMode, PetscViewer, PetscViewerBinaryRead()
9865c6c1daeSBarry Smith @*/
987060da220SMatthew G. Knepley PetscErrorCode PetscViewerBinaryRead(PetscViewer viewer,void *data,PetscInt num,PetscInt *count,PetscDataType dtype)
9885c6c1daeSBarry Smith {
9895c6c1daeSBarry Smith   PetscErrorCode     ierr;
99022a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
9915c6c1daeSBarry Smith 
99222a8f86cSLisandro Dalcin   PetscFunctionBegin;
99322a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer,PETSC_VIEWER_CLASSID,1,PETSCVIEWERBINARY);
99422a8f86cSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer,num,3);
995c98fd787SBarry Smith   ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr);
99622a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary*)viewer->data;
9975c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
998bc196f7cSDave May   if (vbinary->usempiio) {
999060da220SMatthew G. Knepley     ierr = PetscViewerBinaryWriteReadMPIIO(viewer,data,num,count,dtype,PETSC_FALSE);CHKERRQ(ierr);
10005c6c1daeSBarry Smith   } else {
10015c6c1daeSBarry Smith #endif
10029860990eSLisandro Dalcin     ierr = PetscBinarySynchronizedRead(PetscObjectComm((PetscObject)viewer),vbinary->fdes,data,num,count,dtype);CHKERRQ(ierr);
10035c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
10045c6c1daeSBarry Smith   }
10055c6c1daeSBarry Smith #endif
10065c6c1daeSBarry Smith   PetscFunctionReturn(0);
10075c6c1daeSBarry Smith }
10085c6c1daeSBarry Smith 
10095c6c1daeSBarry Smith /*@C
10105c6c1daeSBarry Smith    PetscViewerBinaryWrite - writes to a binary file, only from the first process
10115c6c1daeSBarry Smith 
1012d083f849SBarry Smith    Collective
10135c6c1daeSBarry Smith 
10145c6c1daeSBarry Smith    Input Parameters:
10155c6c1daeSBarry Smith +  viewer - the binary viewer
10165c6c1daeSBarry Smith .  data - location of data
10175824c9d0SPatrick Sanan .  count - number of items of data to write
1018f253e43cSLisandro Dalcin -  dtype - type of data to write
10195c6c1daeSBarry Smith 
10205c6c1daeSBarry Smith    Level: beginner
10215c6c1daeSBarry Smith 
10226a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(),
10235c6c1daeSBarry Smith           VecView(), MatView(), VecLoad(), MatLoad(), PetscViewerBinaryGetDescriptor(), PetscDataType
1024d23dbae5SPatrick Sanan           PetscViewerBinaryGetInfoPointer(), PetscFileMode, PetscViewer, PetscViewerBinaryRead()
10255c6c1daeSBarry Smith @*/
1026f253e43cSLisandro Dalcin PetscErrorCode PetscViewerBinaryWrite(PetscViewer viewer,const void *data,PetscInt count,PetscDataType dtype)
10275c6c1daeSBarry Smith {
10285c6c1daeSBarry Smith   PetscErrorCode     ierr;
102922a8f86cSLisandro Dalcin   PetscViewer_Binary *vbinary;
10305c6c1daeSBarry Smith 
10315c6c1daeSBarry Smith   PetscFunctionBegin;
103222a8f86cSLisandro Dalcin   PetscValidHeaderSpecificType(viewer,PETSC_VIEWER_CLASSID,1,PETSCVIEWERBINARY);
103322a8f86cSLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer,count,3);
1034c98fd787SBarry Smith   ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr);
103522a8f86cSLisandro Dalcin   vbinary = (PetscViewer_Binary*)viewer->data;
10365c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
1037bc196f7cSDave May   if (vbinary->usempiio) {
1038f253e43cSLisandro Dalcin     ierr = PetscViewerBinaryWriteReadMPIIO(viewer,(void*)data,count,NULL,dtype,PETSC_TRUE);CHKERRQ(ierr);
10395c6c1daeSBarry Smith   } else {
10405c6c1daeSBarry Smith #endif
1041f253e43cSLisandro Dalcin     ierr = PetscBinarySynchronizedWrite(PetscObjectComm((PetscObject)viewer),vbinary->fdes,data,count,dtype);CHKERRQ(ierr);
10425c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
10435c6c1daeSBarry Smith   }
10445c6c1daeSBarry Smith #endif
10455c6c1daeSBarry Smith   PetscFunctionReturn(0);
10465c6c1daeSBarry Smith }
10475c6c1daeSBarry Smith 
10485972f5f3SLisandro Dalcin static PetscErrorCode PetscViewerBinaryWriteReadAll(PetscViewer viewer,PetscBool write,void *data,PetscInt count,PetscInt start,PetscInt total,PetscDataType dtype)
10495972f5f3SLisandro Dalcin {
10505972f5f3SLisandro Dalcin   MPI_Comm              comm = PetscObjectComm((PetscObject)viewer);
10515972f5f3SLisandro Dalcin   PetscMPIInt           size,rank;
10525972f5f3SLisandro Dalcin   MPI_Datatype          mdtype;
1053ec4bef21SJose E. Roman   PETSC_UNUSED MPI_Aint lb;
1054ec4bef21SJose E. Roman   MPI_Aint              dsize;
10555972f5f3SLisandro Dalcin   PetscBool             useMPIIO;
10565972f5f3SLisandro Dalcin   PetscErrorCode        ierr;
10575972f5f3SLisandro Dalcin 
10585972f5f3SLisandro Dalcin   PetscFunctionBegin;
10595972f5f3SLisandro Dalcin   PetscValidHeaderSpecificType(viewer,PETSC_VIEWER_CLASSID,1,PETSCVIEWERBINARY);
10605972f5f3SLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer,((start>=0)||(start==PETSC_DETERMINE)),4);
10615972f5f3SLisandro Dalcin   PetscValidLogicalCollectiveBool(viewer,((total>=0)||(total==PETSC_DETERMINE)),5);
10625972f5f3SLisandro Dalcin   PetscValidLogicalCollectiveInt(viewer,total,5);
10635972f5f3SLisandro Dalcin   ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr);
10645972f5f3SLisandro Dalcin 
10655972f5f3SLisandro Dalcin   ierr = PetscDataTypeToMPIDataType(dtype,&mdtype);CHKERRQ(ierr);
1066ffc4695bSBarry Smith   ierr = MPI_Type_get_extent(mdtype,&lb,&dsize);CHKERRMPI(ierr);
1067ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
1068ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
10695972f5f3SLisandro Dalcin 
10705972f5f3SLisandro Dalcin   ierr = PetscViewerBinaryGetUseMPIIO(viewer,&useMPIIO);CHKERRQ(ierr);
10715972f5f3SLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
10725972f5f3SLisandro Dalcin   if (useMPIIO) {
10735972f5f3SLisandro Dalcin     MPI_File       mfdes;
10745972f5f3SLisandro Dalcin     MPI_Offset     off;
10755972f5f3SLisandro Dalcin     PetscMPIInt    cnt;
10765972f5f3SLisandro Dalcin 
10775972f5f3SLisandro Dalcin     if (start == PETSC_DETERMINE) {
1078ffc4695bSBarry Smith       ierr = MPI_Scan(&count,&start,1,MPIU_INT,MPI_SUM,comm);CHKERRMPI(ierr);
10795972f5f3SLisandro Dalcin       start -= count;
10805972f5f3SLisandro Dalcin     }
10815972f5f3SLisandro Dalcin     if (total == PETSC_DETERMINE) {
10825972f5f3SLisandro Dalcin       total = start + count;
1083ffc4695bSBarry Smith       ierr = MPI_Bcast(&total,1,MPIU_INT,size-1,comm);CHKERRMPI(ierr);
10845972f5f3SLisandro Dalcin     }
10855972f5f3SLisandro Dalcin     ierr = PetscMPIIntCast(count,&cnt);CHKERRQ(ierr);
10865972f5f3SLisandro Dalcin     ierr = PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes);CHKERRQ(ierr);
10875972f5f3SLisandro Dalcin     ierr = PetscViewerBinaryGetMPIIOOffset(viewer,&off);CHKERRQ(ierr);
10885972f5f3SLisandro Dalcin     off += (MPI_Offset)(start*dsize);
10895972f5f3SLisandro Dalcin     if (write) {
10905972f5f3SLisandro Dalcin       ierr = MPIU_File_write_at_all(mfdes,off,data,cnt,mdtype,MPI_STATUS_IGNORE);CHKERRQ(ierr);
10915972f5f3SLisandro Dalcin     } else {
10925972f5f3SLisandro Dalcin       ierr = MPIU_File_read_at_all(mfdes,off,data,cnt,mdtype,MPI_STATUS_IGNORE);CHKERRQ(ierr);
10935972f5f3SLisandro Dalcin     }
10945972f5f3SLisandro Dalcin     off  = (MPI_Offset)(total*dsize);
10955972f5f3SLisandro Dalcin     ierr = PetscViewerBinaryAddMPIIOOffset(viewer,off);CHKERRQ(ierr);
10965972f5f3SLisandro Dalcin     PetscFunctionReturn(0);
10975972f5f3SLisandro Dalcin   }
10985972f5f3SLisandro Dalcin #endif
10995972f5f3SLisandro Dalcin   {
11005972f5f3SLisandro Dalcin     int         fdes;
11015972f5f3SLisandro Dalcin     char        *workbuf = NULL;
11025ff7be28SBarry Smith     PetscInt    tcount = !rank ? 0 : count,maxcount=0,message_count,flowcontrolcount;
11035972f5f3SLisandro Dalcin     PetscMPIInt tag,cnt,maxcnt,scnt=0,rcnt=0,j;
11045972f5f3SLisandro Dalcin     MPI_Status  status;
11055972f5f3SLisandro Dalcin 
11065972f5f3SLisandro Dalcin     ierr = PetscCommGetNewTag(comm,&tag);CHKERRQ(ierr);
1107ffc4695bSBarry Smith     ierr = MPI_Reduce(&tcount,&maxcount,1,MPIU_INT,MPI_MAX,0,comm);CHKERRMPI(ierr);
11085972f5f3SLisandro Dalcin     ierr = PetscMPIIntCast(maxcount,&maxcnt);CHKERRQ(ierr);
11095972f5f3SLisandro Dalcin     ierr = PetscMPIIntCast(count,&cnt);CHKERRQ(ierr);
11105972f5f3SLisandro Dalcin 
11115972f5f3SLisandro Dalcin     ierr = PetscViewerBinaryGetDescriptor(viewer,&fdes);CHKERRQ(ierr);
11125972f5f3SLisandro Dalcin     ierr = PetscViewerFlowControlStart(viewer,&message_count,&flowcontrolcount);CHKERRQ(ierr);
11135972f5f3SLisandro Dalcin     if (!rank) {
11145972f5f3SLisandro Dalcin       ierr = PetscMalloc(maxcnt*dsize,&workbuf);CHKERRQ(ierr);
11155972f5f3SLisandro Dalcin       if (write) {
11165972f5f3SLisandro Dalcin         ierr = PetscBinaryWrite(fdes,data,cnt,dtype);CHKERRQ(ierr);
11175972f5f3SLisandro Dalcin       } else {
11185972f5f3SLisandro Dalcin         ierr = PetscBinaryRead(fdes,data,cnt,NULL,dtype);CHKERRQ(ierr);
11195972f5f3SLisandro Dalcin       }
11205972f5f3SLisandro Dalcin       for (j=1; j<size; j++) {
1121*9dddd249SSatish Balay         ierr = PetscViewerFlowControlStepMain(viewer,j,&message_count,flowcontrolcount);CHKERRQ(ierr);
11225972f5f3SLisandro Dalcin         if (write) {
1123ffc4695bSBarry Smith           ierr = MPI_Recv(workbuf,maxcnt,mdtype,j,tag,comm,&status);CHKERRMPI(ierr);
1124ffc4695bSBarry Smith           ierr = MPI_Get_count(&status,mdtype,&rcnt);CHKERRMPI(ierr);
11255972f5f3SLisandro Dalcin           ierr = PetscBinaryWrite(fdes,workbuf,rcnt,dtype);CHKERRQ(ierr);
11265972f5f3SLisandro Dalcin         } else {
1127ffc4695bSBarry Smith           ierr = MPI_Recv(&scnt,1,MPI_INT,j,tag,comm,MPI_STATUS_IGNORE);CHKERRMPI(ierr);
11285972f5f3SLisandro Dalcin           ierr = PetscBinaryRead(fdes,workbuf,scnt,NULL,dtype);CHKERRQ(ierr);
1129ffc4695bSBarry Smith           ierr = MPI_Send(workbuf,scnt,mdtype,j,tag,comm);CHKERRMPI(ierr);
11305972f5f3SLisandro Dalcin         }
11315972f5f3SLisandro Dalcin       }
11325972f5f3SLisandro Dalcin       ierr = PetscFree(workbuf);CHKERRQ(ierr);
1133*9dddd249SSatish Balay       ierr = PetscViewerFlowControlEndMain(viewer,&message_count);CHKERRQ(ierr);
11345972f5f3SLisandro Dalcin     } else {
11355972f5f3SLisandro Dalcin       ierr = PetscViewerFlowControlStepWorker(viewer,rank,&message_count);CHKERRQ(ierr);
11365972f5f3SLisandro Dalcin       if (write) {
1137ffc4695bSBarry Smith         ierr = MPI_Send(data,cnt,mdtype,0,tag,comm);CHKERRMPI(ierr);
11385972f5f3SLisandro Dalcin       } else {
1139ffc4695bSBarry Smith         ierr = MPI_Send(&cnt,1,MPI_INT,0,tag,comm);CHKERRMPI(ierr);
1140ffc4695bSBarry Smith         ierr = MPI_Recv(data,cnt,mdtype,0,tag,comm,MPI_STATUS_IGNORE);CHKERRMPI(ierr);
11415972f5f3SLisandro Dalcin       }
11425972f5f3SLisandro Dalcin       ierr = PetscViewerFlowControlEndWorker(viewer,&message_count);CHKERRQ(ierr);
11435972f5f3SLisandro Dalcin     }
11445972f5f3SLisandro Dalcin   }
11455972f5f3SLisandro Dalcin   PetscFunctionReturn(0);
11465972f5f3SLisandro Dalcin }
11475972f5f3SLisandro Dalcin 
11485972f5f3SLisandro Dalcin /*@C
11495972f5f3SLisandro Dalcin    PetscViewerBinaryReadAll - reads from a binary file from all processes
11505972f5f3SLisandro Dalcin 
11515972f5f3SLisandro Dalcin    Collective
11525972f5f3SLisandro Dalcin 
11535972f5f3SLisandro Dalcin    Input Parameters:
11545972f5f3SLisandro Dalcin +  viewer - the binary viewer
11555972f5f3SLisandro Dalcin .  data - location of data
11565972f5f3SLisandro Dalcin .  count - local number of items of data to read
11575972f5f3SLisandro Dalcin .  start - local start, can be PETSC_DETERMINE
11585972f5f3SLisandro Dalcin .  total - global number of items of data to read, can be PETSC_DETERMINE
11595972f5f3SLisandro Dalcin -  dtype - type of data to read
11605972f5f3SLisandro Dalcin 
11615972f5f3SLisandro Dalcin    Level: advanced
11625972f5f3SLisandro Dalcin 
1163d23dbae5SPatrick Sanan .seealso: PetscViewerBinaryOpen(), PetscViewerBinarySetUseMPIIO(), PetscViewerBinaryRead(), PetscViewerBinaryWriteAll()
11645972f5f3SLisandro Dalcin @*/
11655972f5f3SLisandro Dalcin PetscErrorCode PetscViewerBinaryReadAll(PetscViewer viewer,void *data,PetscInt count,PetscInt start,PetscInt total,PetscDataType dtype)
11665972f5f3SLisandro Dalcin {
11675972f5f3SLisandro Dalcin   PetscErrorCode ierr;
11685972f5f3SLisandro Dalcin   PetscFunctionBegin;
11695972f5f3SLisandro Dalcin   ierr = PetscViewerBinaryWriteReadAll(viewer,PETSC_FALSE,data,count,start,total,dtype);CHKERRQ(ierr);
11705972f5f3SLisandro Dalcin   PetscFunctionReturn(0);
11715972f5f3SLisandro Dalcin }
11725972f5f3SLisandro Dalcin 
11735972f5f3SLisandro Dalcin /*@C
11745972f5f3SLisandro Dalcin    PetscViewerBinaryWriteAll - writes to a binary file from all processes
11755972f5f3SLisandro Dalcin 
11765972f5f3SLisandro Dalcin    Collective
11775972f5f3SLisandro Dalcin 
11785972f5f3SLisandro Dalcin    Input Parameters:
11795972f5f3SLisandro Dalcin +  viewer - the binary viewer
11805972f5f3SLisandro Dalcin .  data - location of data
11815972f5f3SLisandro Dalcin .  count - local number of items of data to write
11825972f5f3SLisandro Dalcin .  start - local start, can be PETSC_DETERMINE
11835972f5f3SLisandro Dalcin .  total - global number of items of data to write, can be PETSC_DETERMINE
11845972f5f3SLisandro Dalcin -  dtype - type of data to write
11855972f5f3SLisandro Dalcin 
11865972f5f3SLisandro Dalcin    Level: advanced
11875972f5f3SLisandro Dalcin 
1188d23dbae5SPatrick Sanan .seealso: PetscViewerBinaryOpen(), PetscViewerBinarySetUseMPIIO(), PetscViewerBinaryWriteAll(), PetscViewerBinaryReadAll()
11895972f5f3SLisandro Dalcin @*/
11905972f5f3SLisandro Dalcin PetscErrorCode PetscViewerBinaryWriteAll(PetscViewer viewer,const void *data,PetscInt count,PetscInt start,PetscInt total,PetscDataType dtype)
11915972f5f3SLisandro Dalcin {
11925972f5f3SLisandro Dalcin   PetscErrorCode ierr;
11935972f5f3SLisandro Dalcin   PetscFunctionBegin;
11945972f5f3SLisandro Dalcin   ierr = PetscViewerBinaryWriteReadAll(viewer,PETSC_TRUE,(void*)data,count,start,total,dtype);CHKERRQ(ierr);
11955972f5f3SLisandro Dalcin   PetscFunctionReturn(0);
11965972f5f3SLisandro Dalcin }
11975972f5f3SLisandro Dalcin 
11985c6c1daeSBarry Smith /*@C
11995c6c1daeSBarry Smith    PetscViewerBinaryWriteStringArray - writes to a binary file, only from the first process an array of strings
12005c6c1daeSBarry Smith 
1201d083f849SBarry Smith    Collective
12025c6c1daeSBarry Smith 
12035c6c1daeSBarry Smith    Input Parameters:
12045c6c1daeSBarry Smith +  viewer - the binary viewer
12055c6c1daeSBarry Smith -  data - location of the array of strings
12065c6c1daeSBarry Smith 
12075c6c1daeSBarry Smith 
12085c6c1daeSBarry Smith    Level: intermediate
12095c6c1daeSBarry Smith 
121095452b02SPatrick Sanan     Notes:
121195452b02SPatrick Sanan     array of strings is null terminated
12125c6c1daeSBarry Smith 
12136a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(),
12145c6c1daeSBarry Smith           VecView(), MatView(), VecLoad(), MatLoad(), PetscViewerBinaryGetDescriptor(),
1215d23dbae5SPatrick Sanan           PetscViewerBinaryGetInfoPointer(), PetscFileMode, PetscViewer, PetscViewerBinaryRead()
12165c6c1daeSBarry Smith @*/
121778fbdcc8SBarry Smith PetscErrorCode PetscViewerBinaryWriteStringArray(PetscViewer viewer,const char * const *data)
12185c6c1daeSBarry Smith {
12195c6c1daeSBarry Smith   PetscErrorCode ierr;
12205c6c1daeSBarry Smith   PetscInt       i,n = 0,*sizes;
1221cc843e7aSLisandro Dalcin   size_t         len;
12225c6c1daeSBarry Smith 
122322a8f86cSLisandro Dalcin   PetscFunctionBegin;
1224c98fd787SBarry Smith   ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr);
12255c6c1daeSBarry Smith   /* count number of strings */
12265c6c1daeSBarry Smith   while (data[n++]);
12275c6c1daeSBarry Smith   n--;
1228854ce69bSBarry Smith   ierr = PetscMalloc1(n+1,&sizes);CHKERRQ(ierr);
12295c6c1daeSBarry Smith   sizes[0] = n;
12305c6c1daeSBarry Smith   for (i=0; i<n; i++) {
1231cc843e7aSLisandro Dalcin     ierr = PetscStrlen(data[i],&len);CHKERRQ(ierr);
1232cc843e7aSLisandro Dalcin     sizes[i+1] = (PetscInt)len + 1; /* size includes space for the null terminator */
12335c6c1daeSBarry Smith   }
1234f253e43cSLisandro Dalcin   ierr = PetscViewerBinaryWrite(viewer,sizes,n+1,PETSC_INT);CHKERRQ(ierr);
12355c6c1daeSBarry Smith   for (i=0; i<n; i++) {
1236f253e43cSLisandro Dalcin     ierr = PetscViewerBinaryWrite(viewer,(void*)data[i],sizes[i+1],PETSC_CHAR);CHKERRQ(ierr);
12375c6c1daeSBarry Smith   }
12385c6c1daeSBarry Smith   ierr = PetscFree(sizes);CHKERRQ(ierr);
12395c6c1daeSBarry Smith   PetscFunctionReturn(0);
12405c6c1daeSBarry Smith }
12415c6c1daeSBarry Smith 
12425c6c1daeSBarry Smith /*@C
12435c6c1daeSBarry Smith    PetscViewerBinaryReadStringArray - reads a binary file an array of strings
12445c6c1daeSBarry Smith 
1245d083f849SBarry Smith    Collective
12465c6c1daeSBarry Smith 
12475c6c1daeSBarry Smith    Input Parameter:
12485c6c1daeSBarry Smith .  viewer - the binary viewer
12495c6c1daeSBarry Smith 
12505c6c1daeSBarry Smith    Output Parameter:
12515c6c1daeSBarry Smith .  data - location of the array of strings
12525c6c1daeSBarry Smith 
12535c6c1daeSBarry Smith    Level: intermediate
12545c6c1daeSBarry Smith 
125595452b02SPatrick Sanan     Notes:
125695452b02SPatrick Sanan     array of strings is null terminated
12575c6c1daeSBarry Smith 
12586a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(),
12595c6c1daeSBarry Smith           VecView(), MatView(), VecLoad(), MatLoad(), PetscViewerBinaryGetDescriptor(),
1260d23dbae5SPatrick Sanan           PetscViewerBinaryGetInfoPointer(), PetscFileMode, PetscViewer, PetscViewerBinaryRead()
12615c6c1daeSBarry Smith @*/
12625c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryReadStringArray(PetscViewer viewer,char ***data)
12635c6c1daeSBarry Smith {
12645c6c1daeSBarry Smith   PetscErrorCode ierr;
1265060da220SMatthew G. Knepley   PetscInt       i,n,*sizes,N = 0;
12665c6c1daeSBarry Smith 
126722a8f86cSLisandro Dalcin   PetscFunctionBegin;
1268c98fd787SBarry Smith   ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr);
12695c6c1daeSBarry Smith   /* count number of strings */
1270060da220SMatthew G. Knepley   ierr = PetscViewerBinaryRead(viewer,&n,1,NULL,PETSC_INT);CHKERRQ(ierr);
1271785e854fSJed Brown   ierr = PetscMalloc1(n,&sizes);CHKERRQ(ierr);
1272060da220SMatthew G. Knepley   ierr = PetscViewerBinaryRead(viewer,sizes,n,NULL,PETSC_INT);CHKERRQ(ierr);
1273a297a907SKarl Rupp   for (i=0; i<n; i++) N += sizes[i];
1274854ce69bSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(char*) + N*sizeof(char),data);CHKERRQ(ierr);
12755c6c1daeSBarry Smith   (*data)[0] = (char*)((*data) + n + 1);
1276a297a907SKarl Rupp   for (i=1; i<n; i++) (*data)[i] = (*data)[i-1] + sizes[i-1];
1277060da220SMatthew G. Knepley   ierr = PetscViewerBinaryRead(viewer,(*data)[0],N,NULL,PETSC_CHAR);CHKERRQ(ierr);
127802c9f0b5SLisandro Dalcin   (*data)[n] = NULL;
12795c6c1daeSBarry Smith   ierr = PetscFree(sizes);CHKERRQ(ierr);
12805c6c1daeSBarry Smith   PetscFunctionReturn(0);
12815c6c1daeSBarry Smith }
12825c6c1daeSBarry Smith 
1283cc843e7aSLisandro Dalcin /*@C
1284cc843e7aSLisandro Dalcin      PetscViewerFileSetMode - Sets the open mode of file
1285cc843e7aSLisandro Dalcin 
1286cc843e7aSLisandro Dalcin     Logically Collective on PetscViewer
1287cc843e7aSLisandro Dalcin 
1288cc843e7aSLisandro Dalcin   Input Parameters:
1289cc843e7aSLisandro Dalcin +  viewer - the PetscViewer; must be a a PETSCVIEWERBINARY, PETSCVIEWERMATLAB, PETSCVIEWERHDF5, or PETSCVIEWERASCII  PetscViewer
1290cc843e7aSLisandro Dalcin -  mode - open mode of file
1291cc843e7aSLisandro Dalcin $    FILE_MODE_WRITE - create new file for output
1292cc843e7aSLisandro Dalcin $    FILE_MODE_READ - open existing file for input
1293cc843e7aSLisandro Dalcin $    FILE_MODE_APPEND - open existing file for output
1294cc843e7aSLisandro Dalcin 
1295cc843e7aSLisandro Dalcin   Level: advanced
1296cc843e7aSLisandro Dalcin 
1297cc843e7aSLisandro Dalcin .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
1298cc843e7aSLisandro Dalcin 
1299cc843e7aSLisandro Dalcin @*/
1300cc843e7aSLisandro Dalcin PetscErrorCode PetscViewerFileSetMode(PetscViewer viewer,PetscFileMode mode)
1301cc843e7aSLisandro Dalcin {
1302cc843e7aSLisandro Dalcin   PetscErrorCode ierr;
1303cc843e7aSLisandro Dalcin 
1304cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1305cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1306cc843e7aSLisandro Dalcin   PetscValidLogicalCollectiveEnum(viewer,mode,2);
1307cc843e7aSLisandro Dalcin   ierr = PetscTryMethod(viewer,"PetscViewerFileSetMode_C",(PetscViewer,PetscFileMode),(viewer,mode));CHKERRQ(ierr);
1308cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1309cc843e7aSLisandro Dalcin }
1310cc843e7aSLisandro Dalcin 
1311cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileSetMode_Binary(PetscViewer viewer,PetscFileMode mode)
1312cc843e7aSLisandro Dalcin {
1313cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data;
1314cc843e7aSLisandro Dalcin 
1315cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1316cc843e7aSLisandro Dalcin   if (viewer->setupcalled && vbinary->filemode != mode) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER,"Cannot change mode to %s after setup",PetscFileModes[mode]);
1317cc843e7aSLisandro Dalcin   vbinary->filemode = mode;
1318cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1319cc843e7aSLisandro Dalcin }
1320cc843e7aSLisandro Dalcin 
1321cc843e7aSLisandro Dalcin /*@C
1322cc843e7aSLisandro Dalcin      PetscViewerFileGetMode - Gets the open mode of file
1323cc843e7aSLisandro Dalcin 
1324cc843e7aSLisandro Dalcin     Not Collective
1325cc843e7aSLisandro Dalcin 
1326cc843e7aSLisandro Dalcin   Input Parameter:
1327cc843e7aSLisandro Dalcin .  viewer - the PetscViewer; must be a PETSCVIEWERBINARY, PETSCVIEWERMATLAB, PETSCVIEWERHDF5, or PETSCVIEWERASCII  PetscViewer
1328cc843e7aSLisandro Dalcin 
1329cc843e7aSLisandro Dalcin   Output Parameter:
1330cc843e7aSLisandro Dalcin .  mode - open mode of file
1331cc843e7aSLisandro Dalcin $    FILE_MODE_WRITE - create new file for binary output
1332cc843e7aSLisandro Dalcin $    FILE_MODE_READ - open existing file for binary input
1333cc843e7aSLisandro Dalcin $    FILE_MODE_APPEND - open existing file for binary output
1334cc843e7aSLisandro Dalcin 
1335cc843e7aSLisandro Dalcin   Level: advanced
1336cc843e7aSLisandro Dalcin 
1337cc843e7aSLisandro Dalcin .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
1338cc843e7aSLisandro Dalcin 
1339cc843e7aSLisandro Dalcin @*/
1340cc843e7aSLisandro Dalcin PetscErrorCode PetscViewerFileGetMode(PetscViewer viewer,PetscFileMode *mode)
1341cc843e7aSLisandro Dalcin {
1342cc843e7aSLisandro Dalcin   PetscErrorCode ierr;
1343cc843e7aSLisandro Dalcin 
1344cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1345cc843e7aSLisandro Dalcin   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1346cc843e7aSLisandro Dalcin   PetscValidPointer(mode,2);
1347cc843e7aSLisandro Dalcin   ierr = PetscUseMethod(viewer,"PetscViewerFileGetMode_C",(PetscViewer,PetscFileMode*),(viewer,mode));CHKERRQ(ierr);
1348cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1349cc843e7aSLisandro Dalcin }
1350cc843e7aSLisandro Dalcin 
1351cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileGetMode_Binary(PetscViewer viewer,PetscFileMode *mode)
1352cc843e7aSLisandro Dalcin {
1353cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data;
1354cc843e7aSLisandro Dalcin 
1355cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1356cc843e7aSLisandro Dalcin   *mode = vbinary->filemode;
1357cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1358cc843e7aSLisandro Dalcin }
1359cc843e7aSLisandro Dalcin 
1360cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileSetName_Binary(PetscViewer viewer,const char name[])
1361cc843e7aSLisandro Dalcin {
1362cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data;
1363cc843e7aSLisandro Dalcin   PetscErrorCode     ierr;
1364cc843e7aSLisandro Dalcin 
1365cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1366cc843e7aSLisandro Dalcin   if (viewer->setupcalled && vbinary->filename) {
1367cc843e7aSLisandro Dalcin     /* gzip can be run after the file with the previous filename has been closed */
1368cc843e7aSLisandro Dalcin     ierr = PetscFree(vbinary->ogzfilename);CHKERRQ(ierr);
1369cc843e7aSLisandro Dalcin     ierr = PetscStrallocpy(vbinary->filename,&vbinary->ogzfilename);CHKERRQ(ierr);
1370cc843e7aSLisandro Dalcin   }
1371cc843e7aSLisandro Dalcin   ierr = PetscFree(vbinary->filename);CHKERRQ(ierr);
1372cc843e7aSLisandro Dalcin   ierr = PetscStrallocpy(name,&vbinary->filename);CHKERRQ(ierr);
1373cc843e7aSLisandro Dalcin   viewer->setupcalled = PETSC_FALSE;
1374cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1375cc843e7aSLisandro Dalcin }
1376cc843e7aSLisandro Dalcin 
137781f0254dSBarry Smith static PetscErrorCode PetscViewerFileGetName_Binary(PetscViewer viewer,const char **name)
13785c6c1daeSBarry Smith {
13795c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data;
13805c6c1daeSBarry Smith 
13815c6c1daeSBarry Smith   PetscFunctionBegin;
13825c6c1daeSBarry Smith   *name = vbinary->filename;
13835c6c1daeSBarry Smith   PetscFunctionReturn(0);
13845c6c1daeSBarry Smith }
13855c6c1daeSBarry Smith 
13865c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
1387e0385b85SDave May static PetscErrorCode PetscViewerFileSetUp_BinaryMPIIO(PetscViewer viewer)
13885c6c1daeSBarry Smith {
13895c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data;
1390e8a65b7dSLisandro Dalcin   int                amode;
1391cc843e7aSLisandro Dalcin   PetscErrorCode     ierr;
13925c6c1daeSBarry Smith 
13935c6c1daeSBarry Smith   PetscFunctionBegin;
1394a297a907SKarl Rupp   vbinary->storecompressed = PETSC_FALSE;
13955c6c1daeSBarry Smith 
1396cc843e7aSLisandro Dalcin   vbinary->moff = 0;
1397cc843e7aSLisandro Dalcin   switch (vbinary->filemode) {
1398e8a65b7dSLisandro Dalcin   case FILE_MODE_READ:   amode = MPI_MODE_RDONLY; break;
1399e8a65b7dSLisandro Dalcin   case FILE_MODE_WRITE:  amode = MPI_MODE_WRONLY | MPI_MODE_CREATE; break;
140022a8f86cSLisandro Dalcin   case FILE_MODE_APPEND: amode = MPI_MODE_WRONLY | MPI_MODE_CREATE | MPI_MODE_APPEND; break;
1401cc843e7aSLisandro Dalcin   default: SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Unsupported file mode %s",PetscFileModes[vbinary->filemode]);
14025c6c1daeSBarry Smith   }
1403ffc4695bSBarry Smith   ierr = MPI_File_open(PetscObjectComm((PetscObject)viewer),vbinary->filename,amode,MPI_INFO_NULL,&vbinary->mfdes);CHKERRMPI(ierr);
140422a8f86cSLisandro Dalcin   /*
140522a8f86cSLisandro Dalcin       The MPI standard does not have MPI_MODE_TRUNCATE. We emulate this behavior by setting the file size to zero.
140622a8f86cSLisandro Dalcin   */
1407ffc4695bSBarry Smith   if (vbinary->filemode == FILE_MODE_WRITE) {ierr = MPI_File_set_size(vbinary->mfdes,0);CHKERRMPI(ierr);}
140822a8f86cSLisandro Dalcin   /*
140922a8f86cSLisandro Dalcin       Initially, all processes view the file as a linear byte stream. Therefore, for files opened with MPI_MODE_APPEND,
141022a8f86cSLisandro Dalcin       MPI_File_get_position[_shared](fh, &offset) returns the absolute byte position at the end of file.
141122a8f86cSLisandro Dalcin       Otherwise, we would need to call MPI_File_get_byte_offset(fh, offset, &byte_offset) to convert
141222a8f86cSLisandro Dalcin       the offset in etype units to an absolute byte position.
141322a8f86cSLisandro Dalcin    */
1414ffc4695bSBarry Smith   if (vbinary->filemode == FILE_MODE_APPEND) {ierr = MPI_File_get_position(vbinary->mfdes,&vbinary->moff);CHKERRMPI(ierr);}
1415cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1416cc843e7aSLisandro Dalcin }
1417cc843e7aSLisandro Dalcin #endif
14185c6c1daeSBarry Smith 
1419cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileSetUp_BinarySTDIO(PetscViewer viewer)
1420cc843e7aSLisandro Dalcin {
1421cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data;
1422cc843e7aSLisandro Dalcin   const char         *fname;
1423cc843e7aSLisandro Dalcin   char               bname[PETSC_MAX_PATH_LEN],*gz;
1424cc843e7aSLisandro Dalcin   PetscBool          found;
1425cc843e7aSLisandro Dalcin   PetscMPIInt        rank;
1426cc843e7aSLisandro Dalcin   PetscErrorCode     ierr;
14275c6c1daeSBarry Smith 
1428cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1429ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank);CHKERRMPI(ierr);
14305c6c1daeSBarry Smith 
1431cc843e7aSLisandro Dalcin   /* if file name ends in .gz strip that off and note user wants file compressed */
1432cc843e7aSLisandro Dalcin   vbinary->storecompressed = PETSC_FALSE;
1433cc843e7aSLisandro Dalcin   if (vbinary->filemode == FILE_MODE_WRITE) {
1434cc843e7aSLisandro Dalcin     ierr = PetscStrstr(vbinary->filename,".gz",&gz);CHKERRQ(ierr);
1435cc843e7aSLisandro Dalcin     if (gz && gz[3] == 0) {*gz = 0; vbinary->storecompressed = PETSC_TRUE;}
1436cc843e7aSLisandro Dalcin   }
1437cc843e7aSLisandro Dalcin #if !defined(PETSC_HAVE_POPEN)
1438cc843e7aSLisandro Dalcin   if (vbinary->storecompressed) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP_SYS,"Cannot run gzip on this machine");
1439cc843e7aSLisandro Dalcin #endif
1440cc843e7aSLisandro Dalcin 
1441cc843e7aSLisandro Dalcin 
1442cc843e7aSLisandro Dalcin   fname = vbinary->filename;
1443cc843e7aSLisandro Dalcin   if (vbinary->filemode == FILE_MODE_READ) { /* possibly get the file from remote site or compressed file */
1444cc843e7aSLisandro Dalcin     ierr  = PetscFileRetrieve(PetscObjectComm((PetscObject)viewer),fname,bname,PETSC_MAX_PATH_LEN,&found);CHKERRQ(ierr);
1445cc843e7aSLisandro Dalcin     if (!found) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_OPEN,"Cannot locate file: %s",fname);
1446cc843e7aSLisandro Dalcin     fname = bname;
14475c6c1daeSBarry Smith   }
14485c6c1daeSBarry Smith 
1449cc843e7aSLisandro Dalcin   vbinary->fdes = -1;
1450cc843e7aSLisandro Dalcin   if (!rank) { /* only first processor opens file*/
1451cc843e7aSLisandro Dalcin     PetscFileMode mode = vbinary->filemode;
1452cc843e7aSLisandro Dalcin     if (mode == FILE_MODE_APPEND) {
1453cc843e7aSLisandro Dalcin       /* check if asked to append to a non-existing file */
1454cc843e7aSLisandro Dalcin       ierr = PetscTestFile(fname,'\0',&found);CHKERRQ(ierr);
1455cc843e7aSLisandro Dalcin       if (!found) mode = FILE_MODE_WRITE;
1456cc843e7aSLisandro Dalcin     }
1457cc843e7aSLisandro Dalcin     ierr = PetscBinaryOpen(fname,mode,&vbinary->fdes);CHKERRQ(ierr);
1458cc843e7aSLisandro Dalcin   }
1459cc843e7aSLisandro Dalcin   PetscFunctionReturn(0);
1460cc843e7aSLisandro Dalcin }
1461cc843e7aSLisandro Dalcin 
1462cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileSetUp_BinaryInfo(PetscViewer viewer)
1463cc843e7aSLisandro Dalcin {
1464cc843e7aSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data;
1465cc843e7aSLisandro Dalcin   PetscMPIInt        rank;
1466cc843e7aSLisandro Dalcin   PetscBool          found;
1467cc843e7aSLisandro Dalcin   PetscErrorCode     ierr;
1468cc843e7aSLisandro Dalcin 
1469cc843e7aSLisandro Dalcin   PetscFunctionBegin;
1470cc843e7aSLisandro Dalcin   vbinary->fdes_info = NULL;
1471ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank);CHKERRMPI(ierr);
1472cc843e7aSLisandro Dalcin   if (!vbinary->skipinfo && (vbinary->filemode == FILE_MODE_READ || !rank)) {
1473cc843e7aSLisandro Dalcin     char infoname[PETSC_MAX_PATH_LEN],iname[PETSC_MAX_PATH_LEN],*gz;
1474cc843e7aSLisandro Dalcin 
1475cc843e7aSLisandro Dalcin     ierr = PetscStrncpy(infoname,vbinary->filename,sizeof(infoname));CHKERRQ(ierr);
1476cc843e7aSLisandro Dalcin     /* remove .gz if it ends file name */
1477cc843e7aSLisandro Dalcin     ierr = PetscStrstr(infoname,".gz",&gz);CHKERRQ(ierr);
1478cc843e7aSLisandro Dalcin     if (gz && gz[3] == 0) *gz = 0;
1479cc843e7aSLisandro Dalcin 
1480a126751eSBarry Smith     ierr = PetscStrlcat(infoname,".info",sizeof(infoname));CHKERRQ(ierr);
1481cc843e7aSLisandro Dalcin     if (vbinary->filemode == FILE_MODE_READ) {
14825c6c1daeSBarry Smith       ierr = PetscFixFilename(infoname,iname);CHKERRQ(ierr);
1483ce94432eSBarry Smith       ierr = PetscFileRetrieve(PetscObjectComm((PetscObject)viewer),iname,infoname,PETSC_MAX_PATH_LEN,&found);CHKERRQ(ierr);
1484cc843e7aSLisandro Dalcin       if (found) {ierr = PetscOptionsInsertFile(PetscObjectComm((PetscObject)viewer),((PetscObject)viewer)->options,infoname,PETSC_FALSE);CHKERRQ(ierr);}
1485cc843e7aSLisandro Dalcin     } else if (!rank) { /* write or append */
1486cc843e7aSLisandro Dalcin       const char *omode = (vbinary->filemode == FILE_MODE_APPEND) ? "a" : "w";
1487cc843e7aSLisandro Dalcin       vbinary->fdes_info = fopen(infoname,omode);
14885c6c1daeSBarry Smith       if (!vbinary->fdes_info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open .info file %s for writing",infoname);
14895c6c1daeSBarry Smith     }
14905c6c1daeSBarry Smith   }
14915c6c1daeSBarry Smith   PetscFunctionReturn(0);
14925c6c1daeSBarry Smith }
14935c6c1daeSBarry Smith 
1494cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerSetUp_Binary(PetscViewer viewer)
14955c6c1daeSBarry Smith {
14965c6c1daeSBarry Smith   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data;
1497cc843e7aSLisandro Dalcin   PetscBool          usempiio;
1498cc843e7aSLisandro Dalcin   PetscErrorCode     ierr;
1499cc843e7aSLisandro Dalcin 
15005c6c1daeSBarry Smith   PetscFunctionBegin;
1501cc843e7aSLisandro Dalcin   if (!vbinary->setfromoptionscalled) {ierr = PetscViewerSetFromOptions(viewer);CHKERRQ(ierr);}
1502cc843e7aSLisandro Dalcin   if (!vbinary->filename) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call PetscViewerFileSetName()");
1503cc843e7aSLisandro Dalcin   if (vbinary->filemode == (PetscFileMode)-1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call PetscViewerFileSetMode()");
1504cc843e7aSLisandro Dalcin   ierr = PetscViewerFileClose_Binary(viewer);CHKERRQ(ierr);
1505cc843e7aSLisandro Dalcin 
1506cc843e7aSLisandro Dalcin   ierr = PetscViewerBinaryGetUseMPIIO(viewer,&usempiio);CHKERRQ(ierr);
1507cc843e7aSLisandro Dalcin   if (usempiio) {
1508cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
1509cc843e7aSLisandro Dalcin     ierr = PetscViewerFileSetUp_BinaryMPIIO(viewer);CHKERRQ(ierr);
1510cc843e7aSLisandro Dalcin #endif
1511cc843e7aSLisandro Dalcin   } else {
1512cc843e7aSLisandro Dalcin     ierr = PetscViewerFileSetUp_BinarySTDIO(viewer);CHKERRQ(ierr);
1513cc843e7aSLisandro Dalcin   }
1514cc843e7aSLisandro Dalcin   ierr = PetscViewerFileSetUp_BinaryInfo(viewer);CHKERRQ(ierr);
1515cc843e7aSLisandro Dalcin 
1516cc843e7aSLisandro Dalcin   ierr = PetscLogObjectState((PetscObject)viewer,"File: %s",vbinary->filename);CHKERRQ(ierr);
15175c6c1daeSBarry Smith   PetscFunctionReturn(0);
15185c6c1daeSBarry Smith }
15195c6c1daeSBarry Smith 
152081f0254dSBarry Smith static PetscErrorCode PetscViewerView_Binary(PetscViewer v,PetscViewer viewer)
15212bf49c77SBarry Smith {
1522cb6ad94fSLisandro Dalcin   PetscViewer_Binary *vbinary = (PetscViewer_Binary*)v->data;
1523cb6ad94fSLisandro Dalcin   const char         *fname = vbinary->filename ? vbinary->filename : "not yet set";
1524cc843e7aSLisandro Dalcin   const char         *fmode = vbinary->filemode != (PetscFileMode) -1 ? PetscFileModes[vbinary->filemode] : "not yet set";
1525cb6ad94fSLisandro Dalcin   PetscBool          usempiio;
1526cc843e7aSLisandro Dalcin   PetscErrorCode     ierr;
15272bf49c77SBarry Smith 
15282bf49c77SBarry Smith   PetscFunctionBegin;
1529cb6ad94fSLisandro Dalcin   ierr = PetscViewerBinaryGetUseMPIIO(v,&usempiio);CHKERRQ(ierr);
1530cb6ad94fSLisandro Dalcin   ierr = PetscViewerASCIIPrintf(viewer,"Filename: %s\n",fname);CHKERRQ(ierr);
1531cb6ad94fSLisandro Dalcin   ierr = PetscViewerASCIIPrintf(viewer,"Mode: %s (%s)\n",fmode,usempiio ? "mpiio" : "stdio");CHKERRQ(ierr);
15322bf49c77SBarry Smith   PetscFunctionReturn(0);
15332bf49c77SBarry Smith }
15342bf49c77SBarry Smith 
153522a8f86cSLisandro Dalcin static PetscErrorCode PetscViewerSetFromOptions_Binary(PetscOptionItems *PetscOptionsObject,PetscViewer viewer)
153603a1d59fSDave May {
153722a8f86cSLisandro Dalcin   PetscViewer_Binary *binary = (PetscViewer_Binary*)viewer->data;
1538d22fd6bcSDave May   char               defaultname[PETSC_MAX_PATH_LEN];
153903a1d59fSDave May   PetscBool          flg;
1540cc843e7aSLisandro Dalcin   PetscErrorCode     ierr;
1541e0385b85SDave May 
154203a1d59fSDave May   PetscFunctionBegin;
154322a8f86cSLisandro Dalcin   if (viewer->setupcalled) PetscFunctionReturn(0);
154403a1d59fSDave May   ierr = PetscOptionsHead(PetscOptionsObject,"Binary PetscViewer Options");CHKERRQ(ierr);
1545d22fd6bcSDave May   ierr = PetscSNPrintf(defaultname,PETSC_MAX_PATH_LEN-1,"binaryoutput");CHKERRQ(ierr);
1546589a23caSBarry Smith   ierr = PetscOptionsString("-viewer_binary_filename","Specify filename","PetscViewerFileSetName",defaultname,defaultname,sizeof(defaultname),&flg);CHKERRQ(ierr);
154722a8f86cSLisandro Dalcin   if (flg) { ierr = PetscViewerFileSetName_Binary(viewer,defaultname);CHKERRQ(ierr); }
154822a8f86cSLisandro Dalcin   ierr = PetscOptionsBool("-viewer_binary_skip_info","Skip writing/reading .info file","PetscViewerBinarySetSkipInfo",binary->skipinfo,&binary->skipinfo,NULL);CHKERRQ(ierr);
1549cc843e7aSLisandro Dalcin   ierr = PetscOptionsBool("-viewer_binary_skip_options","Skip parsing Vec/Mat load options","PetscViewerBinarySetSkipOptions",binary->skipoptions,&binary->skipoptions,NULL);CHKERRQ(ierr);
155022a8f86cSLisandro Dalcin   ierr = PetscOptionsBool("-viewer_binary_skip_header","Skip writing/reading header information","PetscViewerBinarySetSkipHeader",binary->skipheader,&binary->skipheader,NULL);CHKERRQ(ierr);
155103a1d59fSDave May #if defined(PETSC_HAVE_MPIIO)
155222a8f86cSLisandro Dalcin   ierr = PetscOptionsBool("-viewer_binary_mpiio","Use MPI-IO functionality to write/read binary file","PetscViewerBinarySetUseMPIIO",binary->usempiio,&binary->usempiio,NULL);CHKERRQ(ierr);
15535972f5f3SLisandro Dalcin #else
15545972f5f3SLisandro Dalcin   ierr = PetscOptionsBool("-viewer_binary_mpiio","Use MPI-IO functionality to write/read binary file (NOT AVAILABLE)","PetscViewerBinarySetUseMPIIO",PETSC_FALSE,NULL,NULL);CHKERRQ(ierr);
155503a1d59fSDave May #endif
155603a1d59fSDave May   ierr = PetscOptionsTail();CHKERRQ(ierr);
1557bc196f7cSDave May   binary->setfromoptionscalled = PETSC_TRUE;
155803a1d59fSDave May   PetscFunctionReturn(0);
155903a1d59fSDave May }
156003a1d59fSDave May 
15618556b5ebSBarry Smith /*MC
15628556b5ebSBarry Smith    PETSCVIEWERBINARY - A viewer that saves to binary files
15638556b5ebSBarry Smith 
15648556b5ebSBarry Smith 
15658556b5ebSBarry Smith .seealso:  PetscViewerBinaryOpen(), PETSC_VIEWER_STDOUT_(),PETSC_VIEWER_STDOUT_SELF, PETSC_VIEWER_STDOUT_WORLD, PetscViewerCreate(), PetscViewerASCIIOpen(),
15668556b5ebSBarry Smith            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB, PETSCVIEWERDRAW,
156767918a83SBarry Smith            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType(),
156867918a83SBarry Smith            PetscViewerBinaryGetUseMPIIO(), PetscViewerBinarySetUseMPIIO()
15698556b5ebSBarry Smith 
15701b266c99SBarry Smith   Level: beginner
15711b266c99SBarry Smith 
15728556b5ebSBarry Smith M*/
15738556b5ebSBarry Smith 
15748cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_Binary(PetscViewer v)
15755c6c1daeSBarry Smith {
15765c6c1daeSBarry Smith   PetscErrorCode     ierr;
15775c6c1daeSBarry Smith   PetscViewer_Binary *vbinary;
15785c6c1daeSBarry Smith 
15795c6c1daeSBarry Smith   PetscFunctionBegin;
1580b00a9115SJed Brown   ierr    = PetscNewLog(v,&vbinary);CHKERRQ(ierr);
15815c6c1daeSBarry Smith   v->data = (void*)vbinary;
1582cc843e7aSLisandro Dalcin 
158303a1d59fSDave May   v->ops->setfromoptions   = PetscViewerSetFromOptions_Binary;
15845c6c1daeSBarry Smith   v->ops->destroy          = PetscViewerDestroy_Binary;
15852bf49c77SBarry Smith   v->ops->view             = PetscViewerView_Binary;
1586c98fd787SBarry Smith   v->ops->setup            = PetscViewerSetUp_Binary;
1587cc843e7aSLisandro Dalcin   v->ops->flush            = NULL; /* Should we support Flush() ? */
1588cc843e7aSLisandro Dalcin   v->ops->getsubviewer     = PetscViewerGetSubViewer_Binary;
1589cc843e7aSLisandro Dalcin   v->ops->restoresubviewer = PetscViewerRestoreSubViewer_Binary;
1590cc843e7aSLisandro Dalcin   v->ops->read             = PetscViewerBinaryRead;
1591cc843e7aSLisandro Dalcin 
1592cc843e7aSLisandro Dalcin   vbinary->fdes            = -1;
1593e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO)
1594cc843e7aSLisandro Dalcin   vbinary->usempiio        = PETSC_FALSE;
1595e8a65b7dSLisandro Dalcin   vbinary->mfdes           = MPI_FILE_NULL;
1596e8a65b7dSLisandro Dalcin   vbinary->mfsub           = MPI_FILE_NULL;
1597e8a65b7dSLisandro Dalcin #endif
1598cc843e7aSLisandro Dalcin   vbinary->filename        = NULL;
1599cc843e7aSLisandro Dalcin   vbinary->filemode        = (PetscFileMode)-1;
160002c9f0b5SLisandro Dalcin   vbinary->fdes_info       = NULL;
16015c6c1daeSBarry Smith   vbinary->skipinfo        = PETSC_FALSE;
16025c6c1daeSBarry Smith   vbinary->skipoptions     = PETSC_TRUE;
16035c6c1daeSBarry Smith   vbinary->skipheader      = PETSC_FALSE;
16045c6c1daeSBarry Smith   vbinary->storecompressed = PETSC_FALSE;
1605f90597f1SStefano Zampini   vbinary->ogzfilename     = NULL;
16065c6c1daeSBarry Smith   vbinary->flowcontrol     = 256; /* seems a good number for Cray XT-5 */
16075c6c1daeSBarry Smith 
1608cc843e7aSLisandro Dalcin   vbinary->setfromoptionscalled = PETSC_FALSE;
1609cc843e7aSLisandro Dalcin 
1610bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetFlowControl_C",PetscViewerBinaryGetFlowControl_Binary);CHKERRQ(ierr);
1611bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetFlowControl_C",PetscViewerBinarySetFlowControl_Binary);CHKERRQ(ierr);
1612bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetSkipHeader_C",PetscViewerBinaryGetSkipHeader_Binary);CHKERRQ(ierr);
1613cc843e7aSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetSkipHeader_C",PetscViewerBinarySetSkipHeader_Binary);CHKERRQ(ierr);
1614807ea322SDave May   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetSkipOptions_C",PetscViewerBinaryGetSkipOptions_Binary);CHKERRQ(ierr);
1615807ea322SDave May   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetSkipOptions_C",PetscViewerBinarySetSkipOptions_Binary);CHKERRQ(ierr);
1616807ea322SDave May   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetSkipInfo_C",PetscViewerBinaryGetSkipInfo_Binary);CHKERRQ(ierr);
1617807ea322SDave May   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetSkipInfo_C",PetscViewerBinarySetSkipInfo_Binary);CHKERRQ(ierr);
1618bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetInfoPointer_C",PetscViewerBinaryGetInfoPointer_Binary);CHKERRQ(ierr);
1619bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_Binary);CHKERRQ(ierr);
1620cc843e7aSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_Binary);CHKERRQ(ierr);
1621cc843e7aSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_Binary);CHKERRQ(ierr);
1622cc843e7aSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_Binary);CHKERRQ(ierr);
16235c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO)
1624bc196f7cSDave May   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetUseMPIIO_C",PetscViewerBinaryGetUseMPIIO_Binary);CHKERRQ(ierr);
1625bc196f7cSDave May   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetUseMPIIO_C",PetscViewerBinarySetUseMPIIO_Binary);CHKERRQ(ierr);
16265c6c1daeSBarry Smith #endif
16275c6c1daeSBarry Smith   PetscFunctionReturn(0);
16285c6c1daeSBarry Smith }
16295c6c1daeSBarry Smith 
16305c6c1daeSBarry Smith /* ---------------------------------------------------------------------*/
16315c6c1daeSBarry Smith /*
16325c6c1daeSBarry Smith     The variable Petsc_Viewer_Binary_keyval is used to indicate an MPI attribute that
16335c6c1daeSBarry Smith   is attached to a communicator, in this case the attribute is a PetscViewer.
16345c6c1daeSBarry Smith */
1635d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_Binary_keyval = MPI_KEYVAL_INVALID;
16365c6c1daeSBarry Smith 
16375c6c1daeSBarry Smith /*@C
16385c6c1daeSBarry Smith      PETSC_VIEWER_BINARY_ - Creates a binary PetscViewer shared by all processors
16395c6c1daeSBarry Smith                      in a communicator.
16405c6c1daeSBarry Smith 
1641d083f849SBarry Smith      Collective
16425c6c1daeSBarry Smith 
16435c6c1daeSBarry Smith      Input Parameter:
16445c6c1daeSBarry Smith .    comm - the MPI communicator to share the binary PetscViewer
16455c6c1daeSBarry Smith 
16465c6c1daeSBarry Smith      Level: intermediate
16475c6c1daeSBarry Smith 
16485c6c1daeSBarry Smith    Options Database Keys:
16495c6c1daeSBarry Smith +    -viewer_binary_filename <name>
16505c6c1daeSBarry Smith .    -viewer_binary_skip_info
1651a2d7db39SDave May .    -viewer_binary_skip_options
1652a2d7db39SDave May .    -viewer_binary_skip_header
1653a2d7db39SDave May -    -viewer_binary_mpiio
16545c6c1daeSBarry Smith 
16555c6c1daeSBarry Smith    Environmental variables:
16565c6c1daeSBarry Smith -   PETSC_VIEWER_BINARY_FILENAME
16575c6c1daeSBarry Smith 
16585c6c1daeSBarry Smith      Notes:
16595c6c1daeSBarry Smith      Unlike almost all other PETSc routines, PETSC_VIEWER_BINARY_ does not return
16605c6c1daeSBarry Smith      an error code.  The binary PetscViewer is usually used in the form
16615c6c1daeSBarry Smith $       XXXView(XXX object,PETSC_VIEWER_BINARY_(comm));
16625c6c1daeSBarry Smith 
16635c6c1daeSBarry Smith .seealso: PETSC_VIEWER_BINARY_WORLD, PETSC_VIEWER_BINARY_SELF, PetscViewerBinaryOpen(), PetscViewerCreate(),
16645c6c1daeSBarry Smith           PetscViewerDestroy()
16655c6c1daeSBarry Smith @*/
16665c6c1daeSBarry Smith PetscViewer PETSC_VIEWER_BINARY_(MPI_Comm comm)
16675c6c1daeSBarry Smith {
16685c6c1daeSBarry Smith   PetscErrorCode ierr;
16695c6c1daeSBarry Smith   PetscBool      flg;
16705c6c1daeSBarry Smith   PetscViewer    viewer;
16715c6c1daeSBarry Smith   char           fname[PETSC_MAX_PATH_LEN];
16725c6c1daeSBarry Smith   MPI_Comm       ncomm;
16735c6c1daeSBarry Smith 
16745c6c1daeSBarry Smith   PetscFunctionBegin;
167502c9f0b5SLisandro Dalcin   ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
16765c6c1daeSBarry Smith   if (Petsc_Viewer_Binary_keyval == MPI_KEYVAL_INVALID) {
167702c9f0b5SLisandro Dalcin     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_Binary_keyval,NULL);
167802c9f0b5SLisandro Dalcin     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
16795c6c1daeSBarry Smith   }
168047435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_Binary_keyval,(void**)&viewer,(int*)&flg);
168102c9f0b5SLisandro Dalcin   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
16825c6c1daeSBarry Smith   if (!flg) { /* PetscViewer not yet created */
16835c6c1daeSBarry Smith     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_BINARY_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
16842cb5e1ccSBarry Smith     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
16855c6c1daeSBarry Smith     if (!flg) {
16865c6c1daeSBarry Smith       ierr = PetscStrcpy(fname,"binaryoutput");
16872cb5e1ccSBarry Smith       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
16885c6c1daeSBarry Smith     }
16895c6c1daeSBarry Smith     ierr = PetscViewerBinaryOpen(ncomm,fname,FILE_MODE_WRITE,&viewer);
16902cb5e1ccSBarry Smith     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
16915c6c1daeSBarry Smith     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
16922cb5e1ccSBarry Smith     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
169347435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_Binary_keyval,(void*)viewer);
169402c9f0b5SLisandro Dalcin     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
16955c6c1daeSBarry Smith   }
16965c6c1daeSBarry Smith   ierr = PetscCommDestroy(&ncomm);
16972cb5e1ccSBarry Smith   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
16985c6c1daeSBarry Smith   PetscFunctionReturn(viewer);
16995c6c1daeSBarry Smith }
1700