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 29cc843e7aSLisandro Dalcin PetscFunctionBegin; 30cc843e7aSLisandro Dalcin if (vbinary->filemode == FILE_MODE_READ) PetscFunctionReturn(0); 31cc843e7aSLisandro Dalcin if (vbinary->mfsub != MPI_FILE_NULL) { 329566063dSJacob Faibussowitsch PetscCallMPI(MPI_File_sync(vbinary->mfsub)); 33cc843e7aSLisandro Dalcin } 34cc843e7aSLisandro Dalcin if (vbinary->mfdes != MPI_FILE_NULL) { 359566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer))); 369566063dSJacob Faibussowitsch PetscCallMPI(MPI_File_sync(vbinary->mfdes)); 37cc843e7aSLisandro Dalcin } 38cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 39cc843e7aSLisandro Dalcin } 40cc843e7aSLisandro Dalcin #endif 41cc843e7aSLisandro Dalcin 4281f0254dSBarry Smith static PetscErrorCode PetscViewerGetSubViewer_Binary(PetscViewer viewer,MPI_Comm comm,PetscViewer *outviewer) 435c6c1daeSBarry Smith { 44e8a65b7dSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 45cc843e7aSLisandro Dalcin PetscMPIInt rank; 465c6c1daeSBarry Smith 475c6c1daeSBarry Smith PetscFunctionBegin; 489566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer)); 49e8a65b7dSLisandro Dalcin 50e8a65b7dSLisandro Dalcin /* Return subviewer in process zero */ 519566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank)); 52dd400576SPatrick Sanan if (rank == 0) { 53e5afcf28SBarry Smith PetscMPIInt flg; 54e5afcf28SBarry Smith 559566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(PETSC_COMM_SELF,comm,&flg)); 56cc73adaaSBarry Smith PetscCheck(flg == MPI_IDENT || flg == MPI_CONGRUENT,PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscViewerGetSubViewer() for PETSCVIEWERBINARY requires a singleton MPI_Comm"); 579566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(comm,outviewer)); 589566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(*outviewer,PETSCVIEWERBINARY)); 599566063dSJacob Faibussowitsch PetscCall(PetscMemcpy((*outviewer)->data,vbinary,sizeof(PetscViewer_Binary))); 6012f4c3a9SLisandro Dalcin (*outviewer)->setupcalled = PETSC_TRUE; 6196509d17SLisandro Dalcin } else { 6296509d17SLisandro Dalcin *outviewer = NULL; 6396509d17SLisandro Dalcin } 64e8a65b7dSLisandro Dalcin 65e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 66e8a65b7dSLisandro Dalcin if (vbinary->usempiio && *outviewer) { 67e8a65b7dSLisandro Dalcin PetscViewer_Binary *obinary = (PetscViewer_Binary*)(*outviewer)->data; 68e8a65b7dSLisandro Dalcin /* Parent viewer opens a new MPI file handle on PETSC_COMM_SELF and keeps track of it for future reuse */ 69e8a65b7dSLisandro Dalcin if (vbinary->mfsub == MPI_FILE_NULL) { 70e8a65b7dSLisandro Dalcin int amode; 71cc843e7aSLisandro Dalcin switch (vbinary->filemode) { 72e8a65b7dSLisandro Dalcin case FILE_MODE_READ: amode = MPI_MODE_RDONLY; break; 73e8a65b7dSLisandro Dalcin case FILE_MODE_WRITE: amode = MPI_MODE_WRONLY; break; 7422a8f86cSLisandro Dalcin case FILE_MODE_APPEND: amode = MPI_MODE_WRONLY; break; 7598921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported file mode %s",PetscFileModes[vbinary->filemode]); 76e8a65b7dSLisandro Dalcin } 779566063dSJacob Faibussowitsch PetscCallMPI(MPI_File_open(PETSC_COMM_SELF,vbinary->filename,amode,MPI_INFO_NULL,&vbinary->mfsub)); 78e8a65b7dSLisandro Dalcin } 79e8a65b7dSLisandro Dalcin /* Subviewer gets the MPI file handle on PETSC_COMM_SELF */ 80e8a65b7dSLisandro Dalcin obinary->mfdes = vbinary->mfsub; 81e8a65b7dSLisandro Dalcin obinary->mfsub = MPI_FILE_NULL; 82e8a65b7dSLisandro Dalcin obinary->moff = vbinary->moff; 83e8a65b7dSLisandro Dalcin } 84e8a65b7dSLisandro Dalcin #endif 85cc843e7aSLisandro Dalcin 86cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 879566063dSJacob Faibussowitsch PetscCall(PetscViewerBinarySyncMPIIO(viewer)); 88cc843e7aSLisandro Dalcin #endif 895c6c1daeSBarry Smith PetscFunctionReturn(0); 905c6c1daeSBarry Smith } 915c6c1daeSBarry Smith 9281f0254dSBarry Smith static PetscErrorCode PetscViewerRestoreSubViewer_Binary(PetscViewer viewer,MPI_Comm comm,PetscViewer *outviewer) 935c6c1daeSBarry Smith { 94e8a65b7dSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 95cc843e7aSLisandro Dalcin PetscMPIInt rank; 96e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 97e8a65b7dSLisandro Dalcin MPI_Offset moff = 0; 98e8a65b7dSLisandro Dalcin #endif 995c6c1daeSBarry Smith 1005c6c1daeSBarry Smith PetscFunctionBegin; 1019566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank)); 102*c5853193SPierre Jolivet PetscCheck(rank == 0 || !*outviewer,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Subviewer not obtained from viewer"); 103e8a65b7dSLisandro Dalcin 104e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 105e8a65b7dSLisandro Dalcin if (vbinary->usempiio && *outviewer) { 106e8a65b7dSLisandro Dalcin PetscViewer_Binary *obinary = (PetscViewer_Binary*)(*outviewer)->data; 10708401ef6SPierre Jolivet PetscCheck(obinary->mfdes == vbinary->mfsub,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Subviewer not obtained from viewer"); 1089566063dSJacob Faibussowitsch if (obinary->mfsub != MPI_FILE_NULL) PetscCallMPI(MPI_File_close(&obinary->mfsub)); 109e8a65b7dSLisandro Dalcin moff = obinary->moff; 110e8a65b7dSLisandro Dalcin } 111e8a65b7dSLisandro Dalcin #endif 112e8a65b7dSLisandro Dalcin 113e8a65b7dSLisandro Dalcin if (*outviewer) { 114e8a65b7dSLisandro Dalcin PetscViewer_Binary *obinary = (PetscViewer_Binary*)(*outviewer)->data; 11508401ef6SPierre Jolivet PetscCheck(obinary->fdes == vbinary->fdes,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Subviewer not obtained from viewer"); 1169566063dSJacob Faibussowitsch PetscCall(PetscFree((*outviewer)->data)); 1179566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(outviewer)); 1185c6c1daeSBarry Smith } 119e8a65b7dSLisandro Dalcin 120e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 121e8a65b7dSLisandro Dalcin if (vbinary->usempiio) { 122e8a65b7dSLisandro Dalcin PetscInt64 ioff = (PetscInt64)moff; /* We could use MPI_OFFSET datatype (requires MPI 2.2) */ 1239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&ioff,1,MPIU_INT64,0,PetscObjectComm((PetscObject)viewer))); 124e8a65b7dSLisandro Dalcin vbinary->moff = (MPI_Offset)ioff; 125e8a65b7dSLisandro Dalcin } 126e8a65b7dSLisandro Dalcin #endif 127cc843e7aSLisandro Dalcin 128cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 1299566063dSJacob Faibussowitsch PetscCall(PetscViewerBinarySyncMPIIO(viewer)); 130cc843e7aSLisandro Dalcin #endif 1315c6c1daeSBarry Smith PetscFunctionReturn(0); 1325c6c1daeSBarry Smith } 1335c6c1daeSBarry Smith 1345c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO) 1355c6c1daeSBarry Smith /*@C 13685b8072bSPatrick Sanan PetscViewerBinaryGetMPIIOOffset - Gets the current global offset that should be passed to `MPI_File_set_view()` or `MPI_File_{write|read}_at[_all]()` 1375c6c1daeSBarry Smith 1385c6c1daeSBarry Smith Not Collective 1395c6c1daeSBarry Smith 1405c6c1daeSBarry Smith Input Parameter: 14185b8072bSPatrick Sanan . viewer - PetscViewer context, obtained from `PetscViewerBinaryOpen()` 1425c6c1daeSBarry Smith 1435c6c1daeSBarry Smith Output Parameter: 14422a8f86cSLisandro Dalcin . off - the current global offset 1455c6c1daeSBarry Smith 1465c6c1daeSBarry Smith Level: advanced 1475c6c1daeSBarry Smith 1485c6c1daeSBarry Smith Fortran Note: 1495c6c1daeSBarry Smith This routine is not supported in Fortran. 1505c6c1daeSBarry Smith 15185b8072bSPatrick Sanan Use `PetscViewerBinaryAddMPIIOOffset()` to increase this value after you have written a view. 1525c6c1daeSBarry Smith 153db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryAddMPIIOOffset()` 1545c6c1daeSBarry Smith @*/ 1555c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryGetMPIIOOffset(PetscViewer viewer,MPI_Offset *off) 1565c6c1daeSBarry Smith { 15722a8f86cSLisandro Dalcin PetscViewer_Binary *vbinary; 1585c6c1daeSBarry Smith 1595c6c1daeSBarry Smith PetscFunctionBegin; 16022a8f86cSLisandro Dalcin PetscValidHeaderSpecificType(viewer,PETSC_VIEWER_CLASSID,1,PETSCVIEWERBINARY); 16122a8f86cSLisandro Dalcin PetscValidPointer(off,2); 16222a8f86cSLisandro Dalcin vbinary = (PetscViewer_Binary*)viewer->data; 1635c6c1daeSBarry Smith *off = vbinary->moff; 1645c6c1daeSBarry Smith PetscFunctionReturn(0); 1655c6c1daeSBarry Smith } 1665c6c1daeSBarry Smith 1675c6c1daeSBarry Smith /*@C 16822a8f86cSLisandro Dalcin PetscViewerBinaryAddMPIIOOffset - Adds to the current global offset 1695c6c1daeSBarry Smith 17022a8f86cSLisandro Dalcin Logically Collective 1715c6c1daeSBarry Smith 1725c6c1daeSBarry Smith Input Parameters: 1735c6c1daeSBarry Smith + viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 17422a8f86cSLisandro Dalcin - off - the addition to the global offset 1755c6c1daeSBarry Smith 1765c6c1daeSBarry Smith Level: advanced 1775c6c1daeSBarry Smith 1785c6c1daeSBarry Smith Fortran Note: 1795c6c1daeSBarry Smith This routine is not supported in Fortran. 1805c6c1daeSBarry Smith 18185b8072bSPatrick Sanan Use `PetscViewerBinaryGetMPIIOOffset()` to get the value that you should pass to `MPI_File_set_view()` or `MPI_File_{write|read}_at[_all]()` 1825c6c1daeSBarry Smith 183db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen(),PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()` 1845c6c1daeSBarry Smith @*/ 1855c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryAddMPIIOOffset(PetscViewer viewer,MPI_Offset off) 1865c6c1daeSBarry Smith { 18722a8f86cSLisandro Dalcin PetscViewer_Binary *vbinary; 1885c6c1daeSBarry Smith 1895c6c1daeSBarry Smith PetscFunctionBegin; 19022a8f86cSLisandro Dalcin PetscValidHeaderSpecificType(viewer,PETSC_VIEWER_CLASSID,1,PETSCVIEWERBINARY); 19122a8f86cSLisandro Dalcin PetscValidLogicalCollectiveInt(viewer,(PetscInt)off,2); 19222a8f86cSLisandro Dalcin vbinary = (PetscViewer_Binary*)viewer->data; 1935c6c1daeSBarry Smith vbinary->moff += off; 1945c6c1daeSBarry Smith PetscFunctionReturn(0); 1955c6c1daeSBarry Smith } 1965c6c1daeSBarry Smith 1975c6c1daeSBarry Smith /*@C 1985c6c1daeSBarry Smith PetscViewerBinaryGetMPIIODescriptor - Extracts the MPI IO file descriptor from a PetscViewer. 1995c6c1daeSBarry Smith 2005c6c1daeSBarry Smith Not Collective 2015c6c1daeSBarry Smith 2025c6c1daeSBarry Smith Input Parameter: 2035c6c1daeSBarry Smith . viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 2045c6c1daeSBarry Smith 2055c6c1daeSBarry Smith Output Parameter: 2065c6c1daeSBarry Smith . fdes - file descriptor 2075c6c1daeSBarry Smith 2085c6c1daeSBarry Smith Level: advanced 2095c6c1daeSBarry Smith 2105c6c1daeSBarry Smith Fortran Note: 2115c6c1daeSBarry Smith This routine is not supported in Fortran. 2125c6c1daeSBarry Smith 213db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen(),PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()` 2145c6c1daeSBarry Smith @*/ 2155c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryGetMPIIODescriptor(PetscViewer viewer,MPI_File *fdes) 2165c6c1daeSBarry Smith { 21722a8f86cSLisandro Dalcin PetscViewer_Binary *vbinary; 2185c6c1daeSBarry Smith 2195c6c1daeSBarry Smith PetscFunctionBegin; 22022a8f86cSLisandro Dalcin PetscValidHeaderSpecificType(viewer,PETSC_VIEWER_CLASSID,1,PETSCVIEWERBINARY); 22122a8f86cSLisandro Dalcin PetscValidPointer(fdes,2); 2229566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer)); 22322a8f86cSLisandro Dalcin vbinary = (PetscViewer_Binary*)viewer->data; 2245c6c1daeSBarry Smith *fdes = vbinary->mfdes; 2255c6c1daeSBarry Smith PetscFunctionReturn(0); 2265c6c1daeSBarry Smith } 227cc843e7aSLisandro Dalcin #endif 2285c6c1daeSBarry Smith 229cc843e7aSLisandro Dalcin /*@ 230cc843e7aSLisandro Dalcin PetscViewerBinarySetUseMPIIO - Sets a binary viewer to use MPI-IO for reading/writing. Must be called 231cc843e7aSLisandro Dalcin before PetscViewerFileSetName() 232cc843e7aSLisandro Dalcin 233cc843e7aSLisandro Dalcin Logically Collective on PetscViewer 234cc843e7aSLisandro Dalcin 235cc843e7aSLisandro Dalcin Input Parameters: 236cc843e7aSLisandro Dalcin + viewer - the PetscViewer; must be a binary 237cc843e7aSLisandro Dalcin - use - PETSC_TRUE means MPI-IO will be used 238cc843e7aSLisandro Dalcin 239cc843e7aSLisandro Dalcin Options Database: 240cc843e7aSLisandro Dalcin -viewer_binary_mpiio : Flag for using MPI-IO 241cc843e7aSLisandro Dalcin 242cc843e7aSLisandro Dalcin Level: advanced 243cc843e7aSLisandro Dalcin 244db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`, 245db781477SPatrick Sanan `PetscViewerBinaryGetUseMPIIO()` 246cc843e7aSLisandro Dalcin 247cc843e7aSLisandro Dalcin @*/ 248cc843e7aSLisandro Dalcin PetscErrorCode PetscViewerBinarySetUseMPIIO(PetscViewer viewer,PetscBool use) 249a8aae444SDave May { 250a8aae444SDave May PetscFunctionBegin; 251cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 252cc843e7aSLisandro Dalcin PetscValidLogicalCollectiveBool(viewer,use,2); 253cac4c232SBarry Smith PetscTryMethod(viewer,"PetscViewerBinarySetUseMPIIO_C",(PetscViewer,PetscBool),(viewer,use)); 254cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 255cc843e7aSLisandro Dalcin } 256cc843e7aSLisandro Dalcin 257cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 258cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinarySetUseMPIIO_Binary(PetscViewer viewer,PetscBool use) 259cc843e7aSLisandro Dalcin { 260cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 261cc843e7aSLisandro Dalcin PetscFunctionBegin; 262cc73adaaSBarry Smith PetscCheck(!viewer->setupcalled || vbinary->usempiio == use,PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER,"Cannot change MPIIO to %s after setup",PetscBools[use]); 263cc843e7aSLisandro Dalcin vbinary->usempiio = use; 264a8aae444SDave May PetscFunctionReturn(0); 265a8aae444SDave May } 266a8aae444SDave May #endif 267a8aae444SDave May 268cc843e7aSLisandro Dalcin /*@ 269bc196f7cSDave May PetscViewerBinaryGetUseMPIIO - Returns PETSC_TRUE if the binary viewer uses MPI-IO. 2705c6c1daeSBarry Smith 2715c6c1daeSBarry Smith Not Collective 2725c6c1daeSBarry Smith 2735c6c1daeSBarry Smith Input Parameter: 2745c6c1daeSBarry Smith . viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 2755c6c1daeSBarry Smith 2765c6c1daeSBarry Smith Output Parameter: 277d8d19677SJose E. Roman . use - PETSC_TRUE if MPI-IO is being used 2785c6c1daeSBarry Smith 2795c6c1daeSBarry Smith Options Database: 2805c6c1daeSBarry Smith -viewer_binary_mpiio : Flag for using MPI-IO 2815c6c1daeSBarry Smith 2825c6c1daeSBarry Smith Level: advanced 2835c6c1daeSBarry Smith 284bc196f7cSDave May Note: 285bc196f7cSDave May If MPI-IO is not available, this function will always return PETSC_FALSE 286bc196f7cSDave May 2875c6c1daeSBarry Smith Fortran Note: 2885c6c1daeSBarry Smith This routine is not supported in Fortran. 2895c6c1daeSBarry Smith 290db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()` 2915c6c1daeSBarry Smith @*/ 292cc843e7aSLisandro Dalcin PetscErrorCode PetscViewerBinaryGetUseMPIIO(PetscViewer viewer,PetscBool *use) 2935c6c1daeSBarry Smith { 2945c6c1daeSBarry Smith PetscFunctionBegin; 295cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 296cc843e7aSLisandro Dalcin PetscValidBoolPointer(use,2); 297cc843e7aSLisandro Dalcin *use = PETSC_FALSE; 298cac4c232SBarry Smith PetscTryMethod(viewer,"PetscViewerBinaryGetUseMPIIO_C",(PetscViewer,PetscBool*),(viewer,use)); 2995c6c1daeSBarry Smith PetscFunctionReturn(0); 3005c6c1daeSBarry Smith } 3015c6c1daeSBarry Smith 302cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 303cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinaryGetUseMPIIO_Binary(PetscViewer viewer,PetscBool *use) 3045c6c1daeSBarry Smith { 3055c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 3065c6c1daeSBarry Smith 3075c6c1daeSBarry Smith PetscFunctionBegin; 308cc843e7aSLisandro Dalcin *use = vbinary->usempiio; 309cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 310cc843e7aSLisandro Dalcin } 311cc843e7aSLisandro Dalcin #endif 312cc843e7aSLisandro Dalcin 313cc843e7aSLisandro Dalcin /*@ 314cc843e7aSLisandro Dalcin PetscViewerBinarySetFlowControl - Sets how many messages are allowed to outstanding at the same time during parallel IO reads/writes 315cc843e7aSLisandro Dalcin 316cc843e7aSLisandro Dalcin Not Collective 317cc843e7aSLisandro Dalcin 318d8d19677SJose E. Roman Input Parameters: 319cc843e7aSLisandro Dalcin + viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 320cc843e7aSLisandro Dalcin - fc - the number of messages, defaults to 256 if this function was not called 321cc843e7aSLisandro Dalcin 322cc843e7aSLisandro Dalcin Level: advanced 323cc843e7aSLisandro Dalcin 324db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen(),PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinaryGetFlowControl()` 325cc843e7aSLisandro Dalcin 326cc843e7aSLisandro Dalcin @*/ 327cc843e7aSLisandro Dalcin PetscErrorCode PetscViewerBinarySetFlowControl(PetscViewer viewer,PetscInt fc) 328cc843e7aSLisandro Dalcin { 329cc843e7aSLisandro Dalcin PetscFunctionBegin; 330cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 331cc843e7aSLisandro Dalcin PetscValidLogicalCollectiveInt(viewer,fc,2); 332cac4c232SBarry Smith PetscTryMethod(viewer,"PetscViewerBinarySetFlowControl_C",(PetscViewer,PetscInt),(viewer,fc)); 3335c6c1daeSBarry Smith PetscFunctionReturn(0); 3345c6c1daeSBarry Smith } 3355c6c1daeSBarry Smith 336cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinarySetFlowControl_Binary(PetscViewer viewer,PetscInt fc) 337cc843e7aSLisandro Dalcin { 338cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 339cc843e7aSLisandro Dalcin 340cc843e7aSLisandro Dalcin PetscFunctionBegin; 34108401ef6SPierre Jolivet PetscCheck(fc > 1,PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_OUTOFRANGE,"Flow control count must be greater than 1, %" PetscInt_FMT " was set",fc); 342cc843e7aSLisandro Dalcin vbinary->flowcontrol = fc; 343cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 344cc843e7aSLisandro Dalcin } 345cc843e7aSLisandro Dalcin 346cc843e7aSLisandro Dalcin /*@ 3475c6c1daeSBarry Smith PetscViewerBinaryGetFlowControl - Returns how many messages are allowed to outstanding at the same time during parallel IO reads/writes 3485c6c1daeSBarry Smith 3495c6c1daeSBarry Smith Not Collective 3505c6c1daeSBarry Smith 3515c6c1daeSBarry Smith Input Parameter: 3525c6c1daeSBarry Smith . viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 3535c6c1daeSBarry Smith 3545c6c1daeSBarry Smith Output Parameter: 3555c6c1daeSBarry Smith . fc - the number of messages 3565c6c1daeSBarry Smith 3575c6c1daeSBarry Smith Level: advanced 3585c6c1daeSBarry Smith 359db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen(),PetscViewerBinaryGetInfoPointer()`, `PetscViewerBinarySetFlowControl()` 3605c6c1daeSBarry Smith 3615c6c1daeSBarry Smith @*/ 3625c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryGetFlowControl(PetscViewer viewer,PetscInt *fc) 3635c6c1daeSBarry Smith { 3645c6c1daeSBarry Smith PetscFunctionBegin; 365cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 366cc843e7aSLisandro Dalcin PetscValidIntPointer(fc,2); 367cac4c232SBarry Smith PetscUseMethod(viewer,"PetscViewerBinaryGetFlowControl_C",(PetscViewer,PetscInt*),(viewer,fc)); 3685c6c1daeSBarry Smith PetscFunctionReturn(0); 3695c6c1daeSBarry Smith } 3705c6c1daeSBarry Smith 371cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinaryGetFlowControl_Binary(PetscViewer viewer,PetscInt *fc) 3725c6c1daeSBarry Smith { 3735c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 3745c6c1daeSBarry Smith 3755c6c1daeSBarry Smith PetscFunctionBegin; 376cc843e7aSLisandro Dalcin *fc = vbinary->flowcontrol; 3775c6c1daeSBarry Smith PetscFunctionReturn(0); 3785c6c1daeSBarry Smith } 3795c6c1daeSBarry Smith 3805c6c1daeSBarry Smith /*@C 3815c6c1daeSBarry Smith PetscViewerBinaryGetDescriptor - Extracts the file descriptor from a PetscViewer. 3825c6c1daeSBarry Smith 3835872f025SBarry Smith Collective On PetscViewer 3845c6c1daeSBarry Smith 3855c6c1daeSBarry Smith Input Parameter: 3865c6c1daeSBarry Smith . viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 3875c6c1daeSBarry Smith 3885c6c1daeSBarry Smith Output Parameter: 3895c6c1daeSBarry Smith . fdes - file descriptor 3905c6c1daeSBarry Smith 3915c6c1daeSBarry Smith Level: advanced 3925c6c1daeSBarry Smith 3935c6c1daeSBarry Smith Notes: 3945c6c1daeSBarry Smith For writable binary PetscViewers, the descriptor will only be valid for the 3955c6c1daeSBarry Smith first processor in the communicator that shares the PetscViewer. For readable 3965c6c1daeSBarry Smith files it will only be valid on nodes that have the file. If node 0 does not 3975c6c1daeSBarry Smith have the file it generates an error even if another node does have the file. 3985c6c1daeSBarry Smith 3995c6c1daeSBarry Smith Fortran Note: 4005c6c1daeSBarry Smith This routine is not supported in Fortran. 4015c6c1daeSBarry Smith 40295452b02SPatrick Sanan Developer Notes: 40395452b02SPatrick Sanan This must be called on all processes because Dave May changed 4045872f025SBarry Smith the source code that this may be trigger a PetscViewerSetUp() call if it was not previously triggered. 4055872f025SBarry Smith 406db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen(),PetscViewerBinaryGetInfoPointer()` 4075c6c1daeSBarry Smith @*/ 4085c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryGetDescriptor(PetscViewer viewer,int *fdes) 4095c6c1daeSBarry Smith { 41022a8f86cSLisandro Dalcin PetscViewer_Binary *vbinary; 4115c6c1daeSBarry Smith 4125c6c1daeSBarry Smith PetscFunctionBegin; 41322a8f86cSLisandro Dalcin PetscValidHeaderSpecificType(viewer,PETSC_VIEWER_CLASSID,1,PETSCVIEWERBINARY); 41422a8f86cSLisandro Dalcin PetscValidPointer(fdes,2); 4159566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer)); 41622a8f86cSLisandro Dalcin vbinary = (PetscViewer_Binary*)viewer->data; 4175c6c1daeSBarry Smith *fdes = vbinary->fdes; 4185c6c1daeSBarry Smith PetscFunctionReturn(0); 4195c6c1daeSBarry Smith } 4205c6c1daeSBarry Smith 4215c6c1daeSBarry Smith /*@ 4225c6c1daeSBarry Smith PetscViewerBinarySkipInfo - Binary file will not have .info file created with it 4235c6c1daeSBarry Smith 4245c6c1daeSBarry Smith Not Collective 4255c6c1daeSBarry Smith 426fd292e60Sprj- Input Parameter: 4275c6c1daeSBarry Smith . viewer - PetscViewer context, obtained from PetscViewerCreate() 4285c6c1daeSBarry Smith 4295c6c1daeSBarry Smith Options Database Key: 43010699b91SBarry Smith . -viewer_binary_skip_info - true indicates do not generate .info file 4315c6c1daeSBarry Smith 4325c6c1daeSBarry Smith Level: advanced 4335c6c1daeSBarry Smith 43495452b02SPatrick Sanan Notes: 43595452b02SPatrick Sanan This must be called after PetscViewerSetType(). If you use PetscViewerBinaryOpen() then 4365c6c1daeSBarry Smith you can only skip the info file with the -viewer_binary_skip_info flag. To use the function you must open the 437a2d7db39SDave May viewer with PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinarySkipInfo(). 4385c6c1daeSBarry Smith 4395c6c1daeSBarry Smith The .info contains meta information about the data in the binary file, for example the block size if it was 4405c6c1daeSBarry Smith set for a vector or matrix. 4415c6c1daeSBarry Smith 442db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySetSkipOptions()`, 443db781477SPatrick Sanan `PetscViewerBinaryGetSkipOptions()`, `PetscViewerBinaryGetSkipInfo()` 4445c6c1daeSBarry Smith @*/ 4455c6c1daeSBarry Smith PetscErrorCode PetscViewerBinarySkipInfo(PetscViewer viewer) 4465c6c1daeSBarry Smith { 4475c6c1daeSBarry Smith PetscFunctionBegin; 4489566063dSJacob Faibussowitsch PetscCall(PetscViewerBinarySetSkipInfo(viewer,PETSC_TRUE)); 449807ea322SDave May PetscFunctionReturn(0); 450807ea322SDave May } 451807ea322SDave May 452807ea322SDave May /*@ 453807ea322SDave May PetscViewerBinarySetSkipInfo - Binary file will not have .info file created with it 454807ea322SDave May 455807ea322SDave May Not Collective 456807ea322SDave May 457d8d19677SJose E. Roman Input Parameters: 458cc843e7aSLisandro Dalcin + viewer - PetscViewer context, obtained from PetscViewerCreate() 459cc843e7aSLisandro Dalcin - skip - PETSC_TRUE implies the .info file will not be generated 460807ea322SDave May 461807ea322SDave May Options Database Key: 46210699b91SBarry Smith . -viewer_binary_skip_info - true indicates do not generate .info file 463807ea322SDave May 464807ea322SDave May Level: advanced 465807ea322SDave May 466db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySetSkipOptions()`, 467db781477SPatrick Sanan `PetscViewerBinaryGetSkipOptions()`, `PetscViewerBinaryGetSkipInfo()` 468807ea322SDave May @*/ 469807ea322SDave May PetscErrorCode PetscViewerBinarySetSkipInfo(PetscViewer viewer,PetscBool skip) 470807ea322SDave May { 471807ea322SDave May PetscFunctionBegin; 472cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 473cc843e7aSLisandro Dalcin PetscValidLogicalCollectiveBool(viewer,skip,2); 474cac4c232SBarry Smith PetscTryMethod(viewer,"PetscViewerBinarySetSkipInfo_C",(PetscViewer,PetscBool),(viewer,skip)); 475807ea322SDave May PetscFunctionReturn(0); 476807ea322SDave May } 477807ea322SDave May 478cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinarySetSkipInfo_Binary(PetscViewer viewer,PetscBool skip) 479807ea322SDave May { 480807ea322SDave May PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 481807ea322SDave May 482807ea322SDave May PetscFunctionBegin; 483cc843e7aSLisandro Dalcin vbinary->skipinfo = skip; 484807ea322SDave May PetscFunctionReturn(0); 485807ea322SDave May } 486807ea322SDave May 487807ea322SDave May /*@ 488807ea322SDave May PetscViewerBinaryGetSkipInfo - check if viewer wrote a .info file 489807ea322SDave May 490807ea322SDave May Not Collective 491807ea322SDave May 492807ea322SDave May Input Parameter: 493807ea322SDave May . viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 494807ea322SDave May 495807ea322SDave May Output Parameter: 496807ea322SDave May . skip - PETSC_TRUE implies the .info file was not generated 497807ea322SDave May 498807ea322SDave May Level: advanced 499807ea322SDave May 50095452b02SPatrick Sanan Notes: 50195452b02SPatrick Sanan This must be called after PetscViewerSetType() 502807ea322SDave May 503db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`, 504db781477SPatrick Sanan `PetscViewerBinarySetSkipOptions()`, `PetscViewerBinarySetSkipInfo()` 505807ea322SDave May @*/ 506807ea322SDave May PetscErrorCode PetscViewerBinaryGetSkipInfo(PetscViewer viewer,PetscBool *skip) 507807ea322SDave May { 508807ea322SDave May PetscFunctionBegin; 509cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 510cc843e7aSLisandro Dalcin PetscValidBoolPointer(skip,2); 511cac4c232SBarry Smith PetscUseMethod(viewer,"PetscViewerBinaryGetSkipInfo_C",(PetscViewer,PetscBool*),(viewer,skip)); 512807ea322SDave May PetscFunctionReturn(0); 513807ea322SDave May } 514807ea322SDave May 515cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinaryGetSkipInfo_Binary(PetscViewer viewer,PetscBool *skip) 516807ea322SDave May { 517807ea322SDave May PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 518807ea322SDave May 519807ea322SDave May PetscFunctionBegin; 520cc843e7aSLisandro Dalcin *skip = vbinary->skipinfo; 521807ea322SDave May PetscFunctionReturn(0); 522807ea322SDave May } 523807ea322SDave May 5245c6c1daeSBarry Smith /*@ 5255c6c1daeSBarry Smith PetscViewerBinarySetSkipOptions - do not use the PETSc options database when loading objects 5265c6c1daeSBarry Smith 5275c6c1daeSBarry Smith Not Collective 5285c6c1daeSBarry Smith 5295c6c1daeSBarry Smith Input Parameters: 5305c6c1daeSBarry Smith + viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 53110699b91SBarry Smith - skip - PETSC_TRUE means do not use the options from the options database 5325c6c1daeSBarry Smith 5335c6c1daeSBarry Smith Options Database Key: 53410699b91SBarry Smith . -viewer_binary_skip_options - true means do not use the options from the options database 5355c6c1daeSBarry Smith 5365c6c1daeSBarry Smith Level: advanced 5375c6c1daeSBarry Smith 53895452b02SPatrick Sanan Notes: 53995452b02SPatrick Sanan This must be called after PetscViewerSetType() 5405c6c1daeSBarry Smith 541db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`, 542db781477SPatrick Sanan `PetscViewerBinaryGetSkipOptions()` 5435c6c1daeSBarry Smith @*/ 5445c6c1daeSBarry Smith PetscErrorCode PetscViewerBinarySetSkipOptions(PetscViewer viewer,PetscBool skip) 5455c6c1daeSBarry Smith { 5465c6c1daeSBarry Smith PetscFunctionBegin; 547cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 548cc843e7aSLisandro Dalcin PetscValidLogicalCollectiveBool(viewer,skip,2); 549cac4c232SBarry Smith PetscTryMethod(viewer,"PetscViewerBinarySetSkipOptions_C",(PetscViewer,PetscBool),(viewer,skip)); 550807ea322SDave May PetscFunctionReturn(0); 551807ea322SDave May } 552807ea322SDave May 553cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinarySetSkipOptions_Binary(PetscViewer viewer,PetscBool skip) 554807ea322SDave May { 555cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 556807ea322SDave May 557807ea322SDave May PetscFunctionBegin; 558cc843e7aSLisandro Dalcin vbinary->skipoptions = skip; 5595c6c1daeSBarry Smith PetscFunctionReturn(0); 5605c6c1daeSBarry Smith } 5615c6c1daeSBarry Smith 5625c6c1daeSBarry Smith /*@ 5635c6c1daeSBarry Smith PetscViewerBinaryGetSkipOptions - checks if viewer uses the PETSc options database when loading objects 5645c6c1daeSBarry Smith 5655c6c1daeSBarry Smith Not Collective 5665c6c1daeSBarry Smith 5675c6c1daeSBarry Smith Input Parameter: 5685c6c1daeSBarry Smith . viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 5695c6c1daeSBarry Smith 5705c6c1daeSBarry Smith Output Parameter: 5715c6c1daeSBarry Smith . skip - PETSC_TRUE means do not use 5725c6c1daeSBarry Smith 5735c6c1daeSBarry Smith Level: advanced 5745c6c1daeSBarry Smith 57595452b02SPatrick Sanan Notes: 57695452b02SPatrick Sanan This must be called after PetscViewerSetType() 5775c6c1daeSBarry Smith 578db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`, 579db781477SPatrick Sanan `PetscViewerBinarySetSkipOptions()` 5805c6c1daeSBarry Smith @*/ 5815c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryGetSkipOptions(PetscViewer viewer,PetscBool *skip) 5825c6c1daeSBarry Smith { 5835c6c1daeSBarry Smith PetscFunctionBegin; 584cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 585cc843e7aSLisandro Dalcin PetscValidBoolPointer(skip,2); 586cac4c232SBarry Smith PetscUseMethod(viewer,"PetscViewerBinaryGetSkipOptions_C",(PetscViewer,PetscBool*),(viewer,skip)); 5875c6c1daeSBarry Smith PetscFunctionReturn(0); 5885c6c1daeSBarry Smith } 5895c6c1daeSBarry Smith 590cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinaryGetSkipOptions_Binary(PetscViewer viewer,PetscBool *skip) 5915c6c1daeSBarry Smith { 592d21b9a37SPierre Jolivet PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 5935c6c1daeSBarry Smith 5945c6c1daeSBarry Smith PetscFunctionBegin; 595cc843e7aSLisandro Dalcin *skip = vbinary->skipoptions; 5965c6c1daeSBarry Smith PetscFunctionReturn(0); 5975c6c1daeSBarry Smith } 5985c6c1daeSBarry Smith 5995c6c1daeSBarry Smith /*@ 6005c6c1daeSBarry Smith PetscViewerBinarySetSkipHeader - do not write a header with size information on output, just raw data 6015c6c1daeSBarry Smith 6025c6c1daeSBarry Smith Not Collective 6035c6c1daeSBarry Smith 6045c6c1daeSBarry Smith Input Parameters: 6055c6c1daeSBarry Smith + viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 6065c6c1daeSBarry Smith - skip - PETSC_TRUE means do not write header 6075c6c1daeSBarry Smith 6085c6c1daeSBarry Smith Options Database Key: 60910699b91SBarry Smith . -viewer_binary_skip_header - PETSC_TRUE means do not write header 6105c6c1daeSBarry Smith 6115c6c1daeSBarry Smith Level: advanced 6125c6c1daeSBarry Smith 61395452b02SPatrick Sanan Notes: 61495452b02SPatrick Sanan This must be called after PetscViewerSetType() 6155c6c1daeSBarry Smith 61610699b91SBarry Smith Is ignored on anything but a binary viewer 6175c6c1daeSBarry Smith 618db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`, 619db781477SPatrick Sanan `PetscViewerBinaryGetSkipHeader()` 6205c6c1daeSBarry Smith @*/ 6215c6c1daeSBarry Smith PetscErrorCode PetscViewerBinarySetSkipHeader(PetscViewer viewer,PetscBool skip) 6225c6c1daeSBarry Smith { 6235c6c1daeSBarry Smith PetscFunctionBegin; 624cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 625cc843e7aSLisandro Dalcin PetscValidLogicalCollectiveBool(viewer,skip,2); 626cac4c232SBarry Smith PetscTryMethod(viewer,"PetscViewerBinarySetSkipHeader_C",(PetscViewer,PetscBool),(viewer,skip)); 6275c6c1daeSBarry Smith PetscFunctionReturn(0); 6285c6c1daeSBarry Smith } 6295c6c1daeSBarry Smith 630cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinarySetSkipHeader_Binary(PetscViewer viewer,PetscBool skip) 6315c6c1daeSBarry Smith { 6325c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 6335c6c1daeSBarry Smith 6345c6c1daeSBarry Smith PetscFunctionBegin; 635cc843e7aSLisandro Dalcin vbinary->skipheader = skip; 6365c6c1daeSBarry Smith PetscFunctionReturn(0); 6375c6c1daeSBarry Smith } 6385c6c1daeSBarry Smith 6395c6c1daeSBarry Smith /*@ 6405c6c1daeSBarry Smith PetscViewerBinaryGetSkipHeader - checks whether to write a header with size information on output, or just raw data 6415c6c1daeSBarry Smith 6425c6c1daeSBarry Smith Not Collective 6435c6c1daeSBarry Smith 6445c6c1daeSBarry Smith Input Parameter: 6455c6c1daeSBarry Smith . viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 6465c6c1daeSBarry Smith 6475c6c1daeSBarry Smith Output Parameter: 6485c6c1daeSBarry Smith . skip - PETSC_TRUE means do not write header 6495c6c1daeSBarry Smith 6505c6c1daeSBarry Smith Level: advanced 6515c6c1daeSBarry Smith 65295452b02SPatrick Sanan Notes: 65395452b02SPatrick Sanan This must be called after PetscViewerSetType() 6545c6c1daeSBarry Smith 6555c6c1daeSBarry Smith Returns false for PETSCSOCKETVIEWER, you cannot skip the header for it. 6565c6c1daeSBarry Smith 657db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinaryGetDescriptor()`, `PetscViewerBinarySkipInfo()`, 658db781477SPatrick Sanan `PetscViewerBinarySetSkipHeader()` 6595c6c1daeSBarry Smith @*/ 6605c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryGetSkipHeader(PetscViewer viewer,PetscBool *skip) 6615c6c1daeSBarry Smith { 6625c6c1daeSBarry Smith PetscFunctionBegin; 663cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 664cc843e7aSLisandro Dalcin PetscValidBoolPointer(skip,2); 665cac4c232SBarry Smith PetscUseMethod(viewer,"PetscViewerBinaryGetSkipHeader_C",(PetscViewer,PetscBool*),(viewer,skip)); 6665c6c1daeSBarry Smith PetscFunctionReturn(0); 6675c6c1daeSBarry Smith } 6685c6c1daeSBarry Smith 669cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinaryGetSkipHeader_Binary(PetscViewer viewer,PetscBool *skip) 6705c6c1daeSBarry Smith { 671f3b211e4SSatish Balay PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 6725c6c1daeSBarry Smith 6735c6c1daeSBarry Smith PetscFunctionBegin; 674cc843e7aSLisandro Dalcin *skip = vbinary->skipheader; 6755c6c1daeSBarry Smith PetscFunctionReturn(0); 6765c6c1daeSBarry Smith } 6775c6c1daeSBarry Smith 6785c6c1daeSBarry Smith /*@C 6795c6c1daeSBarry Smith PetscViewerBinaryGetInfoPointer - Extracts the file pointer for the ASCII 6805c6c1daeSBarry Smith info file associated with a binary file. 6815c6c1daeSBarry Smith 6825c6c1daeSBarry Smith Not Collective 6835c6c1daeSBarry Smith 6845c6c1daeSBarry Smith Input Parameter: 6855c6c1daeSBarry Smith . viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 6865c6c1daeSBarry Smith 6875c6c1daeSBarry Smith Output Parameter: 6880298fd71SBarry Smith . file - file pointer Always returns NULL if not a binary viewer 6895c6c1daeSBarry Smith 6905c6c1daeSBarry Smith Level: advanced 6915c6c1daeSBarry Smith 6925c6c1daeSBarry Smith Notes: 6935c6c1daeSBarry Smith For writable binary PetscViewers, the descriptor will only be valid for the 6945c6c1daeSBarry Smith first processor in the communicator that shares the PetscViewer. 6955c6c1daeSBarry Smith 6965c6c1daeSBarry Smith Fortran Note: 6975c6c1daeSBarry Smith This routine is not supported in Fortran. 6985c6c1daeSBarry Smith 699db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen(),PetscViewerBinaryGetDescriptor()` 7005c6c1daeSBarry Smith @*/ 7015c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryGetInfoPointer(PetscViewer viewer,FILE **file) 7025c6c1daeSBarry Smith { 7035c6c1daeSBarry Smith PetscFunctionBegin; 704cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 705cc843e7aSLisandro Dalcin PetscValidPointer(file,2); 7060298fd71SBarry Smith *file = NULL; 707cac4c232SBarry Smith PetscTryMethod(viewer,"PetscViewerBinaryGetInfoPointer_C",(PetscViewer,FILE **),(viewer,file)); 7085c6c1daeSBarry Smith PetscFunctionReturn(0); 7095c6c1daeSBarry Smith } 7105c6c1daeSBarry Smith 711cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinaryGetInfoPointer_Binary(PetscViewer viewer,FILE **file) 7125c6c1daeSBarry Smith { 713cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 7145c6c1daeSBarry Smith 7155c6c1daeSBarry Smith PetscFunctionBegin; 7169566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer)); 717cc843e7aSLisandro Dalcin *file = vbinary->fdes_info; 718cc843e7aSLisandro Dalcin if (viewer->format == PETSC_VIEWER_BINARY_MATLAB && !vbinary->matlabheaderwritten) { 7195c6c1daeSBarry Smith if (vbinary->fdes_info) { 720cc843e7aSLisandro Dalcin FILE *info = vbinary->fdes_info; 7219566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,info,"#--- begin code written by PetscViewerBinary for MATLAB format ---#\n")); 7229566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,info,"#$$ Set.filename = '%s';\n",vbinary->filename)); 7239566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,info,"#$$ fd = PetscOpenFile(Set.filename);\n")); 7249566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,info,"#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n")); 725cc843e7aSLisandro Dalcin } 726cc843e7aSLisandro Dalcin vbinary->matlabheaderwritten = PETSC_TRUE; 7275c6c1daeSBarry Smith } 7285c6c1daeSBarry Smith PetscFunctionReturn(0); 7295c6c1daeSBarry Smith } 7305c6c1daeSBarry Smith 731e0385b85SDave May #if defined(PETSC_HAVE_MPIIO) 732e0385b85SDave May static PetscErrorCode PetscViewerFileClose_BinaryMPIIO(PetscViewer v) 733e0385b85SDave May { 734e0385b85SDave May PetscViewer_Binary *vbinary = (PetscViewer_Binary*)v->data; 735e0385b85SDave May 736e0385b85SDave May PetscFunctionBegin; 737e8a65b7dSLisandro Dalcin if (vbinary->mfdes != MPI_FILE_NULL) { 7389566063dSJacob Faibussowitsch PetscCallMPI(MPI_File_close(&vbinary->mfdes)); 739e0385b85SDave May } 740e8a65b7dSLisandro Dalcin if (vbinary->mfsub != MPI_FILE_NULL) { 7419566063dSJacob Faibussowitsch PetscCallMPI(MPI_File_close(&vbinary->mfsub)); 742e8a65b7dSLisandro Dalcin } 743cc843e7aSLisandro Dalcin vbinary->moff = 0; 744e0385b85SDave May PetscFunctionReturn(0); 745e0385b85SDave May } 746e0385b85SDave May #endif 747e0385b85SDave May 748cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileClose_BinarySTDIO(PetscViewer v) 749cc843e7aSLisandro Dalcin { 750cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)v->data; 751cc843e7aSLisandro Dalcin 752cc843e7aSLisandro Dalcin PetscFunctionBegin; 753cc843e7aSLisandro Dalcin if (vbinary->fdes != -1) { 7549566063dSJacob Faibussowitsch PetscCall(PetscBinaryClose(vbinary->fdes)); 755cc843e7aSLisandro Dalcin vbinary->fdes = -1; 756cc843e7aSLisandro Dalcin if (vbinary->storecompressed) { 757cc843e7aSLisandro Dalcin char cmd[8+PETSC_MAX_PATH_LEN],out[64+PETSC_MAX_PATH_LEN] = ""; 758cc843e7aSLisandro Dalcin const char *gzfilename = vbinary->ogzfilename ? vbinary->ogzfilename : vbinary->filename; 759cc843e7aSLisandro Dalcin /* compress the file */ 7609566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(cmd,"gzip -f ",sizeof(cmd))); 7619566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(cmd,gzfilename,sizeof(cmd))); 762cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_POPEN) 763cc843e7aSLisandro Dalcin { 764cc843e7aSLisandro Dalcin FILE *fp; 7659566063dSJacob Faibussowitsch PetscCall(PetscPOpen(PETSC_COMM_SELF,NULL,cmd,"r",&fp)); 766cc73adaaSBarry Smith PetscCheck(!fgets(out,(int)(sizeof(out)-1),fp),PETSC_COMM_SELF,PETSC_ERR_LIB,"Error from command %s\n%s",cmd,out); 7679566063dSJacob Faibussowitsch PetscCall(PetscPClose(PETSC_COMM_SELF,fp)); 768cc843e7aSLisandro Dalcin } 769cc843e7aSLisandro Dalcin #endif 770cc843e7aSLisandro Dalcin } 771cc843e7aSLisandro Dalcin } 7729566063dSJacob Faibussowitsch PetscCall(PetscFree(vbinary->ogzfilename)); 773cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 774cc843e7aSLisandro Dalcin } 775cc843e7aSLisandro Dalcin 776cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileClose_BinaryInfo(PetscViewer v) 777cc843e7aSLisandro Dalcin { 778cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)v->data; 779cc843e7aSLisandro Dalcin 780cc843e7aSLisandro Dalcin PetscFunctionBegin; 781cc843e7aSLisandro Dalcin if (v->format == PETSC_VIEWER_BINARY_MATLAB && vbinary->matlabheaderwritten) { 782cc843e7aSLisandro Dalcin if (vbinary->fdes_info) { 783cc843e7aSLisandro Dalcin FILE *info = vbinary->fdes_info; 7849566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,info,"#--- begin code written by PetscViewerBinary for MATLAB format ---#\n")); 7859566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,info,"#$$ close(fd);\n")); 7869566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,info,"#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n")); 787cc843e7aSLisandro Dalcin } 788cc843e7aSLisandro Dalcin } 789cc843e7aSLisandro Dalcin if (vbinary->fdes_info) { 790cc843e7aSLisandro Dalcin FILE *info = vbinary->fdes_info; 791cc843e7aSLisandro Dalcin vbinary->fdes_info = NULL; 792cc73adaaSBarry Smith PetscCheck(!fclose(info),PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file"); 793cc843e7aSLisandro Dalcin } 794cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 795cc843e7aSLisandro Dalcin } 796cc843e7aSLisandro Dalcin 797cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileClose_Binary(PetscViewer v) 798cc843e7aSLisandro Dalcin { 799cc843e7aSLisandro Dalcin PetscFunctionBegin; 800cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 8019566063dSJacob Faibussowitsch PetscCall(PetscViewerFileClose_BinaryMPIIO(v)); 802cc843e7aSLisandro Dalcin #endif 8039566063dSJacob Faibussowitsch PetscCall(PetscViewerFileClose_BinarySTDIO(v)); 8049566063dSJacob Faibussowitsch PetscCall(PetscViewerFileClose_BinaryInfo(v)); 805cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 806cc843e7aSLisandro Dalcin } 807cc843e7aSLisandro Dalcin 80881f0254dSBarry Smith static PetscErrorCode PetscViewerDestroy_Binary(PetscViewer v) 8095c6c1daeSBarry Smith { 8105c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary*)v->data; 8115c6c1daeSBarry Smith 8125c6c1daeSBarry Smith PetscFunctionBegin; 8139566063dSJacob Faibussowitsch PetscCall(PetscViewerFileClose_Binary(v)); 8149566063dSJacob Faibussowitsch PetscCall(PetscFree(vbinary->filename)); 8159566063dSJacob Faibussowitsch PetscCall(PetscFree(vbinary)); 81622a8f86cSLisandro Dalcin 8179566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetFlowControl_C",NULL)); 8189566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetFlowControl_C",NULL)); 8199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetSkipHeader_C",NULL)); 8209566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetSkipHeader_C",NULL)); 8219566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetSkipOptions_C",NULL)); 8229566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetSkipOptions_C",NULL)); 8239566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetSkipInfo_C",NULL)); 8249566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetSkipInfo_C",NULL)); 8259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetInfoPointer_C",NULL)); 8269566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",NULL)); 8279566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",NULL)); 8289566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",NULL)); 8299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",NULL)); 83022a8f86cSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 8319566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetUseMPIIO_C",NULL)); 8329566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetUseMPIIO_C",NULL)); 83322a8f86cSLisandro Dalcin #endif 834e0385b85SDave May PetscFunctionReturn(0); 835e0385b85SDave May } 8365c6c1daeSBarry Smith 8375c6c1daeSBarry Smith /*@C 8385c6c1daeSBarry Smith PetscViewerBinaryOpen - Opens a file for binary input/output. 8395c6c1daeSBarry Smith 840d083f849SBarry Smith Collective 8415c6c1daeSBarry Smith 8425c6c1daeSBarry Smith Input Parameters: 8435c6c1daeSBarry Smith + comm - MPI communicator 8445c6c1daeSBarry Smith . name - name of file 845cc843e7aSLisandro Dalcin - mode - open mode of file 8465c6c1daeSBarry Smith $ FILE_MODE_WRITE - create new file for binary output 8475c6c1daeSBarry Smith $ FILE_MODE_READ - open existing file for binary input 8485c6c1daeSBarry Smith $ FILE_MODE_APPEND - open existing file for binary output 8495c6c1daeSBarry Smith 8505c6c1daeSBarry Smith Output Parameter: 851cc843e7aSLisandro Dalcin . viewer - PetscViewer for binary input/output to use with the specified file 8525c6c1daeSBarry Smith 8535c6c1daeSBarry Smith Options Database Keys: 85463c55180SPatrick Sanan + -viewer_binary_filename <name> - 85563c55180SPatrick Sanan . -viewer_binary_skip_info - 85663c55180SPatrick Sanan . -viewer_binary_skip_options - 85763c55180SPatrick Sanan . -viewer_binary_skip_header - 85863c55180SPatrick Sanan - -viewer_binary_mpiio - 8595c6c1daeSBarry Smith 8605c6c1daeSBarry Smith Level: beginner 8615c6c1daeSBarry Smith 8625c6c1daeSBarry Smith Note: 8635c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 8645c6c1daeSBarry Smith 8655c6c1daeSBarry Smith For reading files, the filename may begin with ftp:// or http:// and/or 8665c6c1daeSBarry Smith end with .gz; in this case file is brought over and uncompressed. 8675c6c1daeSBarry Smith 8685c6c1daeSBarry Smith For creating files, if the file name ends with .gz it is automatically 8695c6c1daeSBarry Smith compressed when closed. 8705c6c1daeSBarry Smith 871db781477SPatrick Sanan .seealso: `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`, 872db781477SPatrick Sanan `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`, 873db781477SPatrick Sanan `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()`, `PetscViewerBinarySetUseMPIIO()`, 874db781477SPatrick Sanan `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinaryGetMPIIOOffset()` 8755c6c1daeSBarry Smith @*/ 876cc843e7aSLisandro Dalcin PetscErrorCode PetscViewerBinaryOpen(MPI_Comm comm,const char name[],PetscFileMode mode,PetscViewer *viewer) 8775c6c1daeSBarry Smith { 8785c6c1daeSBarry Smith PetscFunctionBegin; 8799566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(comm,viewer)); 8809566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(*viewer,PETSCVIEWERBINARY)); 8819566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(*viewer,mode)); 8829566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(*viewer,name)); 8839566063dSJacob Faibussowitsch PetscCall(PetscViewerSetFromOptions(*viewer)); 8845c6c1daeSBarry Smith PetscFunctionReturn(0); 8855c6c1daeSBarry Smith } 8865c6c1daeSBarry Smith 8875c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO) 888060da220SMatthew G. Knepley static PetscErrorCode PetscViewerBinaryWriteReadMPIIO(PetscViewer viewer,void *data,PetscInt num,PetscInt *count,PetscDataType dtype,PetscBool write) 8895c6c1daeSBarry Smith { 89022a8f86cSLisandro Dalcin MPI_Comm comm = PetscObjectComm((PetscObject)viewer); 8915c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 89222a8f86cSLisandro Dalcin MPI_File mfdes = vbinary->mfdes; 8935c6c1daeSBarry Smith MPI_Datatype mdtype; 89469e10bbaSLisandro Dalcin PetscMPIInt rank,cnt; 8955c6c1daeSBarry Smith MPI_Status status; 8965c6c1daeSBarry Smith MPI_Aint ul,dsize; 8975c6c1daeSBarry Smith 8985c6c1daeSBarry Smith PetscFunctionBegin; 8999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 9009566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(num,&cnt)); 9019566063dSJacob Faibussowitsch PetscCall(PetscDataTypeToMPIDataType(dtype,&mdtype)); 9025c6c1daeSBarry Smith if (write) { 903dd400576SPatrick Sanan if (rank == 0) { 9049566063dSJacob Faibussowitsch PetscCall(MPIU_File_write_at(mfdes,vbinary->moff,data,cnt,mdtype,&status)); 90569e10bbaSLisandro Dalcin } 9065c6c1daeSBarry Smith } else { 907dd400576SPatrick Sanan if (rank == 0) { 9089566063dSJacob Faibussowitsch PetscCall(MPIU_File_read_at(mfdes,vbinary->moff,data,cnt,mdtype,&status)); 9099566063dSJacob Faibussowitsch if (cnt > 0) PetscCallMPI(MPI_Get_count(&status,mdtype,&cnt)); 9105c6c1daeSBarry Smith } 9119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&cnt,1,MPI_INT,0,comm)); 9129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(data,cnt,mdtype,0,comm)); 91369e10bbaSLisandro Dalcin } 9149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_get_extent(mdtype,&ul,&dsize)); 9155c6c1daeSBarry Smith vbinary->moff += dsize*cnt; 9169860990eSLisandro Dalcin if (count) *count = cnt; 9175c6c1daeSBarry Smith PetscFunctionReturn(0); 9185c6c1daeSBarry Smith } 9195c6c1daeSBarry Smith #endif 9205c6c1daeSBarry Smith 9215c6c1daeSBarry Smith /*@C 9225c6c1daeSBarry Smith PetscViewerBinaryRead - Reads from a binary file, all processors get the same result 9235c6c1daeSBarry Smith 924d083f849SBarry Smith Collective 9255c6c1daeSBarry Smith 9265c6c1daeSBarry Smith Input Parameters: 9275c6c1daeSBarry Smith + viewer - the binary viewer 928907376e6SBarry Smith . data - location of the data to be written 929060da220SMatthew G. Knepley . num - number of items of data to read 930907376e6SBarry Smith - dtype - type of data to read 9315c6c1daeSBarry Smith 932f8e4bde8SMatthew G. Knepley Output Parameters: 9335972f5f3SLisandro Dalcin . count - number of items of data actually read, or NULL. 934f8e4bde8SMatthew G. Knepley 9355c6c1daeSBarry Smith Level: beginner 9365c6c1daeSBarry Smith 937db781477SPatrick Sanan .seealso: `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`, 938db781477SPatrick Sanan `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`, 939db781477SPatrick Sanan `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()` 9405c6c1daeSBarry Smith @*/ 941060da220SMatthew G. Knepley PetscErrorCode PetscViewerBinaryRead(PetscViewer viewer,void *data,PetscInt num,PetscInt *count,PetscDataType dtype) 9425c6c1daeSBarry Smith { 94322a8f86cSLisandro Dalcin PetscViewer_Binary *vbinary; 9445c6c1daeSBarry Smith 94522a8f86cSLisandro Dalcin PetscFunctionBegin; 94622a8f86cSLisandro Dalcin PetscValidHeaderSpecificType(viewer,PETSC_VIEWER_CLASSID,1,PETSCVIEWERBINARY); 94722a8f86cSLisandro Dalcin PetscValidLogicalCollectiveInt(viewer,num,3); 9489566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer)); 94922a8f86cSLisandro Dalcin vbinary = (PetscViewer_Binary*)viewer->data; 9505c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO) 951bc196f7cSDave May if (vbinary->usempiio) { 9529566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWriteReadMPIIO(viewer,data,num,count,dtype,PETSC_FALSE)); 9535c6c1daeSBarry Smith } else { 9545c6c1daeSBarry Smith #endif 9559566063dSJacob Faibussowitsch PetscCall(PetscBinarySynchronizedRead(PetscObjectComm((PetscObject)viewer),vbinary->fdes,data,num,count,dtype)); 9565c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO) 9575c6c1daeSBarry Smith } 9585c6c1daeSBarry Smith #endif 9595c6c1daeSBarry Smith PetscFunctionReturn(0); 9605c6c1daeSBarry Smith } 9615c6c1daeSBarry Smith 9625c6c1daeSBarry Smith /*@C 9635c6c1daeSBarry Smith PetscViewerBinaryWrite - writes to a binary file, only from the first process 9645c6c1daeSBarry Smith 965d083f849SBarry Smith Collective 9665c6c1daeSBarry Smith 9675c6c1daeSBarry Smith Input Parameters: 9685c6c1daeSBarry Smith + viewer - the binary viewer 9695c6c1daeSBarry Smith . data - location of data 9705824c9d0SPatrick Sanan . count - number of items of data to write 971f253e43cSLisandro Dalcin - dtype - type of data to write 9725c6c1daeSBarry Smith 9735c6c1daeSBarry Smith Level: beginner 9745c6c1daeSBarry Smith 975db781477SPatrick Sanan .seealso: `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`, 976db781477SPatrick Sanan `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`, `PetscDataType` 977db781477SPatrick Sanan `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()` 9785c6c1daeSBarry Smith @*/ 979f253e43cSLisandro Dalcin PetscErrorCode PetscViewerBinaryWrite(PetscViewer viewer,const void *data,PetscInt count,PetscDataType dtype) 9805c6c1daeSBarry Smith { 98122a8f86cSLisandro Dalcin PetscViewer_Binary *vbinary; 9825c6c1daeSBarry Smith 9835c6c1daeSBarry Smith PetscFunctionBegin; 98422a8f86cSLisandro Dalcin PetscValidHeaderSpecificType(viewer,PETSC_VIEWER_CLASSID,1,PETSCVIEWERBINARY); 98522a8f86cSLisandro Dalcin PetscValidLogicalCollectiveInt(viewer,count,3); 9869566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer)); 98722a8f86cSLisandro Dalcin vbinary = (PetscViewer_Binary*)viewer->data; 9885c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO) 989bc196f7cSDave May if (vbinary->usempiio) { 9909566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWriteReadMPIIO(viewer,(void*)data,count,NULL,dtype,PETSC_TRUE)); 9915c6c1daeSBarry Smith } else { 9925c6c1daeSBarry Smith #endif 9939566063dSJacob Faibussowitsch PetscCall(PetscBinarySynchronizedWrite(PetscObjectComm((PetscObject)viewer),vbinary->fdes,data,count,dtype)); 9945c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO) 9955c6c1daeSBarry Smith } 9965c6c1daeSBarry Smith #endif 9975c6c1daeSBarry Smith PetscFunctionReturn(0); 9985c6c1daeSBarry Smith } 9995c6c1daeSBarry Smith 10005972f5f3SLisandro Dalcin static PetscErrorCode PetscViewerBinaryWriteReadAll(PetscViewer viewer,PetscBool write,void *data,PetscInt count,PetscInt start,PetscInt total,PetscDataType dtype) 10015972f5f3SLisandro Dalcin { 10025972f5f3SLisandro Dalcin MPI_Comm comm = PetscObjectComm((PetscObject)viewer); 10035972f5f3SLisandro Dalcin PetscMPIInt size,rank; 10045972f5f3SLisandro Dalcin MPI_Datatype mdtype; 1005ec4bef21SJose E. Roman PETSC_UNUSED MPI_Aint lb; 1006ec4bef21SJose E. Roman MPI_Aint dsize; 10075972f5f3SLisandro Dalcin PetscBool useMPIIO; 10085972f5f3SLisandro Dalcin 10095972f5f3SLisandro Dalcin PetscFunctionBegin; 10105972f5f3SLisandro Dalcin PetscValidHeaderSpecificType(viewer,PETSC_VIEWER_CLASSID,1,PETSCVIEWERBINARY); 1011064a246eSJacob Faibussowitsch PetscValidLogicalCollectiveBool(viewer,((start>=0)||(start==PETSC_DETERMINE)),5); 1012064a246eSJacob Faibussowitsch PetscValidLogicalCollectiveBool(viewer,((total>=0)||(total==PETSC_DETERMINE)),6); 1013064a246eSJacob Faibussowitsch PetscValidLogicalCollectiveInt(viewer,total,6); 10149566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer)); 10155972f5f3SLisandro Dalcin 10169566063dSJacob Faibussowitsch PetscCall(PetscDataTypeToMPIDataType(dtype,&mdtype)); 10179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_get_extent(mdtype,&lb,&dsize)); 10189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 10199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 10205972f5f3SLisandro Dalcin 10219566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetUseMPIIO(viewer,&useMPIIO)); 10225972f5f3SLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 10235972f5f3SLisandro Dalcin if (useMPIIO) { 10245972f5f3SLisandro Dalcin MPI_File mfdes; 10255972f5f3SLisandro Dalcin MPI_Offset off; 10265972f5f3SLisandro Dalcin PetscMPIInt cnt; 10275972f5f3SLisandro Dalcin 10285972f5f3SLisandro Dalcin if (start == PETSC_DETERMINE) { 10299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&count,&start,1,MPIU_INT,MPI_SUM,comm)); 10305972f5f3SLisandro Dalcin start -= count; 10315972f5f3SLisandro Dalcin } 10325972f5f3SLisandro Dalcin if (total == PETSC_DETERMINE) { 10335972f5f3SLisandro Dalcin total = start + count; 10349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&total,1,MPIU_INT,size-1,comm)); 10355972f5f3SLisandro Dalcin } 10369566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(count,&cnt)); 10379566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes)); 10389566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetMPIIOOffset(viewer,&off)); 10395972f5f3SLisandro Dalcin off += (MPI_Offset)(start*dsize); 10405972f5f3SLisandro Dalcin if (write) { 10419566063dSJacob Faibussowitsch PetscCall(MPIU_File_write_at_all(mfdes,off,data,cnt,mdtype,MPI_STATUS_IGNORE)); 10425972f5f3SLisandro Dalcin } else { 10439566063dSJacob Faibussowitsch PetscCall(MPIU_File_read_at_all(mfdes,off,data,cnt,mdtype,MPI_STATUS_IGNORE)); 10445972f5f3SLisandro Dalcin } 10455972f5f3SLisandro Dalcin off = (MPI_Offset)(total*dsize); 10469566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryAddMPIIOOffset(viewer,off)); 10475972f5f3SLisandro Dalcin PetscFunctionReturn(0); 10485972f5f3SLisandro Dalcin } 10495972f5f3SLisandro Dalcin #endif 10505972f5f3SLisandro Dalcin { 10515972f5f3SLisandro Dalcin int fdes; 10525972f5f3SLisandro Dalcin char *workbuf = NULL; 1053dd400576SPatrick Sanan PetscInt tcount = rank == 0 ? 0 : count,maxcount=0,message_count,flowcontrolcount; 10545972f5f3SLisandro Dalcin PetscMPIInt tag,cnt,maxcnt,scnt=0,rcnt=0,j; 10555972f5f3SLisandro Dalcin MPI_Status status; 10565972f5f3SLisandro Dalcin 10579566063dSJacob Faibussowitsch PetscCall(PetscCommGetNewTag(comm,&tag)); 10589566063dSJacob Faibussowitsch PetscCallMPI(MPI_Reduce(&tcount,&maxcount,1,MPIU_INT,MPI_MAX,0,comm)); 10599566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(maxcount,&maxcnt)); 10609566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(count,&cnt)); 10615972f5f3SLisandro Dalcin 10629566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetDescriptor(viewer,&fdes)); 10639566063dSJacob Faibussowitsch PetscCall(PetscViewerFlowControlStart(viewer,&message_count,&flowcontrolcount)); 1064dd400576SPatrick Sanan if (rank == 0) { 10659566063dSJacob Faibussowitsch PetscCall(PetscMalloc(maxcnt*dsize,&workbuf)); 10665972f5f3SLisandro Dalcin if (write) { 10679566063dSJacob Faibussowitsch PetscCall(PetscBinaryWrite(fdes,data,cnt,dtype)); 10685972f5f3SLisandro Dalcin } else { 10699566063dSJacob Faibussowitsch PetscCall(PetscBinaryRead(fdes,data,cnt,NULL,dtype)); 10705972f5f3SLisandro Dalcin } 10715972f5f3SLisandro Dalcin for (j=1; j<size; j++) { 10729566063dSJacob Faibussowitsch PetscCall(PetscViewerFlowControlStepMain(viewer,j,&message_count,flowcontrolcount)); 10735972f5f3SLisandro Dalcin if (write) { 10749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(workbuf,maxcnt,mdtype,j,tag,comm,&status)); 10759566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status,mdtype,&rcnt)); 10769566063dSJacob Faibussowitsch PetscCall(PetscBinaryWrite(fdes,workbuf,rcnt,dtype)); 10775972f5f3SLisandro Dalcin } else { 10789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(&scnt,1,MPI_INT,j,tag,comm,MPI_STATUS_IGNORE)); 10799566063dSJacob Faibussowitsch PetscCall(PetscBinaryRead(fdes,workbuf,scnt,NULL,dtype)); 10809566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(workbuf,scnt,mdtype,j,tag,comm)); 10815972f5f3SLisandro Dalcin } 10825972f5f3SLisandro Dalcin } 10839566063dSJacob Faibussowitsch PetscCall(PetscFree(workbuf)); 10849566063dSJacob Faibussowitsch PetscCall(PetscViewerFlowControlEndMain(viewer,&message_count)); 10855972f5f3SLisandro Dalcin } else { 10869566063dSJacob Faibussowitsch PetscCall(PetscViewerFlowControlStepWorker(viewer,rank,&message_count)); 10875972f5f3SLisandro Dalcin if (write) { 10889566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(data,cnt,mdtype,0,tag,comm)); 10895972f5f3SLisandro Dalcin } else { 10909566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(&cnt,1,MPI_INT,0,tag,comm)); 10919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(data,cnt,mdtype,0,tag,comm,MPI_STATUS_IGNORE)); 10925972f5f3SLisandro Dalcin } 10939566063dSJacob Faibussowitsch PetscCall(PetscViewerFlowControlEndWorker(viewer,&message_count)); 10945972f5f3SLisandro Dalcin } 10955972f5f3SLisandro Dalcin } 10965972f5f3SLisandro Dalcin PetscFunctionReturn(0); 10975972f5f3SLisandro Dalcin } 10985972f5f3SLisandro Dalcin 10995972f5f3SLisandro Dalcin /*@C 11005972f5f3SLisandro Dalcin PetscViewerBinaryReadAll - reads from a binary file from all processes 11015972f5f3SLisandro Dalcin 11025972f5f3SLisandro Dalcin Collective 11035972f5f3SLisandro Dalcin 11045972f5f3SLisandro Dalcin Input Parameters: 11055972f5f3SLisandro Dalcin + viewer - the binary viewer 11065972f5f3SLisandro Dalcin . data - location of data 11075972f5f3SLisandro Dalcin . count - local number of items of data to read 11085972f5f3SLisandro Dalcin . start - local start, can be PETSC_DETERMINE 11095972f5f3SLisandro Dalcin . total - global number of items of data to read, can be PETSC_DETERMINE 11105972f5f3SLisandro Dalcin - dtype - type of data to read 11115972f5f3SLisandro Dalcin 11125972f5f3SLisandro Dalcin Level: advanced 11135972f5f3SLisandro Dalcin 1114db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryRead()`, `PetscViewerBinaryWriteAll()` 11155972f5f3SLisandro Dalcin @*/ 11165972f5f3SLisandro Dalcin PetscErrorCode PetscViewerBinaryReadAll(PetscViewer viewer,void *data,PetscInt count,PetscInt start,PetscInt total,PetscDataType dtype) 11175972f5f3SLisandro Dalcin { 11185972f5f3SLisandro Dalcin PetscFunctionBegin; 11199566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWriteReadAll(viewer,PETSC_FALSE,data,count,start,total,dtype)); 11205972f5f3SLisandro Dalcin PetscFunctionReturn(0); 11215972f5f3SLisandro Dalcin } 11225972f5f3SLisandro Dalcin 11235972f5f3SLisandro Dalcin /*@C 11245972f5f3SLisandro Dalcin PetscViewerBinaryWriteAll - writes to a binary file from all processes 11255972f5f3SLisandro Dalcin 11265972f5f3SLisandro Dalcin Collective 11275972f5f3SLisandro Dalcin 11285972f5f3SLisandro Dalcin Input Parameters: 11295972f5f3SLisandro Dalcin + viewer - the binary viewer 11305972f5f3SLisandro Dalcin . data - location of data 11315972f5f3SLisandro Dalcin . count - local number of items of data to write 11325972f5f3SLisandro Dalcin . start - local start, can be PETSC_DETERMINE 11335972f5f3SLisandro Dalcin . total - global number of items of data to write, can be PETSC_DETERMINE 11345972f5f3SLisandro Dalcin - dtype - type of data to write 11355972f5f3SLisandro Dalcin 11365972f5f3SLisandro Dalcin Level: advanced 11375972f5f3SLisandro Dalcin 1138db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PetscViewerBinarySetUseMPIIO()`, `PetscViewerBinaryWriteAll()`, `PetscViewerBinaryReadAll()` 11395972f5f3SLisandro Dalcin @*/ 11405972f5f3SLisandro Dalcin PetscErrorCode PetscViewerBinaryWriteAll(PetscViewer viewer,const void *data,PetscInt count,PetscInt start,PetscInt total,PetscDataType dtype) 11415972f5f3SLisandro Dalcin { 11425972f5f3SLisandro Dalcin PetscFunctionBegin; 11439566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWriteReadAll(viewer,PETSC_TRUE,(void*)data,count,start,total,dtype)); 11445972f5f3SLisandro Dalcin PetscFunctionReturn(0); 11455972f5f3SLisandro Dalcin } 11465972f5f3SLisandro Dalcin 11475c6c1daeSBarry Smith /*@C 11485c6c1daeSBarry Smith PetscViewerBinaryWriteStringArray - writes to a binary file, only from the first process an array of strings 11495c6c1daeSBarry Smith 1150d083f849SBarry Smith Collective 11515c6c1daeSBarry Smith 11525c6c1daeSBarry Smith Input Parameters: 11535c6c1daeSBarry Smith + viewer - the binary viewer 11545c6c1daeSBarry Smith - data - location of the array of strings 11555c6c1daeSBarry Smith 11565c6c1daeSBarry Smith Level: intermediate 11575c6c1daeSBarry Smith 115895452b02SPatrick Sanan Notes: 115995452b02SPatrick Sanan array of strings is null terminated 11605c6c1daeSBarry Smith 1161db781477SPatrick Sanan .seealso: `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`, 1162db781477SPatrick Sanan `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`, 1163db781477SPatrick Sanan `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()` 11645c6c1daeSBarry Smith @*/ 116578fbdcc8SBarry Smith PetscErrorCode PetscViewerBinaryWriteStringArray(PetscViewer viewer,const char * const *data) 11665c6c1daeSBarry Smith { 11675c6c1daeSBarry Smith PetscInt i,n = 0,*sizes; 1168cc843e7aSLisandro Dalcin size_t len; 11695c6c1daeSBarry Smith 117022a8f86cSLisandro Dalcin PetscFunctionBegin; 11719566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer)); 11725c6c1daeSBarry Smith /* count number of strings */ 11735c6c1daeSBarry Smith while (data[n++]); 11745c6c1daeSBarry Smith n--; 11759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n+1,&sizes)); 11765c6c1daeSBarry Smith sizes[0] = n; 11775c6c1daeSBarry Smith for (i=0; i<n; i++) { 11789566063dSJacob Faibussowitsch PetscCall(PetscStrlen(data[i],&len)); 1179cc843e7aSLisandro Dalcin sizes[i+1] = (PetscInt)len + 1; /* size includes space for the null terminator */ 11805c6c1daeSBarry Smith } 11819566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWrite(viewer,sizes,n+1,PETSC_INT)); 11825c6c1daeSBarry Smith for (i=0; i<n; i++) { 11839566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWrite(viewer,(void*)data[i],sizes[i+1],PETSC_CHAR)); 11845c6c1daeSBarry Smith } 11859566063dSJacob Faibussowitsch PetscCall(PetscFree(sizes)); 11865c6c1daeSBarry Smith PetscFunctionReturn(0); 11875c6c1daeSBarry Smith } 11885c6c1daeSBarry Smith 11895c6c1daeSBarry Smith /*@C 11905c6c1daeSBarry Smith PetscViewerBinaryReadStringArray - reads a binary file an array of strings 11915c6c1daeSBarry Smith 1192d083f849SBarry Smith Collective 11935c6c1daeSBarry Smith 11945c6c1daeSBarry Smith Input Parameter: 11955c6c1daeSBarry Smith . viewer - the binary viewer 11965c6c1daeSBarry Smith 11975c6c1daeSBarry Smith Output Parameter: 11985c6c1daeSBarry Smith . data - location of the array of strings 11995c6c1daeSBarry Smith 12005c6c1daeSBarry Smith Level: intermediate 12015c6c1daeSBarry Smith 120295452b02SPatrick Sanan Notes: 120395452b02SPatrick Sanan array of strings is null terminated 12045c6c1daeSBarry Smith 1205db781477SPatrick Sanan .seealso: `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`, 1206db781477SPatrick Sanan `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`, `PetscViewerBinaryGetDescriptor()`, 1207db781477SPatrick Sanan `PetscViewerBinaryGetInfoPointer()`, `PetscFileMode`, `PetscViewer`, `PetscViewerBinaryRead()` 12085c6c1daeSBarry Smith @*/ 12095c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryReadStringArray(PetscViewer viewer,char ***data) 12105c6c1daeSBarry Smith { 1211060da220SMatthew G. Knepley PetscInt i,n,*sizes,N = 0; 12125c6c1daeSBarry Smith 121322a8f86cSLisandro Dalcin PetscFunctionBegin; 12149566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer)); 12155c6c1daeSBarry Smith /* count number of strings */ 12169566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer,&n,1,NULL,PETSC_INT)); 12179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&sizes)); 12189566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer,sizes,n,NULL,PETSC_INT)); 1219a297a907SKarl Rupp for (i=0; i<n; i++) N += sizes[i]; 12209566063dSJacob Faibussowitsch PetscCall(PetscMalloc((n+1)*sizeof(char*) + N*sizeof(char),data)); 12215c6c1daeSBarry Smith (*data)[0] = (char*)((*data) + n + 1); 1222a297a907SKarl Rupp for (i=1; i<n; i++) (*data)[i] = (*data)[i-1] + sizes[i-1]; 12239566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer,(*data)[0],N,NULL,PETSC_CHAR)); 122402c9f0b5SLisandro Dalcin (*data)[n] = NULL; 12259566063dSJacob Faibussowitsch PetscCall(PetscFree(sizes)); 12265c6c1daeSBarry Smith PetscFunctionReturn(0); 12275c6c1daeSBarry Smith } 12285c6c1daeSBarry Smith 1229cc843e7aSLisandro Dalcin /*@C 1230cc843e7aSLisandro Dalcin PetscViewerFileSetMode - Sets the open mode of file 1231cc843e7aSLisandro Dalcin 1232cc843e7aSLisandro Dalcin Logically Collective on PetscViewer 1233cc843e7aSLisandro Dalcin 1234cc843e7aSLisandro Dalcin Input Parameters: 1235cc843e7aSLisandro Dalcin + viewer - the PetscViewer; must be a a PETSCVIEWERBINARY, PETSCVIEWERMATLAB, PETSCVIEWERHDF5, or PETSCVIEWERASCII PetscViewer 1236cc843e7aSLisandro Dalcin - mode - open mode of file 1237cc843e7aSLisandro Dalcin $ FILE_MODE_WRITE - create new file for output 1238cc843e7aSLisandro Dalcin $ FILE_MODE_READ - open existing file for input 1239cc843e7aSLisandro Dalcin $ FILE_MODE_APPEND - open existing file for output 1240cc843e7aSLisandro Dalcin 1241cc843e7aSLisandro Dalcin Level: advanced 1242cc843e7aSLisandro Dalcin 1243db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()` 1244cc843e7aSLisandro Dalcin 1245cc843e7aSLisandro Dalcin @*/ 1246cc843e7aSLisandro Dalcin PetscErrorCode PetscViewerFileSetMode(PetscViewer viewer,PetscFileMode mode) 1247cc843e7aSLisandro Dalcin { 1248cc843e7aSLisandro Dalcin PetscFunctionBegin; 1249cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1250cc843e7aSLisandro Dalcin PetscValidLogicalCollectiveEnum(viewer,mode,2); 125108401ef6SPierre Jolivet PetscCheck(mode != FILE_MODE_UNDEFINED,PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Cannot set FILE_MODE_UNDEFINED"); 1252cc73adaaSBarry Smith else PetscCheck(mode >= FILE_MODE_UNDEFINED && mode <= FILE_MODE_APPEND_UPDATE,PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_OUTOFRANGE,"Invalid file mode %d",(int)mode); 1253cac4c232SBarry Smith PetscTryMethod(viewer,"PetscViewerFileSetMode_C",(PetscViewer,PetscFileMode),(viewer,mode)); 1254cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 1255cc843e7aSLisandro Dalcin } 1256cc843e7aSLisandro Dalcin 1257cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileSetMode_Binary(PetscViewer viewer,PetscFileMode mode) 1258cc843e7aSLisandro Dalcin { 1259cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 1260cc843e7aSLisandro Dalcin 1261cc843e7aSLisandro Dalcin PetscFunctionBegin; 1262cc73adaaSBarry Smith PetscCheck(!viewer->setupcalled || vbinary->filemode == mode,PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER,"Cannot change mode to %s after setup",PetscFileModes[mode]); 1263cc843e7aSLisandro Dalcin vbinary->filemode = mode; 1264cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 1265cc843e7aSLisandro Dalcin } 1266cc843e7aSLisandro Dalcin 1267cc843e7aSLisandro Dalcin /*@C 1268cc843e7aSLisandro Dalcin PetscViewerFileGetMode - Gets the open mode of file 1269cc843e7aSLisandro Dalcin 1270cc843e7aSLisandro Dalcin Not Collective 1271cc843e7aSLisandro Dalcin 1272cc843e7aSLisandro Dalcin Input Parameter: 1273cc843e7aSLisandro Dalcin . viewer - the PetscViewer; must be a PETSCVIEWERBINARY, PETSCVIEWERMATLAB, PETSCVIEWERHDF5, or PETSCVIEWERASCII PetscViewer 1274cc843e7aSLisandro Dalcin 1275cc843e7aSLisandro Dalcin Output Parameter: 1276cc843e7aSLisandro Dalcin . mode - open mode of file 1277cc843e7aSLisandro Dalcin $ FILE_MODE_WRITE - create new file for binary output 1278cc843e7aSLisandro Dalcin $ FILE_MODE_READ - open existing file for binary input 1279cc843e7aSLisandro Dalcin $ FILE_MODE_APPEND - open existing file for binary output 1280cc843e7aSLisandro Dalcin 1281cc843e7aSLisandro Dalcin Level: advanced 1282cc843e7aSLisandro Dalcin 1283db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()` 1284cc843e7aSLisandro Dalcin 1285cc843e7aSLisandro Dalcin @*/ 1286cc843e7aSLisandro Dalcin PetscErrorCode PetscViewerFileGetMode(PetscViewer viewer,PetscFileMode *mode) 1287cc843e7aSLisandro Dalcin { 1288cc843e7aSLisandro Dalcin PetscFunctionBegin; 1289cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1290cc843e7aSLisandro Dalcin PetscValidPointer(mode,2); 1291cac4c232SBarry Smith PetscUseMethod(viewer,"PetscViewerFileGetMode_C",(PetscViewer,PetscFileMode*),(viewer,mode)); 1292cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 1293cc843e7aSLisandro Dalcin } 1294cc843e7aSLisandro Dalcin 1295cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileGetMode_Binary(PetscViewer viewer,PetscFileMode *mode) 1296cc843e7aSLisandro Dalcin { 1297cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 1298cc843e7aSLisandro Dalcin 1299cc843e7aSLisandro Dalcin PetscFunctionBegin; 1300cc843e7aSLisandro Dalcin *mode = vbinary->filemode; 1301cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 1302cc843e7aSLisandro Dalcin } 1303cc843e7aSLisandro Dalcin 1304cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileSetName_Binary(PetscViewer viewer,const char name[]) 1305cc843e7aSLisandro Dalcin { 1306cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 1307cc843e7aSLisandro Dalcin 1308cc843e7aSLisandro Dalcin PetscFunctionBegin; 1309cc843e7aSLisandro Dalcin if (viewer->setupcalled && vbinary->filename) { 1310cc843e7aSLisandro Dalcin /* gzip can be run after the file with the previous filename has been closed */ 13119566063dSJacob Faibussowitsch PetscCall(PetscFree(vbinary->ogzfilename)); 13129566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(vbinary->filename,&vbinary->ogzfilename)); 1313cc843e7aSLisandro Dalcin } 13149566063dSJacob Faibussowitsch PetscCall(PetscFree(vbinary->filename)); 13159566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name,&vbinary->filename)); 1316cc843e7aSLisandro Dalcin viewer->setupcalled = PETSC_FALSE; 1317cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 1318cc843e7aSLisandro Dalcin } 1319cc843e7aSLisandro Dalcin 132081f0254dSBarry Smith static PetscErrorCode PetscViewerFileGetName_Binary(PetscViewer viewer,const char **name) 13215c6c1daeSBarry Smith { 13225c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 13235c6c1daeSBarry Smith 13245c6c1daeSBarry Smith PetscFunctionBegin; 13255c6c1daeSBarry Smith *name = vbinary->filename; 13265c6c1daeSBarry Smith PetscFunctionReturn(0); 13275c6c1daeSBarry Smith } 13285c6c1daeSBarry Smith 13295c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO) 1330e0385b85SDave May static PetscErrorCode PetscViewerFileSetUp_BinaryMPIIO(PetscViewer viewer) 13315c6c1daeSBarry Smith { 13325c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 1333e8a65b7dSLisandro Dalcin int amode; 13345c6c1daeSBarry Smith 13355c6c1daeSBarry Smith PetscFunctionBegin; 1336a297a907SKarl Rupp vbinary->storecompressed = PETSC_FALSE; 13375c6c1daeSBarry Smith 1338cc843e7aSLisandro Dalcin vbinary->moff = 0; 1339cc843e7aSLisandro Dalcin switch (vbinary->filemode) { 1340e8a65b7dSLisandro Dalcin case FILE_MODE_READ: amode = MPI_MODE_RDONLY; break; 1341e8a65b7dSLisandro Dalcin case FILE_MODE_WRITE: amode = MPI_MODE_WRONLY | MPI_MODE_CREATE; break; 134222a8f86cSLisandro Dalcin case FILE_MODE_APPEND: amode = MPI_MODE_WRONLY | MPI_MODE_CREATE | MPI_MODE_APPEND; break; 13437e4fd573SVaclav Hapla case FILE_MODE_UNDEFINED: SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerSetUp()"); 134498921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Unsupported file mode %s",PetscFileModes[vbinary->filemode]); 13455c6c1daeSBarry Smith } 13469566063dSJacob Faibussowitsch PetscCallMPI(MPI_File_open(PetscObjectComm((PetscObject)viewer),vbinary->filename,amode,MPI_INFO_NULL,&vbinary->mfdes)); 134722a8f86cSLisandro Dalcin /* 134822a8f86cSLisandro Dalcin The MPI standard does not have MPI_MODE_TRUNCATE. We emulate this behavior by setting the file size to zero. 134922a8f86cSLisandro Dalcin */ 13509566063dSJacob Faibussowitsch if (vbinary->filemode == FILE_MODE_WRITE) PetscCallMPI(MPI_File_set_size(vbinary->mfdes,0)); 135122a8f86cSLisandro Dalcin /* 135222a8f86cSLisandro Dalcin Initially, all processes view the file as a linear byte stream. Therefore, for files opened with MPI_MODE_APPEND, 135322a8f86cSLisandro Dalcin MPI_File_get_position[_shared](fh, &offset) returns the absolute byte position at the end of file. 135422a8f86cSLisandro Dalcin Otherwise, we would need to call MPI_File_get_byte_offset(fh, offset, &byte_offset) to convert 135522a8f86cSLisandro Dalcin the offset in etype units to an absolute byte position. 135622a8f86cSLisandro Dalcin */ 13579566063dSJacob Faibussowitsch if (vbinary->filemode == FILE_MODE_APPEND) PetscCallMPI(MPI_File_get_position(vbinary->mfdes,&vbinary->moff)); 1358cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 1359cc843e7aSLisandro Dalcin } 1360cc843e7aSLisandro Dalcin #endif 13615c6c1daeSBarry Smith 1362cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileSetUp_BinarySTDIO(PetscViewer viewer) 1363cc843e7aSLisandro Dalcin { 1364cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 1365cc843e7aSLisandro Dalcin const char *fname; 1366cc843e7aSLisandro Dalcin char bname[PETSC_MAX_PATH_LEN],*gz; 1367cc843e7aSLisandro Dalcin PetscBool found; 1368cc843e7aSLisandro Dalcin PetscMPIInt rank; 13695c6c1daeSBarry Smith 1370cc843e7aSLisandro Dalcin PetscFunctionBegin; 13719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank)); 13725c6c1daeSBarry Smith 1373cc843e7aSLisandro Dalcin /* if file name ends in .gz strip that off and note user wants file compressed */ 1374cc843e7aSLisandro Dalcin vbinary->storecompressed = PETSC_FALSE; 1375cc843e7aSLisandro Dalcin if (vbinary->filemode == FILE_MODE_WRITE) { 13769566063dSJacob Faibussowitsch PetscCall(PetscStrstr(vbinary->filename,".gz",&gz)); 1377cc843e7aSLisandro Dalcin if (gz && gz[3] == 0) {*gz = 0; vbinary->storecompressed = PETSC_TRUE;} 1378cc843e7aSLisandro Dalcin } 1379cc843e7aSLisandro Dalcin #if !defined(PETSC_HAVE_POPEN) 138028b400f6SJacob Faibussowitsch PetscCheck(!vbinary->storecompressed,PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP_SYS,"Cannot run gzip on this machine"); 1381cc843e7aSLisandro Dalcin #endif 1382cc843e7aSLisandro Dalcin 1383cc843e7aSLisandro Dalcin fname = vbinary->filename; 1384cc843e7aSLisandro Dalcin if (vbinary->filemode == FILE_MODE_READ) { /* possibly get the file from remote site or compressed file */ 13859566063dSJacob Faibussowitsch PetscCall(PetscFileRetrieve(PetscObjectComm((PetscObject)viewer),fname,bname,PETSC_MAX_PATH_LEN,&found)); 138628b400f6SJacob Faibussowitsch PetscCheck(found,PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_OPEN,"Cannot locate file: %s",fname); 1387cc843e7aSLisandro Dalcin fname = bname; 13885c6c1daeSBarry Smith } 13895c6c1daeSBarry Smith 1390cc843e7aSLisandro Dalcin vbinary->fdes = -1; 1391dd400576SPatrick Sanan if (rank == 0) { /* only first processor opens file*/ 1392cc843e7aSLisandro Dalcin PetscFileMode mode = vbinary->filemode; 1393cc843e7aSLisandro Dalcin if (mode == FILE_MODE_APPEND) { 1394cc843e7aSLisandro Dalcin /* check if asked to append to a non-existing file */ 13959566063dSJacob Faibussowitsch PetscCall(PetscTestFile(fname,'\0',&found)); 1396cc843e7aSLisandro Dalcin if (!found) mode = FILE_MODE_WRITE; 1397cc843e7aSLisandro Dalcin } 13989566063dSJacob Faibussowitsch PetscCall(PetscBinaryOpen(fname,mode,&vbinary->fdes)); 1399cc843e7aSLisandro Dalcin } 1400cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 1401cc843e7aSLisandro Dalcin } 1402cc843e7aSLisandro Dalcin 1403cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileSetUp_BinaryInfo(PetscViewer viewer) 1404cc843e7aSLisandro Dalcin { 1405cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 1406cc843e7aSLisandro Dalcin PetscMPIInt rank; 1407cc843e7aSLisandro Dalcin PetscBool found; 1408cc843e7aSLisandro Dalcin 1409cc843e7aSLisandro Dalcin PetscFunctionBegin; 1410cc843e7aSLisandro Dalcin vbinary->fdes_info = NULL; 14119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank)); 1412dd400576SPatrick Sanan if (!vbinary->skipinfo && (vbinary->filemode == FILE_MODE_READ || rank == 0)) { 1413cc843e7aSLisandro Dalcin char infoname[PETSC_MAX_PATH_LEN],iname[PETSC_MAX_PATH_LEN],*gz; 1414cc843e7aSLisandro Dalcin 14159566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(infoname,vbinary->filename,sizeof(infoname))); 1416cc843e7aSLisandro Dalcin /* remove .gz if it ends file name */ 14179566063dSJacob Faibussowitsch PetscCall(PetscStrstr(infoname,".gz",&gz)); 1418cc843e7aSLisandro Dalcin if (gz && gz[3] == 0) *gz = 0; 1419cc843e7aSLisandro Dalcin 14209566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(infoname,".info",sizeof(infoname))); 1421cc843e7aSLisandro Dalcin if (vbinary->filemode == FILE_MODE_READ) { 14229566063dSJacob Faibussowitsch PetscCall(PetscFixFilename(infoname,iname)); 14239566063dSJacob Faibussowitsch PetscCall(PetscFileRetrieve(PetscObjectComm((PetscObject)viewer),iname,infoname,PETSC_MAX_PATH_LEN,&found)); 14249566063dSJacob Faibussowitsch if (found) PetscCall(PetscOptionsInsertFile(PetscObjectComm((PetscObject)viewer),((PetscObject)viewer)->options,infoname,PETSC_FALSE)); 1425dd400576SPatrick Sanan } else if (rank == 0) { /* write or append */ 1426cc843e7aSLisandro Dalcin const char *omode = (vbinary->filemode == FILE_MODE_APPEND) ? "a" : "w"; 1427cc843e7aSLisandro Dalcin vbinary->fdes_info = fopen(infoname,omode); 142828b400f6SJacob Faibussowitsch PetscCheck(vbinary->fdes_info,PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open .info file %s for writing",infoname); 14295c6c1daeSBarry Smith } 14305c6c1daeSBarry Smith } 14315c6c1daeSBarry Smith PetscFunctionReturn(0); 14325c6c1daeSBarry Smith } 14335c6c1daeSBarry Smith 1434cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerSetUp_Binary(PetscViewer viewer) 14355c6c1daeSBarry Smith { 14365c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 1437cc843e7aSLisandro Dalcin PetscBool usempiio; 1438cc843e7aSLisandro Dalcin 14395c6c1daeSBarry Smith PetscFunctionBegin; 14409566063dSJacob Faibussowitsch if (!vbinary->setfromoptionscalled) PetscCall(PetscViewerSetFromOptions(viewer)); 144128b400f6SJacob Faibussowitsch PetscCheck(vbinary->filename,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call PetscViewerFileSetName()"); 144208401ef6SPierre Jolivet PetscCheck(vbinary->filemode != (PetscFileMode)-1,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call PetscViewerFileSetMode()"); 14439566063dSJacob Faibussowitsch PetscCall(PetscViewerFileClose_Binary(viewer)); 1444cc843e7aSLisandro Dalcin 14459566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetUseMPIIO(viewer,&usempiio)); 1446cc843e7aSLisandro Dalcin if (usempiio) { 1447cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 14489566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetUp_BinaryMPIIO(viewer)); 1449cc843e7aSLisandro Dalcin #endif 1450cc843e7aSLisandro Dalcin } else { 14519566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetUp_BinarySTDIO(viewer)); 1452cc843e7aSLisandro Dalcin } 14539566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetUp_BinaryInfo(viewer)); 1454cc843e7aSLisandro Dalcin 14559566063dSJacob Faibussowitsch PetscCall(PetscLogObjectState((PetscObject)viewer,"File: %s",vbinary->filename)); 14565c6c1daeSBarry Smith PetscFunctionReturn(0); 14575c6c1daeSBarry Smith } 14585c6c1daeSBarry Smith 145981f0254dSBarry Smith static PetscErrorCode PetscViewerView_Binary(PetscViewer v,PetscViewer viewer) 14602bf49c77SBarry Smith { 1461cb6ad94fSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)v->data; 1462cb6ad94fSLisandro Dalcin const char *fname = vbinary->filename ? vbinary->filename : "not yet set"; 1463cc843e7aSLisandro Dalcin const char *fmode = vbinary->filemode != (PetscFileMode) -1 ? PetscFileModes[vbinary->filemode] : "not yet set"; 1464cb6ad94fSLisandro Dalcin PetscBool usempiio; 14652bf49c77SBarry Smith 14662bf49c77SBarry Smith PetscFunctionBegin; 14679566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetUseMPIIO(v,&usempiio)); 14689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Filename: %s\n",fname)); 14699566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Mode: %s (%s)\n",fmode,usempiio ? "mpiio" : "stdio")); 14702bf49c77SBarry Smith PetscFunctionReturn(0); 14712bf49c77SBarry Smith } 14722bf49c77SBarry Smith 147322a8f86cSLisandro Dalcin static PetscErrorCode PetscViewerSetFromOptions_Binary(PetscOptionItems *PetscOptionsObject,PetscViewer viewer) 147403a1d59fSDave May { 147522a8f86cSLisandro Dalcin PetscViewer_Binary *binary = (PetscViewer_Binary*)viewer->data; 1476d22fd6bcSDave May char defaultname[PETSC_MAX_PATH_LEN]; 147703a1d59fSDave May PetscBool flg; 1478e0385b85SDave May 147903a1d59fSDave May PetscFunctionBegin; 148022a8f86cSLisandro Dalcin if (viewer->setupcalled) PetscFunctionReturn(0); 1481d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"Binary PetscViewer Options"); 14829566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(defaultname,PETSC_MAX_PATH_LEN-1,"binaryoutput")); 14839566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-viewer_binary_filename","Specify filename","PetscViewerFileSetName",defaultname,defaultname,sizeof(defaultname),&flg)); 14849566063dSJacob Faibussowitsch if (flg) PetscCall(PetscViewerFileSetName_Binary(viewer,defaultname)); 14859566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_binary_skip_info","Skip writing/reading .info file","PetscViewerBinarySetSkipInfo",binary->skipinfo,&binary->skipinfo,NULL)); 14869566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_binary_skip_options","Skip parsing Vec/Mat load options","PetscViewerBinarySetSkipOptions",binary->skipoptions,&binary->skipoptions,NULL)); 14879566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_binary_skip_header","Skip writing/reading header information","PetscViewerBinarySetSkipHeader",binary->skipheader,&binary->skipheader,NULL)); 148803a1d59fSDave May #if defined(PETSC_HAVE_MPIIO) 14899566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_binary_mpiio","Use MPI-IO functionality to write/read binary file","PetscViewerBinarySetUseMPIIO",binary->usempiio,&binary->usempiio,NULL)); 14905972f5f3SLisandro Dalcin #else 14919566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_binary_mpiio","Use MPI-IO functionality to write/read binary file (NOT AVAILABLE)","PetscViewerBinarySetUseMPIIO",PETSC_FALSE,NULL,NULL)); 149203a1d59fSDave May #endif 1493d0609cedSBarry Smith PetscOptionsHeadEnd(); 1494bc196f7cSDave May binary->setfromoptionscalled = PETSC_TRUE; 149503a1d59fSDave May PetscFunctionReturn(0); 149603a1d59fSDave May } 149703a1d59fSDave May 14988556b5ebSBarry Smith /*MC 14998556b5ebSBarry Smith PETSCVIEWERBINARY - A viewer that saves to binary files 15008556b5ebSBarry Smith 1501db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `PETSC_VIEWER_STDOUT_(),PETSC_VIEWER_STDOUT_SELF`, `PETSC_VIEWER_STDOUT_WORLD`, `PetscViewerCreate()`, `PetscViewerASCIIOpen()`, 1502db781477SPatrick Sanan `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`, `PETSCVIEWERDRAW`, 1503db781477SPatrick Sanan `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()`, 1504db781477SPatrick Sanan `PetscViewerBinaryGetUseMPIIO()`, `PetscViewerBinarySetUseMPIIO()` 15058556b5ebSBarry Smith 15061b266c99SBarry Smith Level: beginner 15071b266c99SBarry Smith 15088556b5ebSBarry Smith M*/ 15098556b5ebSBarry Smith 15108cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_Binary(PetscViewer v) 15115c6c1daeSBarry Smith { 15125c6c1daeSBarry Smith PetscViewer_Binary *vbinary; 15135c6c1daeSBarry Smith 15145c6c1daeSBarry Smith PetscFunctionBegin; 15159566063dSJacob Faibussowitsch PetscCall(PetscNewLog(v,&vbinary)); 15165c6c1daeSBarry Smith v->data = (void*)vbinary; 1517cc843e7aSLisandro Dalcin 151803a1d59fSDave May v->ops->setfromoptions = PetscViewerSetFromOptions_Binary; 15195c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_Binary; 15202bf49c77SBarry Smith v->ops->view = PetscViewerView_Binary; 1521c98fd787SBarry Smith v->ops->setup = PetscViewerSetUp_Binary; 1522cc843e7aSLisandro Dalcin v->ops->flush = NULL; /* Should we support Flush() ? */ 1523cc843e7aSLisandro Dalcin v->ops->getsubviewer = PetscViewerGetSubViewer_Binary; 1524cc843e7aSLisandro Dalcin v->ops->restoresubviewer = PetscViewerRestoreSubViewer_Binary; 1525cc843e7aSLisandro Dalcin v->ops->read = PetscViewerBinaryRead; 1526cc843e7aSLisandro Dalcin 1527cc843e7aSLisandro Dalcin vbinary->fdes = -1; 1528e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 1529cc843e7aSLisandro Dalcin vbinary->usempiio = PETSC_FALSE; 1530e8a65b7dSLisandro Dalcin vbinary->mfdes = MPI_FILE_NULL; 1531e8a65b7dSLisandro Dalcin vbinary->mfsub = MPI_FILE_NULL; 1532e8a65b7dSLisandro Dalcin #endif 1533cc843e7aSLisandro Dalcin vbinary->filename = NULL; 15347e4fd573SVaclav Hapla vbinary->filemode = FILE_MODE_UNDEFINED; 153502c9f0b5SLisandro Dalcin vbinary->fdes_info = NULL; 15365c6c1daeSBarry Smith vbinary->skipinfo = PETSC_FALSE; 15375c6c1daeSBarry Smith vbinary->skipoptions = PETSC_TRUE; 15385c6c1daeSBarry Smith vbinary->skipheader = PETSC_FALSE; 15395c6c1daeSBarry Smith vbinary->storecompressed = PETSC_FALSE; 1540f90597f1SStefano Zampini vbinary->ogzfilename = NULL; 15415c6c1daeSBarry Smith vbinary->flowcontrol = 256; /* seems a good number for Cray XT-5 */ 15425c6c1daeSBarry Smith 1543cc843e7aSLisandro Dalcin vbinary->setfromoptionscalled = PETSC_FALSE; 1544cc843e7aSLisandro Dalcin 15459566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetFlowControl_C",PetscViewerBinaryGetFlowControl_Binary)); 15469566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetFlowControl_C",PetscViewerBinarySetFlowControl_Binary)); 15479566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetSkipHeader_C",PetscViewerBinaryGetSkipHeader_Binary)); 15489566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetSkipHeader_C",PetscViewerBinarySetSkipHeader_Binary)); 15499566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetSkipOptions_C",PetscViewerBinaryGetSkipOptions_Binary)); 15509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetSkipOptions_C",PetscViewerBinarySetSkipOptions_Binary)); 15519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetSkipInfo_C",PetscViewerBinaryGetSkipInfo_Binary)); 15529566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetSkipInfo_C",PetscViewerBinarySetSkipInfo_Binary)); 15539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetInfoPointer_C",PetscViewerBinaryGetInfoPointer_Binary)); 15549566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_Binary)); 15559566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_Binary)); 15569566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_Binary)); 15579566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_Binary)); 15585c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO) 15599566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetUseMPIIO_C",PetscViewerBinaryGetUseMPIIO_Binary)); 15609566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetUseMPIIO_C",PetscViewerBinarySetUseMPIIO_Binary)); 15615c6c1daeSBarry Smith #endif 15625c6c1daeSBarry Smith PetscFunctionReturn(0); 15635c6c1daeSBarry Smith } 15645c6c1daeSBarry Smith 15655c6c1daeSBarry Smith /* ---------------------------------------------------------------------*/ 15665c6c1daeSBarry Smith /* 15675c6c1daeSBarry Smith The variable Petsc_Viewer_Binary_keyval is used to indicate an MPI attribute that 15685c6c1daeSBarry Smith is attached to a communicator, in this case the attribute is a PetscViewer. 15695c6c1daeSBarry Smith */ 1570d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_Binary_keyval = MPI_KEYVAL_INVALID; 15715c6c1daeSBarry Smith 15725c6c1daeSBarry Smith /*@C 15735c6c1daeSBarry Smith PETSC_VIEWER_BINARY_ - Creates a binary PetscViewer shared by all processors 15745c6c1daeSBarry Smith in a communicator. 15755c6c1daeSBarry Smith 1576d083f849SBarry Smith Collective 15775c6c1daeSBarry Smith 15785c6c1daeSBarry Smith Input Parameter: 15795c6c1daeSBarry Smith . comm - the MPI communicator to share the binary PetscViewer 15805c6c1daeSBarry Smith 15815c6c1daeSBarry Smith Level: intermediate 15825c6c1daeSBarry Smith 15835c6c1daeSBarry Smith Options Database Keys: 158410699b91SBarry Smith + -viewer_binary_filename <name> - filename in which to store the binary data, defaults to binaryoutput 158510699b91SBarry Smith . -viewer_binary_skip_info - true means do not create .info file for this viewer 158610699b91SBarry Smith . -viewer_binary_skip_options - true means do not use the options database for this viewer 158710699b91SBarry Smith . -viewer_binary_skip_header - true means do not store the usual header information in the binary file 158810699b91SBarry Smith - -viewer_binary_mpiio - true means use the file via MPI-IO, maybe faster for large files and many MPI ranks 15895c6c1daeSBarry Smith 15905c6c1daeSBarry Smith Environmental variables: 159110699b91SBarry Smith - PETSC_VIEWER_BINARY_FILENAME - filename in which to store the binary data, defaults to binaryoutput 15925c6c1daeSBarry Smith 15935c6c1daeSBarry Smith Notes: 15945c6c1daeSBarry Smith Unlike almost all other PETSc routines, PETSC_VIEWER_BINARY_ does not return 15955c6c1daeSBarry Smith an error code. The binary PetscViewer is usually used in the form 15965c6c1daeSBarry Smith $ XXXView(XXX object,PETSC_VIEWER_BINARY_(comm)); 15975c6c1daeSBarry Smith 1598db781477SPatrick Sanan .seealso: `PETSC_VIEWER_BINARY_WORLD`, `PETSC_VIEWER_BINARY_SELF`, `PetscViewerBinaryOpen()`, `PetscViewerCreate()`, 1599db781477SPatrick Sanan `PetscViewerDestroy()` 16005c6c1daeSBarry Smith @*/ 16015c6c1daeSBarry Smith PetscViewer PETSC_VIEWER_BINARY_(MPI_Comm comm) 16025c6c1daeSBarry Smith { 16035c6c1daeSBarry Smith PetscErrorCode ierr; 16045c6c1daeSBarry Smith PetscBool flg; 16055c6c1daeSBarry Smith PetscViewer viewer; 16065c6c1daeSBarry Smith char fname[PETSC_MAX_PATH_LEN]; 16075c6c1daeSBarry Smith MPI_Comm ncomm; 16085c6c1daeSBarry Smith 16095c6c1daeSBarry Smith PetscFunctionBegin; 161002c9f0b5SLisandro 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);} 16115c6c1daeSBarry Smith if (Petsc_Viewer_Binary_keyval == MPI_KEYVAL_INVALID) { 161202c9f0b5SLisandro Dalcin ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_Binary_keyval,NULL); 161302c9f0b5SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 16145c6c1daeSBarry Smith } 161547435625SJed Brown ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_Binary_keyval,(void**)&viewer,(int*)&flg); 161602c9f0b5SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 16175c6c1daeSBarry Smith if (!flg) { /* PetscViewer not yet created */ 16185c6c1daeSBarry Smith ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_BINARY_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 16192cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 16205c6c1daeSBarry Smith if (!flg) { 16215c6c1daeSBarry Smith ierr = PetscStrcpy(fname,"binaryoutput"); 16222cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 16235c6c1daeSBarry Smith } 16245c6c1daeSBarry Smith ierr = PetscViewerBinaryOpen(ncomm,fname,FILE_MODE_WRITE,&viewer); 16252cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 16265c6c1daeSBarry Smith ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 16272cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 162847435625SJed Brown ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_Binary_keyval,(void*)viewer); 162902c9f0b5SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 16305c6c1daeSBarry Smith } 16315c6c1daeSBarry Smith ierr = PetscCommDestroy(&ncomm); 16322cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 16335c6c1daeSBarry Smith PetscFunctionReturn(viewer); 16345c6c1daeSBarry Smith } 1635