1af0996ceSBarry Smith #include <petsc/private/viewerimpl.h> /*I "petscviewer.h" I*/ 25c6c1daeSBarry Smith 35c6c1daeSBarry Smith typedef struct { 45c6c1daeSBarry Smith int fdes; /* file descriptor, ignored if using MPI IO */ 55c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO) 6bc196f7cSDave May PetscBool usempiio; 75c6c1daeSBarry Smith MPI_File mfdes; /* ignored unless using MPI IO */ 8e8a65b7dSLisandro Dalcin MPI_File mfsub; /* subviewer support */ 95c6c1daeSBarry Smith MPI_Offset moff; 105c6c1daeSBarry Smith #endif 11cc843e7aSLisandro Dalcin char *filename; /* file name */ 12cc843e7aSLisandro Dalcin PetscFileMode filemode; /* read/write/append mode */ 135c6c1daeSBarry Smith FILE *fdes_info; /* optional file containing info on binary file*/ 145c6c1daeSBarry Smith PetscBool storecompressed; /* gzip the write binary file when closing it*/ 15f90597f1SStefano Zampini char *ogzfilename; /* gzip can be run after the filename has been updated */ 165c6c1daeSBarry Smith PetscBool skipinfo; /* Don't create info file for writing; don't use for reading */ 175c6c1daeSBarry Smith PetscBool skipoptions; /* don't use PETSc options database when loading */ 185c6c1daeSBarry Smith PetscInt flowcontrol; /* allow only <flowcontrol> messages outstanding at a time while doing IO */ 195c6c1daeSBarry Smith PetscBool skipheader; /* don't write header, only raw data */ 20a261c58fSBarry Smith PetscBool matlabheaderwritten; /* if format is PETSC_VIEWER_BINARY_MATLAB has the MATLAB .info header been written yet */ 21c98fd787SBarry Smith PetscBool setfromoptionscalled; 225c6c1daeSBarry Smith } PetscViewer_Binary; 235c6c1daeSBarry Smith 24cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 25cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinarySyncMPIIO(PetscViewer viewer) 26cc843e7aSLisandro Dalcin { 27cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 28cc843e7aSLisandro Dalcin PetscErrorCode ierr; 29cc843e7aSLisandro Dalcin 30cc843e7aSLisandro Dalcin PetscFunctionBegin; 31cc843e7aSLisandro Dalcin if (vbinary->filemode == FILE_MODE_READ) PetscFunctionReturn(0); 32cc843e7aSLisandro Dalcin if (vbinary->mfsub != MPI_FILE_NULL) { 33ffc4695bSBarry Smith ierr = MPI_File_sync(vbinary->mfsub);CHKERRMPI(ierr); 34cc843e7aSLisandro Dalcin } 35cc843e7aSLisandro Dalcin if (vbinary->mfdes != MPI_FILE_NULL) { 36ffc4695bSBarry Smith ierr = MPI_Barrier(PetscObjectComm((PetscObject)viewer));CHKERRMPI(ierr); 37ffc4695bSBarry Smith ierr = MPI_File_sync(vbinary->mfdes);CHKERRMPI(ierr); 38cc843e7aSLisandro Dalcin } 39cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 40cc843e7aSLisandro Dalcin } 41cc843e7aSLisandro Dalcin #endif 42cc843e7aSLisandro Dalcin 4381f0254dSBarry Smith static PetscErrorCode PetscViewerGetSubViewer_Binary(PetscViewer viewer,MPI_Comm comm,PetscViewer *outviewer) 445c6c1daeSBarry Smith { 45e8a65b7dSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 46cc843e7aSLisandro Dalcin PetscMPIInt rank; 47cc843e7aSLisandro Dalcin PetscErrorCode ierr; 485c6c1daeSBarry Smith 495c6c1daeSBarry Smith PetscFunctionBegin; 50c98fd787SBarry Smith ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr); 51e8a65b7dSLisandro Dalcin 52e8a65b7dSLisandro Dalcin /* Return subviewer in process zero */ 53ffc4695bSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank);CHKERRMPI(ierr); 54dd400576SPatrick Sanan if (rank == 0) { 55e5afcf28SBarry Smith PetscMPIInt flg; 56e5afcf28SBarry Smith 57ffc4695bSBarry Smith ierr = MPI_Comm_compare(PETSC_COMM_SELF,comm,&flg);CHKERRMPI(ierr); 58e5afcf28SBarry Smith if (flg != MPI_IDENT && flg != MPI_CONGRUENT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscViewerGetSubViewer() for PETSCVIEWERBINARY requires a singleton MPI_Comm"); 59e5afcf28SBarry Smith ierr = PetscViewerCreate(comm,outviewer);CHKERRQ(ierr); 605c6c1daeSBarry Smith ierr = PetscViewerSetType(*outviewer,PETSCVIEWERBINARY);CHKERRQ(ierr); 61e8a65b7dSLisandro Dalcin ierr = PetscMemcpy((*outviewer)->data,vbinary,sizeof(PetscViewer_Binary));CHKERRQ(ierr); 6212f4c3a9SLisandro Dalcin (*outviewer)->setupcalled = PETSC_TRUE; 6396509d17SLisandro Dalcin } else { 6496509d17SLisandro Dalcin *outviewer = NULL; 6596509d17SLisandro Dalcin } 66e8a65b7dSLisandro Dalcin 67e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 68e8a65b7dSLisandro Dalcin if (vbinary->usempiio && *outviewer) { 69e8a65b7dSLisandro Dalcin PetscViewer_Binary *obinary = (PetscViewer_Binary*)(*outviewer)->data; 70e8a65b7dSLisandro Dalcin /* Parent viewer opens a new MPI file handle on PETSC_COMM_SELF and keeps track of it for future reuse */ 71e8a65b7dSLisandro Dalcin if (vbinary->mfsub == MPI_FILE_NULL) { 72e8a65b7dSLisandro Dalcin int amode; 73cc843e7aSLisandro Dalcin switch (vbinary->filemode) { 74e8a65b7dSLisandro Dalcin case FILE_MODE_READ: amode = MPI_MODE_RDONLY; break; 75e8a65b7dSLisandro Dalcin case FILE_MODE_WRITE: amode = MPI_MODE_WRONLY; break; 7622a8f86cSLisandro Dalcin case FILE_MODE_APPEND: amode = MPI_MODE_WRONLY; break; 77*98921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported file mode %s",PetscFileModes[vbinary->filemode]); 78e8a65b7dSLisandro Dalcin } 79ffc4695bSBarry Smith ierr = MPI_File_open(PETSC_COMM_SELF,vbinary->filename,amode,MPI_INFO_NULL,&vbinary->mfsub);CHKERRMPI(ierr); 80e8a65b7dSLisandro Dalcin } 81e8a65b7dSLisandro Dalcin /* Subviewer gets the MPI file handle on PETSC_COMM_SELF */ 82e8a65b7dSLisandro Dalcin obinary->mfdes = vbinary->mfsub; 83e8a65b7dSLisandro Dalcin obinary->mfsub = MPI_FILE_NULL; 84e8a65b7dSLisandro Dalcin obinary->moff = vbinary->moff; 85e8a65b7dSLisandro Dalcin } 86e8a65b7dSLisandro Dalcin #endif 87cc843e7aSLisandro Dalcin 88cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 89cc843e7aSLisandro Dalcin ierr = PetscViewerBinarySyncMPIIO(viewer);CHKERRQ(ierr); 90cc843e7aSLisandro Dalcin #endif 915c6c1daeSBarry Smith PetscFunctionReturn(0); 925c6c1daeSBarry Smith } 935c6c1daeSBarry Smith 9481f0254dSBarry Smith static PetscErrorCode PetscViewerRestoreSubViewer_Binary(PetscViewer viewer,MPI_Comm comm,PetscViewer *outviewer) 955c6c1daeSBarry Smith { 96e8a65b7dSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 97cc843e7aSLisandro Dalcin PetscMPIInt rank; 98e8a65b7dSLisandro Dalcin PetscErrorCode ierr; 99e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 100e8a65b7dSLisandro Dalcin MPI_Offset moff = 0; 101e8a65b7dSLisandro Dalcin #endif 1025c6c1daeSBarry Smith 1035c6c1daeSBarry Smith PetscFunctionBegin; 104ffc4695bSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank);CHKERRMPI(ierr); 105e8a65b7dSLisandro Dalcin if (rank && *outviewer) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Subviewer not obtained from viewer"); 106e8a65b7dSLisandro Dalcin 107e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 108e8a65b7dSLisandro Dalcin if (vbinary->usempiio && *outviewer) { 109e8a65b7dSLisandro Dalcin PetscViewer_Binary *obinary = (PetscViewer_Binary*)(*outviewer)->data; 110e8a65b7dSLisandro Dalcin if (obinary->mfdes != vbinary->mfsub) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Subviewer not obtained from viewer"); 111ffc4695bSBarry Smith if (obinary->mfsub != MPI_FILE_NULL) {ierr = MPI_File_close(&obinary->mfsub);CHKERRMPI(ierr);} 112e8a65b7dSLisandro Dalcin moff = obinary->moff; 113e8a65b7dSLisandro Dalcin } 114e8a65b7dSLisandro Dalcin #endif 115e8a65b7dSLisandro Dalcin 116e8a65b7dSLisandro Dalcin if (*outviewer) { 117e8a65b7dSLisandro Dalcin PetscViewer_Binary *obinary = (PetscViewer_Binary*)(*outviewer)->data; 118e8a65b7dSLisandro Dalcin if (obinary->fdes != vbinary->fdes) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Subviewer not obtained from viewer"); 1195c6c1daeSBarry Smith ierr = PetscFree((*outviewer)->data);CHKERRQ(ierr); 1205c6c1daeSBarry Smith ierr = PetscHeaderDestroy(outviewer);CHKERRQ(ierr); 1215c6c1daeSBarry Smith } 122e8a65b7dSLisandro Dalcin 123e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 124e8a65b7dSLisandro Dalcin if (vbinary->usempiio) { 125e8a65b7dSLisandro Dalcin PetscInt64 ioff = (PetscInt64)moff; /* We could use MPI_OFFSET datatype (requires MPI 2.2) */ 126ffc4695bSBarry Smith ierr = MPI_Bcast(&ioff,1,MPIU_INT64,0,PetscObjectComm((PetscObject)viewer));CHKERRMPI(ierr); 127e8a65b7dSLisandro Dalcin vbinary->moff = (MPI_Offset)ioff; 128e8a65b7dSLisandro Dalcin } 129e8a65b7dSLisandro Dalcin #endif 130cc843e7aSLisandro Dalcin 131cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 132cc843e7aSLisandro Dalcin ierr = PetscViewerBinarySyncMPIIO(viewer);CHKERRQ(ierr); 133cc843e7aSLisandro Dalcin #endif 1345c6c1daeSBarry Smith PetscFunctionReturn(0); 1355c6c1daeSBarry Smith } 1365c6c1daeSBarry Smith 1375c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO) 1385c6c1daeSBarry Smith /*@C 13922a8f86cSLisandro Dalcin PetscViewerBinaryGetMPIIOOffset - Gets the current global offset that should be passed to MPI_File_set_view() or MPI_File_{write|read}_at[_all]() 1405c6c1daeSBarry Smith 1415c6c1daeSBarry Smith Not Collective 1425c6c1daeSBarry Smith 1435c6c1daeSBarry Smith Input Parameter: 1445c6c1daeSBarry Smith . viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 1455c6c1daeSBarry Smith 1465c6c1daeSBarry Smith Output Parameter: 14722a8f86cSLisandro Dalcin . off - the current global offset 1485c6c1daeSBarry Smith 1495c6c1daeSBarry Smith Level: advanced 1505c6c1daeSBarry Smith 1515c6c1daeSBarry Smith Fortran Note: 1525c6c1daeSBarry Smith This routine is not supported in Fortran. 1535c6c1daeSBarry Smith 1545c6c1daeSBarry Smith Use PetscViewerBinaryAddMPIIOOffset() to increase this value after you have written a view. 1555c6c1daeSBarry Smith 15667918a83SBarry Smith .seealso: PetscViewerBinaryOpen(), PetscViewerBinaryGetInfoPointer(), PetscViewerBinaryGetUseMPIIO(), PetscViewerBinarySetUseMPIIO(), PetscViewerBinaryAddMPIIOOffset() 1575c6c1daeSBarry Smith @*/ 1585c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryGetMPIIOOffset(PetscViewer viewer,MPI_Offset *off) 1595c6c1daeSBarry Smith { 16022a8f86cSLisandro Dalcin PetscViewer_Binary *vbinary; 1615c6c1daeSBarry Smith 1625c6c1daeSBarry Smith PetscFunctionBegin; 16322a8f86cSLisandro Dalcin PetscValidHeaderSpecificType(viewer,PETSC_VIEWER_CLASSID,1,PETSCVIEWERBINARY); 16422a8f86cSLisandro Dalcin PetscValidPointer(off,2); 16522a8f86cSLisandro Dalcin vbinary = (PetscViewer_Binary*)viewer->data; 1665c6c1daeSBarry Smith *off = vbinary->moff; 1675c6c1daeSBarry Smith PetscFunctionReturn(0); 1685c6c1daeSBarry Smith } 1695c6c1daeSBarry Smith 1705c6c1daeSBarry Smith /*@C 17122a8f86cSLisandro Dalcin PetscViewerBinaryAddMPIIOOffset - Adds to the current global offset 1725c6c1daeSBarry Smith 17322a8f86cSLisandro Dalcin Logically Collective 1745c6c1daeSBarry Smith 1755c6c1daeSBarry Smith Input Parameters: 1765c6c1daeSBarry Smith + viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 17722a8f86cSLisandro Dalcin - off - the addition to the global offset 1785c6c1daeSBarry Smith 1795c6c1daeSBarry Smith Level: advanced 1805c6c1daeSBarry Smith 1815c6c1daeSBarry Smith Fortran Note: 1825c6c1daeSBarry Smith This routine is not supported in Fortran. 1835c6c1daeSBarry Smith 18422a8f86cSLisandro Dalcin Use PetscViewerBinaryGetMPIIOOffset() to get the value that you should pass to MPI_File_set_view() or MPI_File_{write|read}_at[_all]() 1855c6c1daeSBarry Smith 18667918a83SBarry Smith .seealso: PetscViewerBinaryOpen(),PetscViewerBinaryGetInfoPointer(), PetscViewerBinaryGetUseMPIIO(), PetscViewerBinarySetUseMPIIO(), PetscViewerBinaryGetMPIIOOffset() 1875c6c1daeSBarry Smith @*/ 1885c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryAddMPIIOOffset(PetscViewer viewer,MPI_Offset off) 1895c6c1daeSBarry Smith { 19022a8f86cSLisandro Dalcin PetscViewer_Binary *vbinary; 1915c6c1daeSBarry Smith 1925c6c1daeSBarry Smith PetscFunctionBegin; 19322a8f86cSLisandro Dalcin PetscValidHeaderSpecificType(viewer,PETSC_VIEWER_CLASSID,1,PETSCVIEWERBINARY); 19422a8f86cSLisandro Dalcin PetscValidLogicalCollectiveInt(viewer,(PetscInt)off,2); 19522a8f86cSLisandro Dalcin vbinary = (PetscViewer_Binary*)viewer->data; 1965c6c1daeSBarry Smith vbinary->moff += off; 1975c6c1daeSBarry Smith PetscFunctionReturn(0); 1985c6c1daeSBarry Smith } 1995c6c1daeSBarry Smith 2005c6c1daeSBarry Smith /*@C 2015c6c1daeSBarry Smith PetscViewerBinaryGetMPIIODescriptor - Extracts the MPI IO file descriptor from a PetscViewer. 2025c6c1daeSBarry Smith 2035c6c1daeSBarry Smith Not Collective 2045c6c1daeSBarry Smith 2055c6c1daeSBarry Smith Input Parameter: 2065c6c1daeSBarry Smith . viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 2075c6c1daeSBarry Smith 2085c6c1daeSBarry Smith Output Parameter: 2095c6c1daeSBarry Smith . fdes - file descriptor 2105c6c1daeSBarry Smith 2115c6c1daeSBarry Smith Level: advanced 2125c6c1daeSBarry Smith 2135c6c1daeSBarry Smith Fortran Note: 2145c6c1daeSBarry Smith This routine is not supported in Fortran. 2155c6c1daeSBarry Smith 21667918a83SBarry Smith .seealso: PetscViewerBinaryOpen(),PetscViewerBinaryGetInfoPointer(), PetscViewerBinaryGetUseMPIIO(), PetscViewerBinarySetUseMPIIO(), PetscViewerBinaryGetMPIIOOffset() 2175c6c1daeSBarry Smith @*/ 2185c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryGetMPIIODescriptor(PetscViewer viewer,MPI_File *fdes) 2195c6c1daeSBarry Smith { 22003a1d59fSDave May PetscErrorCode ierr; 22122a8f86cSLisandro Dalcin PetscViewer_Binary *vbinary; 2225c6c1daeSBarry Smith 2235c6c1daeSBarry Smith PetscFunctionBegin; 22422a8f86cSLisandro Dalcin PetscValidHeaderSpecificType(viewer,PETSC_VIEWER_CLASSID,1,PETSCVIEWERBINARY); 22522a8f86cSLisandro Dalcin PetscValidPointer(fdes,2); 226c98fd787SBarry Smith ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr); 22722a8f86cSLisandro Dalcin vbinary = (PetscViewer_Binary*)viewer->data; 2285c6c1daeSBarry Smith *fdes = vbinary->mfdes; 2295c6c1daeSBarry Smith PetscFunctionReturn(0); 2305c6c1daeSBarry Smith } 231cc843e7aSLisandro Dalcin #endif 2325c6c1daeSBarry Smith 233cc843e7aSLisandro Dalcin /*@ 234cc843e7aSLisandro Dalcin PetscViewerBinarySetUseMPIIO - Sets a binary viewer to use MPI-IO for reading/writing. Must be called 235cc843e7aSLisandro Dalcin before PetscViewerFileSetName() 236cc843e7aSLisandro Dalcin 237cc843e7aSLisandro Dalcin Logically Collective on PetscViewer 238cc843e7aSLisandro Dalcin 239cc843e7aSLisandro Dalcin Input Parameters: 240cc843e7aSLisandro Dalcin + viewer - the PetscViewer; must be a binary 241cc843e7aSLisandro Dalcin - use - PETSC_TRUE means MPI-IO will be used 242cc843e7aSLisandro Dalcin 243cc843e7aSLisandro Dalcin Options Database: 244cc843e7aSLisandro Dalcin -viewer_binary_mpiio : Flag for using MPI-IO 245cc843e7aSLisandro Dalcin 246cc843e7aSLisandro Dalcin Level: advanced 247cc843e7aSLisandro Dalcin 248cc843e7aSLisandro Dalcin .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 249cc843e7aSLisandro Dalcin PetscViewerBinaryGetUseMPIIO() 250cc843e7aSLisandro Dalcin 251cc843e7aSLisandro Dalcin @*/ 252cc843e7aSLisandro Dalcin PetscErrorCode PetscViewerBinarySetUseMPIIO(PetscViewer viewer,PetscBool use) 253a8aae444SDave May { 254cc843e7aSLisandro Dalcin PetscErrorCode ierr; 255a8aae444SDave May 256a8aae444SDave May PetscFunctionBegin; 257cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 258cc843e7aSLisandro Dalcin PetscValidLogicalCollectiveBool(viewer,use,2); 259cc843e7aSLisandro Dalcin ierr = PetscTryMethod(viewer,"PetscViewerBinarySetUseMPIIO_C",(PetscViewer,PetscBool),(viewer,use));CHKERRQ(ierr); 260cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 261cc843e7aSLisandro Dalcin } 262cc843e7aSLisandro Dalcin 263cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 264cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinarySetUseMPIIO_Binary(PetscViewer viewer,PetscBool use) 265cc843e7aSLisandro Dalcin { 266cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 267cc843e7aSLisandro Dalcin PetscFunctionBegin; 268*98921bdaSJacob Faibussowitsch if (viewer->setupcalled && vbinary->usempiio != use) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER,"Cannot change MPIIO to %s after setup",PetscBools[use]); 269cc843e7aSLisandro Dalcin vbinary->usempiio = use; 270a8aae444SDave May PetscFunctionReturn(0); 271a8aae444SDave May } 272a8aae444SDave May #endif 273a8aae444SDave May 274cc843e7aSLisandro Dalcin /*@ 275bc196f7cSDave May PetscViewerBinaryGetUseMPIIO - Returns PETSC_TRUE if the binary viewer uses MPI-IO. 2765c6c1daeSBarry Smith 2775c6c1daeSBarry Smith Not Collective 2785c6c1daeSBarry Smith 2795c6c1daeSBarry Smith Input Parameter: 2805c6c1daeSBarry Smith . viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 2815c6c1daeSBarry Smith 2825c6c1daeSBarry Smith Output Parameter: 283d8d19677SJose E. Roman . use - PETSC_TRUE if MPI-IO is being used 2845c6c1daeSBarry Smith 2855c6c1daeSBarry Smith Options Database: 2865c6c1daeSBarry Smith -viewer_binary_mpiio : Flag for using MPI-IO 2875c6c1daeSBarry Smith 2885c6c1daeSBarry Smith Level: advanced 2895c6c1daeSBarry Smith 290bc196f7cSDave May Note: 291bc196f7cSDave May If MPI-IO is not available, this function will always return PETSC_FALSE 292bc196f7cSDave May 2935c6c1daeSBarry Smith Fortran Note: 2945c6c1daeSBarry Smith This routine is not supported in Fortran. 2955c6c1daeSBarry Smith 29667918a83SBarry Smith .seealso: PetscViewerBinaryOpen(), PetscViewerBinaryGetInfoPointer(), PetscViewerBinarySetUseMPIIO(), PetscViewerBinaryGetMPIIOOffset() 2975c6c1daeSBarry Smith @*/ 298cc843e7aSLisandro Dalcin PetscErrorCode PetscViewerBinaryGetUseMPIIO(PetscViewer viewer,PetscBool *use) 2995c6c1daeSBarry Smith { 300a8aae444SDave May PetscErrorCode ierr; 3015c6c1daeSBarry Smith 3025c6c1daeSBarry Smith PetscFunctionBegin; 303cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 304cc843e7aSLisandro Dalcin PetscValidBoolPointer(use,2); 305cc843e7aSLisandro Dalcin *use = PETSC_FALSE; 306cc843e7aSLisandro Dalcin ierr = PetscTryMethod(viewer,"PetscViewerBinaryGetUseMPIIO_C",(PetscViewer,PetscBool*),(viewer,use));CHKERRQ(ierr); 3075c6c1daeSBarry Smith PetscFunctionReturn(0); 3085c6c1daeSBarry Smith } 3095c6c1daeSBarry Smith 310cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 311cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinaryGetUseMPIIO_Binary(PetscViewer viewer,PetscBool *use) 3125c6c1daeSBarry Smith { 3135c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 3145c6c1daeSBarry Smith 3155c6c1daeSBarry Smith PetscFunctionBegin; 316cc843e7aSLisandro Dalcin *use = vbinary->usempiio; 317cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 318cc843e7aSLisandro Dalcin } 319cc843e7aSLisandro Dalcin #endif 320cc843e7aSLisandro Dalcin 321cc843e7aSLisandro Dalcin /*@ 322cc843e7aSLisandro Dalcin PetscViewerBinarySetFlowControl - Sets how many messages are allowed to outstanding at the same time during parallel IO reads/writes 323cc843e7aSLisandro Dalcin 324cc843e7aSLisandro Dalcin Not Collective 325cc843e7aSLisandro Dalcin 326d8d19677SJose E. Roman Input Parameters: 327cc843e7aSLisandro Dalcin + viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 328cc843e7aSLisandro Dalcin - fc - the number of messages, defaults to 256 if this function was not called 329cc843e7aSLisandro Dalcin 330cc843e7aSLisandro Dalcin Level: advanced 331cc843e7aSLisandro Dalcin 332cc843e7aSLisandro Dalcin .seealso: PetscViewerBinaryOpen(),PetscViewerBinaryGetInfoPointer(), PetscViewerBinaryGetFlowControl() 333cc843e7aSLisandro Dalcin 334cc843e7aSLisandro Dalcin @*/ 335cc843e7aSLisandro Dalcin PetscErrorCode PetscViewerBinarySetFlowControl(PetscViewer viewer,PetscInt fc) 336cc843e7aSLisandro Dalcin { 337cc843e7aSLisandro Dalcin PetscErrorCode ierr; 338cc843e7aSLisandro Dalcin 339cc843e7aSLisandro Dalcin PetscFunctionBegin; 340cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 341cc843e7aSLisandro Dalcin PetscValidLogicalCollectiveInt(viewer,fc,2); 342cc843e7aSLisandro Dalcin ierr = PetscTryMethod(viewer,"PetscViewerBinarySetFlowControl_C",(PetscViewer,PetscInt),(viewer,fc));CHKERRQ(ierr); 3435c6c1daeSBarry Smith PetscFunctionReturn(0); 3445c6c1daeSBarry Smith } 3455c6c1daeSBarry Smith 346cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinarySetFlowControl_Binary(PetscViewer viewer,PetscInt fc) 347cc843e7aSLisandro Dalcin { 348cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 349cc843e7aSLisandro Dalcin 350cc843e7aSLisandro Dalcin PetscFunctionBegin; 351*98921bdaSJacob Faibussowitsch if (fc <= 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_OUTOFRANGE,"Flow control count must be greater than 1, %" PetscInt_FMT " was set",fc); 352cc843e7aSLisandro Dalcin vbinary->flowcontrol = fc; 353cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 354cc843e7aSLisandro Dalcin } 355cc843e7aSLisandro Dalcin 356cc843e7aSLisandro Dalcin /*@ 3575c6c1daeSBarry Smith PetscViewerBinaryGetFlowControl - Returns how many messages are allowed to outstanding at the same time during parallel IO reads/writes 3585c6c1daeSBarry Smith 3595c6c1daeSBarry Smith Not Collective 3605c6c1daeSBarry Smith 3615c6c1daeSBarry Smith Input Parameter: 3625c6c1daeSBarry Smith . viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 3635c6c1daeSBarry Smith 3645c6c1daeSBarry Smith Output Parameter: 3655c6c1daeSBarry Smith . fc - the number of messages 3665c6c1daeSBarry Smith 3675c6c1daeSBarry Smith Level: advanced 3685c6c1daeSBarry Smith 3695c6c1daeSBarry Smith .seealso: PetscViewerBinaryOpen(),PetscViewerBinaryGetInfoPointer(), PetscViewerBinarySetFlowControl() 3705c6c1daeSBarry Smith 3715c6c1daeSBarry Smith @*/ 3725c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryGetFlowControl(PetscViewer viewer,PetscInt *fc) 3735c6c1daeSBarry Smith { 3745c6c1daeSBarry Smith PetscErrorCode ierr; 3755c6c1daeSBarry Smith 3765c6c1daeSBarry Smith PetscFunctionBegin; 377cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 378cc843e7aSLisandro Dalcin PetscValidIntPointer(fc,2); 379163d334eSBarry Smith ierr = PetscUseMethod(viewer,"PetscViewerBinaryGetFlowControl_C",(PetscViewer,PetscInt*),(viewer,fc));CHKERRQ(ierr); 3805c6c1daeSBarry Smith PetscFunctionReturn(0); 3815c6c1daeSBarry Smith } 3825c6c1daeSBarry Smith 383cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinaryGetFlowControl_Binary(PetscViewer viewer,PetscInt *fc) 3845c6c1daeSBarry Smith { 3855c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 3865c6c1daeSBarry Smith 3875c6c1daeSBarry Smith PetscFunctionBegin; 388cc843e7aSLisandro Dalcin *fc = vbinary->flowcontrol; 3895c6c1daeSBarry Smith PetscFunctionReturn(0); 3905c6c1daeSBarry Smith } 3915c6c1daeSBarry Smith 3925c6c1daeSBarry Smith /*@C 3935c6c1daeSBarry Smith PetscViewerBinaryGetDescriptor - Extracts the file descriptor from a PetscViewer. 3945c6c1daeSBarry Smith 3955872f025SBarry Smith Collective On PetscViewer 3965c6c1daeSBarry Smith 3975c6c1daeSBarry Smith Input Parameter: 3985c6c1daeSBarry Smith . viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 3995c6c1daeSBarry Smith 4005c6c1daeSBarry Smith Output Parameter: 4015c6c1daeSBarry Smith . fdes - file descriptor 4025c6c1daeSBarry Smith 4035c6c1daeSBarry Smith Level: advanced 4045c6c1daeSBarry Smith 4055c6c1daeSBarry Smith Notes: 4065c6c1daeSBarry Smith For writable binary PetscViewers, the descriptor will only be valid for the 4075c6c1daeSBarry Smith first processor in the communicator that shares the PetscViewer. For readable 4085c6c1daeSBarry Smith files it will only be valid on nodes that have the file. If node 0 does not 4095c6c1daeSBarry Smith have the file it generates an error even if another node does have the file. 4105c6c1daeSBarry Smith 4115c6c1daeSBarry Smith Fortran Note: 4125c6c1daeSBarry Smith This routine is not supported in Fortran. 4135c6c1daeSBarry Smith 41495452b02SPatrick Sanan Developer Notes: 41595452b02SPatrick Sanan This must be called on all processes because Dave May changed 4165872f025SBarry Smith the source code that this may be trigger a PetscViewerSetUp() call if it was not previously triggered. 4175872f025SBarry Smith 4185c6c1daeSBarry Smith .seealso: PetscViewerBinaryOpen(),PetscViewerBinaryGetInfoPointer() 4195c6c1daeSBarry Smith @*/ 4205c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryGetDescriptor(PetscViewer viewer,int *fdes) 4215c6c1daeSBarry Smith { 42203a1d59fSDave May PetscErrorCode ierr; 42322a8f86cSLisandro Dalcin PetscViewer_Binary *vbinary; 4245c6c1daeSBarry Smith 4255c6c1daeSBarry Smith PetscFunctionBegin; 42622a8f86cSLisandro Dalcin PetscValidHeaderSpecificType(viewer,PETSC_VIEWER_CLASSID,1,PETSCVIEWERBINARY); 42722a8f86cSLisandro Dalcin PetscValidPointer(fdes,2); 428c98fd787SBarry Smith ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr); 42922a8f86cSLisandro Dalcin vbinary = (PetscViewer_Binary*)viewer->data; 4305c6c1daeSBarry Smith *fdes = vbinary->fdes; 4315c6c1daeSBarry Smith PetscFunctionReturn(0); 4325c6c1daeSBarry Smith } 4335c6c1daeSBarry Smith 4345c6c1daeSBarry Smith /*@ 4355c6c1daeSBarry Smith PetscViewerBinarySkipInfo - Binary file will not have .info file created with it 4365c6c1daeSBarry Smith 4375c6c1daeSBarry Smith Not Collective 4385c6c1daeSBarry Smith 439fd292e60Sprj- Input Parameter: 4405c6c1daeSBarry Smith . viewer - PetscViewer context, obtained from PetscViewerCreate() 4415c6c1daeSBarry Smith 4425c6c1daeSBarry Smith Options Database Key: 44310699b91SBarry Smith . -viewer_binary_skip_info - true indicates do not generate .info file 4445c6c1daeSBarry Smith 4455c6c1daeSBarry Smith Level: advanced 4465c6c1daeSBarry Smith 44795452b02SPatrick Sanan Notes: 44895452b02SPatrick Sanan This must be called after PetscViewerSetType(). If you use PetscViewerBinaryOpen() then 4495c6c1daeSBarry Smith you can only skip the info file with the -viewer_binary_skip_info flag. To use the function you must open the 450a2d7db39SDave May viewer with PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinarySkipInfo(). 4515c6c1daeSBarry Smith 4525c6c1daeSBarry Smith The .info contains meta information about the data in the binary file, for example the block size if it was 4535c6c1daeSBarry Smith set for a vector or matrix. 4545c6c1daeSBarry Smith 4555c6c1daeSBarry Smith .seealso: PetscViewerBinaryOpen(), PetscViewerBinaryGetDescriptor(), PetscViewerBinarySetSkipOptions(), 456807ea322SDave May PetscViewerBinaryGetSkipOptions(), PetscViewerBinaryGetSkipInfo() 4575c6c1daeSBarry Smith @*/ 4585c6c1daeSBarry Smith PetscErrorCode PetscViewerBinarySkipInfo(PetscViewer viewer) 4595c6c1daeSBarry Smith { 460cc843e7aSLisandro Dalcin PetscErrorCode ierr; 4615c6c1daeSBarry Smith 4625c6c1daeSBarry Smith PetscFunctionBegin; 463cc843e7aSLisandro Dalcin ierr = PetscViewerBinarySetSkipInfo(viewer,PETSC_TRUE);CHKERRQ(ierr); 464807ea322SDave May PetscFunctionReturn(0); 465807ea322SDave May } 466807ea322SDave May 467807ea322SDave May /*@ 468807ea322SDave May PetscViewerBinarySetSkipInfo - Binary file will not have .info file created with it 469807ea322SDave May 470807ea322SDave May Not Collective 471807ea322SDave May 472d8d19677SJose E. Roman Input Parameters: 473cc843e7aSLisandro Dalcin + viewer - PetscViewer context, obtained from PetscViewerCreate() 474cc843e7aSLisandro Dalcin - skip - PETSC_TRUE implies the .info file will not be generated 475807ea322SDave May 476807ea322SDave May Options Database Key: 47710699b91SBarry Smith . -viewer_binary_skip_info - true indicates do not generate .info file 478807ea322SDave May 479807ea322SDave May Level: advanced 480807ea322SDave May 481807ea322SDave May .seealso: PetscViewerBinaryOpen(), PetscViewerBinaryGetDescriptor(), PetscViewerBinarySetSkipOptions(), 482807ea322SDave May PetscViewerBinaryGetSkipOptions(), PetscViewerBinaryGetSkipInfo() 483807ea322SDave May @*/ 484807ea322SDave May PetscErrorCode PetscViewerBinarySetSkipInfo(PetscViewer viewer,PetscBool skip) 485807ea322SDave May { 486807ea322SDave May PetscErrorCode ierr; 487807ea322SDave May 488807ea322SDave May PetscFunctionBegin; 489cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 490cc843e7aSLisandro Dalcin PetscValidLogicalCollectiveBool(viewer,skip,2); 491cc843e7aSLisandro Dalcin ierr = PetscTryMethod(viewer,"PetscViewerBinarySetSkipInfo_C",(PetscViewer,PetscBool),(viewer,skip));CHKERRQ(ierr); 492807ea322SDave May PetscFunctionReturn(0); 493807ea322SDave May } 494807ea322SDave May 495cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinarySetSkipInfo_Binary(PetscViewer viewer,PetscBool skip) 496807ea322SDave May { 497807ea322SDave May PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 498807ea322SDave May 499807ea322SDave May PetscFunctionBegin; 500cc843e7aSLisandro Dalcin vbinary->skipinfo = skip; 501807ea322SDave May PetscFunctionReturn(0); 502807ea322SDave May } 503807ea322SDave May 504807ea322SDave May /*@ 505807ea322SDave May PetscViewerBinaryGetSkipInfo - check if viewer wrote a .info file 506807ea322SDave May 507807ea322SDave May Not Collective 508807ea322SDave May 509807ea322SDave May Input Parameter: 510807ea322SDave May . viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 511807ea322SDave May 512807ea322SDave May Output Parameter: 513807ea322SDave May . skip - PETSC_TRUE implies the .info file was not generated 514807ea322SDave May 515807ea322SDave May Level: advanced 516807ea322SDave May 51795452b02SPatrick Sanan Notes: 51895452b02SPatrick Sanan This must be called after PetscViewerSetType() 519807ea322SDave May 520807ea322SDave May .seealso: PetscViewerBinaryOpen(), PetscViewerBinaryGetDescriptor(), PetscViewerBinarySkipInfo(), 521807ea322SDave May PetscViewerBinarySetSkipOptions(), PetscViewerBinarySetSkipInfo() 522807ea322SDave May @*/ 523807ea322SDave May PetscErrorCode PetscViewerBinaryGetSkipInfo(PetscViewer viewer,PetscBool *skip) 524807ea322SDave May { 525807ea322SDave May PetscErrorCode ierr; 526807ea322SDave May 527807ea322SDave May PetscFunctionBegin; 528cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 529cc843e7aSLisandro Dalcin PetscValidBoolPointer(skip,2); 530807ea322SDave May ierr = PetscUseMethod(viewer,"PetscViewerBinaryGetSkipInfo_C",(PetscViewer,PetscBool*),(viewer,skip));CHKERRQ(ierr); 531807ea322SDave May PetscFunctionReturn(0); 532807ea322SDave May } 533807ea322SDave May 534cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinaryGetSkipInfo_Binary(PetscViewer viewer,PetscBool *skip) 535807ea322SDave May { 536807ea322SDave May PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 537807ea322SDave May 538807ea322SDave May PetscFunctionBegin; 539cc843e7aSLisandro Dalcin *skip = vbinary->skipinfo; 540807ea322SDave May PetscFunctionReturn(0); 541807ea322SDave May } 542807ea322SDave May 5435c6c1daeSBarry Smith /*@ 5445c6c1daeSBarry Smith PetscViewerBinarySetSkipOptions - do not use the PETSc options database when loading objects 5455c6c1daeSBarry Smith 5465c6c1daeSBarry Smith Not Collective 5475c6c1daeSBarry Smith 5485c6c1daeSBarry Smith Input Parameters: 5495c6c1daeSBarry Smith + viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 55010699b91SBarry Smith - skip - PETSC_TRUE means do not use the options from the options database 5515c6c1daeSBarry Smith 5525c6c1daeSBarry Smith Options Database Key: 55310699b91SBarry Smith . -viewer_binary_skip_options - true means do not use the options from the options database 5545c6c1daeSBarry Smith 5555c6c1daeSBarry Smith Level: advanced 5565c6c1daeSBarry Smith 55795452b02SPatrick Sanan Notes: 55895452b02SPatrick Sanan This must be called after PetscViewerSetType() 5595c6c1daeSBarry Smith 5605c6c1daeSBarry Smith .seealso: PetscViewerBinaryOpen(), PetscViewerBinaryGetDescriptor(), PetscViewerBinarySkipInfo(), 5615c6c1daeSBarry Smith PetscViewerBinaryGetSkipOptions() 5625c6c1daeSBarry Smith @*/ 5635c6c1daeSBarry Smith PetscErrorCode PetscViewerBinarySetSkipOptions(PetscViewer viewer,PetscBool skip) 5645c6c1daeSBarry Smith { 565807ea322SDave May PetscErrorCode ierr; 5665c6c1daeSBarry Smith 5675c6c1daeSBarry Smith PetscFunctionBegin; 568cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 569cc843e7aSLisandro Dalcin PetscValidLogicalCollectiveBool(viewer,skip,2); 570cc843e7aSLisandro Dalcin ierr = PetscTryMethod(viewer,"PetscViewerBinarySetSkipOptions_C",(PetscViewer,PetscBool),(viewer,skip));CHKERRQ(ierr); 571807ea322SDave May PetscFunctionReturn(0); 572807ea322SDave May } 573807ea322SDave May 574cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinarySetSkipOptions_Binary(PetscViewer viewer,PetscBool skip) 575807ea322SDave May { 576cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 577807ea322SDave May 578807ea322SDave May PetscFunctionBegin; 579cc843e7aSLisandro Dalcin vbinary->skipoptions = skip; 5805c6c1daeSBarry Smith PetscFunctionReturn(0); 5815c6c1daeSBarry Smith } 5825c6c1daeSBarry Smith 5835c6c1daeSBarry Smith /*@ 5845c6c1daeSBarry Smith PetscViewerBinaryGetSkipOptions - checks if viewer uses the PETSc options database when loading objects 5855c6c1daeSBarry Smith 5865c6c1daeSBarry Smith Not Collective 5875c6c1daeSBarry Smith 5885c6c1daeSBarry Smith Input Parameter: 5895c6c1daeSBarry Smith . viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 5905c6c1daeSBarry Smith 5915c6c1daeSBarry Smith Output Parameter: 5925c6c1daeSBarry Smith . skip - PETSC_TRUE means do not use 5935c6c1daeSBarry Smith 5945c6c1daeSBarry Smith Level: advanced 5955c6c1daeSBarry Smith 59695452b02SPatrick Sanan Notes: 59795452b02SPatrick Sanan This must be called after PetscViewerSetType() 5985c6c1daeSBarry Smith 5995c6c1daeSBarry Smith .seealso: PetscViewerBinaryOpen(), PetscViewerBinaryGetDescriptor(), PetscViewerBinarySkipInfo(), 6005c6c1daeSBarry Smith PetscViewerBinarySetSkipOptions() 6015c6c1daeSBarry Smith @*/ 6025c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryGetSkipOptions(PetscViewer viewer,PetscBool *skip) 6035c6c1daeSBarry Smith { 604807ea322SDave May PetscErrorCode ierr; 6055c6c1daeSBarry Smith 6065c6c1daeSBarry Smith PetscFunctionBegin; 607cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 608cc843e7aSLisandro Dalcin PetscValidBoolPointer(skip,2); 609807ea322SDave May ierr = PetscUseMethod(viewer,"PetscViewerBinaryGetSkipOptions_C",(PetscViewer,PetscBool*),(viewer,skip));CHKERRQ(ierr); 6105c6c1daeSBarry Smith PetscFunctionReturn(0); 6115c6c1daeSBarry Smith } 6125c6c1daeSBarry Smith 613cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinaryGetSkipOptions_Binary(PetscViewer viewer,PetscBool *skip) 6145c6c1daeSBarry Smith { 615d21b9a37SPierre Jolivet PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 6165c6c1daeSBarry Smith 6175c6c1daeSBarry Smith PetscFunctionBegin; 618cc843e7aSLisandro Dalcin *skip = vbinary->skipoptions; 6195c6c1daeSBarry Smith PetscFunctionReturn(0); 6205c6c1daeSBarry Smith } 6215c6c1daeSBarry Smith 6225c6c1daeSBarry Smith /*@ 6235c6c1daeSBarry Smith PetscViewerBinarySetSkipHeader - do not write a header with size information on output, just raw data 6245c6c1daeSBarry Smith 6255c6c1daeSBarry Smith Not Collective 6265c6c1daeSBarry Smith 6275c6c1daeSBarry Smith Input Parameters: 6285c6c1daeSBarry Smith + viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 6295c6c1daeSBarry Smith - skip - PETSC_TRUE means do not write header 6305c6c1daeSBarry Smith 6315c6c1daeSBarry Smith Options Database Key: 63210699b91SBarry Smith . -viewer_binary_skip_header - PETSC_TRUE means do not write header 6335c6c1daeSBarry Smith 6345c6c1daeSBarry Smith Level: advanced 6355c6c1daeSBarry Smith 63695452b02SPatrick Sanan Notes: 63795452b02SPatrick Sanan This must be called after PetscViewerSetType() 6385c6c1daeSBarry Smith 63910699b91SBarry Smith Is ignored on anything but a binary viewer 6405c6c1daeSBarry Smith 6415c6c1daeSBarry Smith .seealso: PetscViewerBinaryOpen(), PetscViewerBinaryGetDescriptor(), PetscViewerBinarySkipInfo(), 6425c6c1daeSBarry Smith PetscViewerBinaryGetSkipHeader() 6435c6c1daeSBarry Smith @*/ 6445c6c1daeSBarry Smith PetscErrorCode PetscViewerBinarySetSkipHeader(PetscViewer viewer,PetscBool skip) 6455c6c1daeSBarry Smith { 6465c6c1daeSBarry Smith PetscErrorCode ierr; 6475c6c1daeSBarry Smith 6485c6c1daeSBarry Smith PetscFunctionBegin; 649cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 650cc843e7aSLisandro Dalcin PetscValidLogicalCollectiveBool(viewer,skip,2); 651cc843e7aSLisandro Dalcin ierr = PetscTryMethod(viewer,"PetscViewerBinarySetSkipHeader_C",(PetscViewer,PetscBool),(viewer,skip));CHKERRQ(ierr); 6525c6c1daeSBarry Smith PetscFunctionReturn(0); 6535c6c1daeSBarry Smith } 6545c6c1daeSBarry Smith 655cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinarySetSkipHeader_Binary(PetscViewer viewer,PetscBool skip) 6565c6c1daeSBarry Smith { 6575c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 6585c6c1daeSBarry Smith 6595c6c1daeSBarry Smith PetscFunctionBegin; 660cc843e7aSLisandro Dalcin vbinary->skipheader = skip; 6615c6c1daeSBarry Smith PetscFunctionReturn(0); 6625c6c1daeSBarry Smith } 6635c6c1daeSBarry Smith 6645c6c1daeSBarry Smith /*@ 6655c6c1daeSBarry Smith PetscViewerBinaryGetSkipHeader - checks whether to write a header with size information on output, or just raw data 6665c6c1daeSBarry Smith 6675c6c1daeSBarry Smith Not Collective 6685c6c1daeSBarry Smith 6695c6c1daeSBarry Smith Input Parameter: 6705c6c1daeSBarry Smith . viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 6715c6c1daeSBarry Smith 6725c6c1daeSBarry Smith Output Parameter: 6735c6c1daeSBarry Smith . skip - PETSC_TRUE means do not write header 6745c6c1daeSBarry Smith 6755c6c1daeSBarry Smith Level: advanced 6765c6c1daeSBarry Smith 67795452b02SPatrick Sanan Notes: 67895452b02SPatrick Sanan This must be called after PetscViewerSetType() 6795c6c1daeSBarry Smith 6805c6c1daeSBarry Smith Returns false for PETSCSOCKETVIEWER, you cannot skip the header for it. 6815c6c1daeSBarry Smith 6825c6c1daeSBarry Smith .seealso: PetscViewerBinaryOpen(), PetscViewerBinaryGetDescriptor(), PetscViewerBinarySkipInfo(), 6835c6c1daeSBarry Smith PetscViewerBinarySetSkipHeader() 6845c6c1daeSBarry Smith @*/ 6855c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryGetSkipHeader(PetscViewer viewer,PetscBool *skip) 6865c6c1daeSBarry Smith { 6875c6c1daeSBarry Smith PetscErrorCode ierr; 6885c6c1daeSBarry Smith 6895c6c1daeSBarry Smith PetscFunctionBegin; 690cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 691cc843e7aSLisandro Dalcin PetscValidBoolPointer(skip,2); 692163d334eSBarry Smith ierr = PetscUseMethod(viewer,"PetscViewerBinaryGetSkipHeader_C",(PetscViewer,PetscBool*),(viewer,skip));CHKERRQ(ierr); 6935c6c1daeSBarry Smith PetscFunctionReturn(0); 6945c6c1daeSBarry Smith } 6955c6c1daeSBarry Smith 696cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinaryGetSkipHeader_Binary(PetscViewer viewer,PetscBool *skip) 6975c6c1daeSBarry Smith { 698f3b211e4SSatish Balay PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 6995c6c1daeSBarry Smith 7005c6c1daeSBarry Smith PetscFunctionBegin; 701cc843e7aSLisandro Dalcin *skip = vbinary->skipheader; 7025c6c1daeSBarry Smith PetscFunctionReturn(0); 7035c6c1daeSBarry Smith } 7045c6c1daeSBarry Smith 7055c6c1daeSBarry Smith /*@C 7065c6c1daeSBarry Smith PetscViewerBinaryGetInfoPointer - Extracts the file pointer for the ASCII 7075c6c1daeSBarry Smith info file associated with a binary file. 7085c6c1daeSBarry Smith 7095c6c1daeSBarry Smith Not Collective 7105c6c1daeSBarry Smith 7115c6c1daeSBarry Smith Input Parameter: 7125c6c1daeSBarry Smith . viewer - PetscViewer context, obtained from PetscViewerBinaryOpen() 7135c6c1daeSBarry Smith 7145c6c1daeSBarry Smith Output Parameter: 7150298fd71SBarry Smith . file - file pointer Always returns NULL if not a binary viewer 7165c6c1daeSBarry Smith 7175c6c1daeSBarry Smith Level: advanced 7185c6c1daeSBarry Smith 7195c6c1daeSBarry Smith Notes: 7205c6c1daeSBarry Smith For writable binary PetscViewers, the descriptor will only be valid for the 7215c6c1daeSBarry Smith first processor in the communicator that shares the PetscViewer. 7225c6c1daeSBarry Smith 7235c6c1daeSBarry Smith Fortran Note: 7245c6c1daeSBarry Smith This routine is not supported in Fortran. 7255c6c1daeSBarry Smith 7265c6c1daeSBarry Smith .seealso: PetscViewerBinaryOpen(),PetscViewerBinaryGetDescriptor() 7275c6c1daeSBarry Smith @*/ 7285c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryGetInfoPointer(PetscViewer viewer,FILE **file) 7295c6c1daeSBarry Smith { 7305c6c1daeSBarry Smith PetscErrorCode ierr; 7315c6c1daeSBarry Smith 7325c6c1daeSBarry Smith PetscFunctionBegin; 733cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 734cc843e7aSLisandro Dalcin PetscValidPointer(file,2); 7350298fd71SBarry Smith *file = NULL; 7365c6c1daeSBarry Smith ierr = PetscTryMethod(viewer,"PetscViewerBinaryGetInfoPointer_C",(PetscViewer,FILE **),(viewer,file));CHKERRQ(ierr); 7375c6c1daeSBarry Smith PetscFunctionReturn(0); 7385c6c1daeSBarry Smith } 7395c6c1daeSBarry Smith 740cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerBinaryGetInfoPointer_Binary(PetscViewer viewer,FILE **file) 7415c6c1daeSBarry Smith { 742cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 7435c6c1daeSBarry Smith PetscErrorCode ierr; 7445c6c1daeSBarry Smith 7455c6c1daeSBarry Smith PetscFunctionBegin; 746cc843e7aSLisandro Dalcin ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr); 747cc843e7aSLisandro Dalcin *file = vbinary->fdes_info; 748cc843e7aSLisandro Dalcin if (viewer->format == PETSC_VIEWER_BINARY_MATLAB && !vbinary->matlabheaderwritten) { 7495c6c1daeSBarry Smith if (vbinary->fdes_info) { 750cc843e7aSLisandro Dalcin FILE *info = vbinary->fdes_info; 751cc843e7aSLisandro Dalcin ierr = PetscFPrintf(PETSC_COMM_SELF,info,"#--- begin code written by PetscViewerBinary for MATLAB format ---#\n");CHKERRQ(ierr); 752cc843e7aSLisandro Dalcin ierr = PetscFPrintf(PETSC_COMM_SELF,info,"#$$ Set.filename = '%s';\n",vbinary->filename);CHKERRQ(ierr); 753cc843e7aSLisandro Dalcin ierr = PetscFPrintf(PETSC_COMM_SELF,info,"#$$ fd = PetscOpenFile(Set.filename);\n");CHKERRQ(ierr); 754cc843e7aSLisandro Dalcin ierr = PetscFPrintf(PETSC_COMM_SELF,info,"#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n");CHKERRQ(ierr); 755cc843e7aSLisandro Dalcin } 756cc843e7aSLisandro Dalcin vbinary->matlabheaderwritten = PETSC_TRUE; 7575c6c1daeSBarry Smith } 7585c6c1daeSBarry Smith PetscFunctionReturn(0); 7595c6c1daeSBarry Smith } 7605c6c1daeSBarry Smith 761e0385b85SDave May #if defined(PETSC_HAVE_MPIIO) 762e0385b85SDave May static PetscErrorCode PetscViewerFileClose_BinaryMPIIO(PetscViewer v) 763e0385b85SDave May { 764e0385b85SDave May PetscViewer_Binary *vbinary = (PetscViewer_Binary*)v->data; 765e0385b85SDave May PetscErrorCode ierr; 766e0385b85SDave May 767e0385b85SDave May PetscFunctionBegin; 768e8a65b7dSLisandro Dalcin if (vbinary->mfdes != MPI_FILE_NULL) { 769ffc4695bSBarry Smith ierr = MPI_File_close(&vbinary->mfdes);CHKERRMPI(ierr); 770e0385b85SDave May } 771e8a65b7dSLisandro Dalcin if (vbinary->mfsub != MPI_FILE_NULL) { 772ffc4695bSBarry Smith ierr = MPI_File_close(&vbinary->mfsub);CHKERRMPI(ierr); 773e8a65b7dSLisandro Dalcin } 774cc843e7aSLisandro Dalcin vbinary->moff = 0; 775e0385b85SDave May PetscFunctionReturn(0); 776e0385b85SDave May } 777e0385b85SDave May #endif 778e0385b85SDave May 779cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileClose_BinarySTDIO(PetscViewer v) 780cc843e7aSLisandro Dalcin { 781cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)v->data; 782cc843e7aSLisandro Dalcin PetscErrorCode ierr; 783cc843e7aSLisandro Dalcin 784cc843e7aSLisandro Dalcin PetscFunctionBegin; 785cc843e7aSLisandro Dalcin if (vbinary->fdes != -1) { 786cc843e7aSLisandro Dalcin ierr = PetscBinaryClose(vbinary->fdes);CHKERRQ(ierr); 787cc843e7aSLisandro Dalcin vbinary->fdes = -1; 788cc843e7aSLisandro Dalcin if (vbinary->storecompressed) { 789cc843e7aSLisandro Dalcin char cmd[8+PETSC_MAX_PATH_LEN],out[64+PETSC_MAX_PATH_LEN] = ""; 790cc843e7aSLisandro Dalcin const char *gzfilename = vbinary->ogzfilename ? vbinary->ogzfilename : vbinary->filename; 791cc843e7aSLisandro Dalcin /* compress the file */ 792cc843e7aSLisandro Dalcin ierr = PetscStrncpy(cmd,"gzip -f ",sizeof(cmd));CHKERRQ(ierr); 793cc843e7aSLisandro Dalcin ierr = PetscStrlcat(cmd,gzfilename,sizeof(cmd));CHKERRQ(ierr); 794cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_POPEN) 795cc843e7aSLisandro Dalcin { 796cc843e7aSLisandro Dalcin FILE *fp; 797cc843e7aSLisandro Dalcin ierr = PetscPOpen(PETSC_COMM_SELF,NULL,cmd,"r",&fp);CHKERRQ(ierr); 798*98921bdaSJacob Faibussowitsch if (fgets(out,(int)(sizeof(out)-1),fp)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error from command %s\n%s",cmd,out); 799cc843e7aSLisandro Dalcin ierr = PetscPClose(PETSC_COMM_SELF,fp);CHKERRQ(ierr); 800cc843e7aSLisandro Dalcin } 801cc843e7aSLisandro Dalcin #endif 802cc843e7aSLisandro Dalcin } 803cc843e7aSLisandro Dalcin } 804cc843e7aSLisandro Dalcin ierr = PetscFree(vbinary->ogzfilename);CHKERRQ(ierr); 805cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 806cc843e7aSLisandro Dalcin } 807cc843e7aSLisandro Dalcin 808cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileClose_BinaryInfo(PetscViewer v) 809cc843e7aSLisandro Dalcin { 810cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)v->data; 811cc843e7aSLisandro Dalcin PetscErrorCode ierr; 812cc843e7aSLisandro Dalcin 813cc843e7aSLisandro Dalcin PetscFunctionBegin; 814cc843e7aSLisandro Dalcin if (v->format == PETSC_VIEWER_BINARY_MATLAB && vbinary->matlabheaderwritten) { 815cc843e7aSLisandro Dalcin if (vbinary->fdes_info) { 816cc843e7aSLisandro Dalcin FILE *info = vbinary->fdes_info; 817cc843e7aSLisandro Dalcin ierr = PetscFPrintf(PETSC_COMM_SELF,info,"#--- begin code written by PetscViewerBinary for MATLAB format ---#\n");CHKERRQ(ierr); 818cc843e7aSLisandro Dalcin ierr = PetscFPrintf(PETSC_COMM_SELF,info,"#$$ close(fd);\n");CHKERRQ(ierr); 819cc843e7aSLisandro Dalcin ierr = PetscFPrintf(PETSC_COMM_SELF,info,"#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n");CHKERRQ(ierr); 820cc843e7aSLisandro Dalcin } 821cc843e7aSLisandro Dalcin } 822cc843e7aSLisandro Dalcin if (vbinary->fdes_info) { 823cc843e7aSLisandro Dalcin FILE *info = vbinary->fdes_info; 824cc843e7aSLisandro Dalcin vbinary->fdes_info = NULL; 825cc843e7aSLisandro Dalcin if (fclose(info)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file"); 826cc843e7aSLisandro Dalcin } 827cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 828cc843e7aSLisandro Dalcin } 829cc843e7aSLisandro Dalcin 830cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileClose_Binary(PetscViewer v) 831cc843e7aSLisandro Dalcin { 832cc843e7aSLisandro Dalcin PetscErrorCode ierr; 833cc843e7aSLisandro Dalcin 834cc843e7aSLisandro Dalcin PetscFunctionBegin; 835cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 836cc843e7aSLisandro Dalcin ierr = PetscViewerFileClose_BinaryMPIIO(v);CHKERRQ(ierr); 837cc843e7aSLisandro Dalcin #endif 838cc843e7aSLisandro Dalcin ierr = PetscViewerFileClose_BinarySTDIO(v);CHKERRQ(ierr); 839cc843e7aSLisandro Dalcin ierr = PetscViewerFileClose_BinaryInfo(v);CHKERRQ(ierr); 840cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 841cc843e7aSLisandro Dalcin } 842cc843e7aSLisandro Dalcin 84381f0254dSBarry Smith static PetscErrorCode PetscViewerDestroy_Binary(PetscViewer v) 8445c6c1daeSBarry Smith { 8455c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary*)v->data; 8465c6c1daeSBarry Smith PetscErrorCode ierr; 8475c6c1daeSBarry Smith 8485c6c1daeSBarry Smith PetscFunctionBegin; 849e0385b85SDave May ierr = PetscViewerFileClose_Binary(v);CHKERRQ(ierr); 850f90597f1SStefano Zampini ierr = PetscFree(vbinary->filename);CHKERRQ(ierr); 851e0385b85SDave May ierr = PetscFree(vbinary);CHKERRQ(ierr); 85222a8f86cSLisandro Dalcin 85322a8f86cSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetFlowControl_C",NULL);CHKERRQ(ierr); 85422a8f86cSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetFlowControl_C",NULL);CHKERRQ(ierr); 85522a8f86cSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetSkipHeader_C",NULL);CHKERRQ(ierr); 856cc843e7aSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetSkipHeader_C",NULL);CHKERRQ(ierr); 85722a8f86cSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetSkipOptions_C",NULL);CHKERRQ(ierr); 85822a8f86cSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetSkipOptions_C",NULL);CHKERRQ(ierr); 85922a8f86cSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetSkipInfo_C",NULL);CHKERRQ(ierr); 86022a8f86cSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetSkipInfo_C",NULL);CHKERRQ(ierr); 86122a8f86cSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetInfoPointer_C",NULL);CHKERRQ(ierr); 86222a8f86cSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr); 863cc843e7aSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr); 864cc843e7aSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",NULL);CHKERRQ(ierr); 865cc843e7aSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr); 86622a8f86cSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 86722a8f86cSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetUseMPIIO_C",NULL);CHKERRQ(ierr); 86822a8f86cSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetUseMPIIO_C",NULL);CHKERRQ(ierr); 86922a8f86cSLisandro Dalcin #endif 870e0385b85SDave May PetscFunctionReturn(0); 871e0385b85SDave May } 8725c6c1daeSBarry Smith 8735c6c1daeSBarry Smith /*@C 8745c6c1daeSBarry Smith PetscViewerBinaryOpen - Opens a file for binary input/output. 8755c6c1daeSBarry Smith 876d083f849SBarry Smith Collective 8775c6c1daeSBarry Smith 8785c6c1daeSBarry Smith Input Parameters: 8795c6c1daeSBarry Smith + comm - MPI communicator 8805c6c1daeSBarry Smith . name - name of file 881cc843e7aSLisandro Dalcin - mode - open mode of file 8825c6c1daeSBarry Smith $ FILE_MODE_WRITE - create new file for binary output 8835c6c1daeSBarry Smith $ FILE_MODE_READ - open existing file for binary input 8845c6c1daeSBarry Smith $ FILE_MODE_APPEND - open existing file for binary output 8855c6c1daeSBarry Smith 8865c6c1daeSBarry Smith Output Parameter: 887cc843e7aSLisandro Dalcin . viewer - PetscViewer for binary input/output to use with the specified file 8885c6c1daeSBarry Smith 8895c6c1daeSBarry Smith Options Database Keys: 89063c55180SPatrick Sanan + -viewer_binary_filename <name> - 89163c55180SPatrick Sanan . -viewer_binary_skip_info - 89263c55180SPatrick Sanan . -viewer_binary_skip_options - 89363c55180SPatrick Sanan . -viewer_binary_skip_header - 89463c55180SPatrick Sanan - -viewer_binary_mpiio - 8955c6c1daeSBarry Smith 8965c6c1daeSBarry Smith Level: beginner 8975c6c1daeSBarry Smith 8985c6c1daeSBarry Smith Note: 8995c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 9005c6c1daeSBarry Smith 9015c6c1daeSBarry Smith For reading files, the filename may begin with ftp:// or http:// and/or 9025c6c1daeSBarry Smith end with .gz; in this case file is brought over and uncompressed. 9035c6c1daeSBarry Smith 9045c6c1daeSBarry Smith For creating files, if the file name ends with .gz it is automatically 9055c6c1daeSBarry Smith compressed when closed. 9065c6c1daeSBarry Smith 9076a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), 9085c6c1daeSBarry Smith VecView(), MatView(), VecLoad(), MatLoad(), PetscViewerBinaryGetDescriptor(), 90967918a83SBarry Smith PetscViewerBinaryGetInfoPointer(), PetscFileMode, PetscViewer, PetscViewerBinaryRead(), PetscViewerBinarySetUseMPIIO(), 91067918a83SBarry Smith PetscViewerBinaryGetUseMPIIO(), PetscViewerBinaryGetMPIIOOffset() 9115c6c1daeSBarry Smith @*/ 912cc843e7aSLisandro Dalcin PetscErrorCode PetscViewerBinaryOpen(MPI_Comm comm,const char name[],PetscFileMode mode,PetscViewer *viewer) 9135c6c1daeSBarry Smith { 9145c6c1daeSBarry Smith PetscErrorCode ierr; 9155c6c1daeSBarry Smith 9165c6c1daeSBarry Smith PetscFunctionBegin; 917cc843e7aSLisandro Dalcin ierr = PetscViewerCreate(comm,viewer);CHKERRQ(ierr); 918cc843e7aSLisandro Dalcin ierr = PetscViewerSetType(*viewer,PETSCVIEWERBINARY);CHKERRQ(ierr); 919cc843e7aSLisandro Dalcin ierr = PetscViewerFileSetMode(*viewer,mode);CHKERRQ(ierr); 920cc843e7aSLisandro Dalcin ierr = PetscViewerFileSetName(*viewer,name);CHKERRQ(ierr); 921cc843e7aSLisandro Dalcin ierr = PetscViewerSetFromOptions(*viewer);CHKERRQ(ierr); 9225c6c1daeSBarry Smith PetscFunctionReturn(0); 9235c6c1daeSBarry Smith } 9245c6c1daeSBarry Smith 9255c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO) 926060da220SMatthew G. Knepley static PetscErrorCode PetscViewerBinaryWriteReadMPIIO(PetscViewer viewer,void *data,PetscInt num,PetscInt *count,PetscDataType dtype,PetscBool write) 9275c6c1daeSBarry Smith { 92822a8f86cSLisandro Dalcin MPI_Comm comm = PetscObjectComm((PetscObject)viewer); 9295c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 93022a8f86cSLisandro Dalcin MPI_File mfdes = vbinary->mfdes; 9315c6c1daeSBarry Smith PetscErrorCode ierr; 9325c6c1daeSBarry Smith MPI_Datatype mdtype; 93369e10bbaSLisandro Dalcin PetscMPIInt rank,cnt; 9345c6c1daeSBarry Smith MPI_Status status; 9355c6c1daeSBarry Smith MPI_Aint ul,dsize; 9365c6c1daeSBarry Smith 9375c6c1daeSBarry Smith PetscFunctionBegin; 938ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 939060da220SMatthew G. Knepley ierr = PetscMPIIntCast(num,&cnt);CHKERRQ(ierr); 9405c6c1daeSBarry Smith ierr = PetscDataTypeToMPIDataType(dtype,&mdtype);CHKERRQ(ierr); 9415c6c1daeSBarry Smith if (write) { 942dd400576SPatrick Sanan if (rank == 0) { 94369e10bbaSLisandro Dalcin ierr = MPIU_File_write_at(mfdes,vbinary->moff,data,cnt,mdtype,&status);CHKERRQ(ierr); 94469e10bbaSLisandro Dalcin } 9455c6c1daeSBarry Smith } else { 946dd400576SPatrick Sanan if (rank == 0) { 94769e10bbaSLisandro Dalcin ierr = MPIU_File_read_at(mfdes,vbinary->moff,data,cnt,mdtype,&status);CHKERRQ(ierr); 948ffc4695bSBarry Smith if (cnt > 0) {ierr = MPI_Get_count(&status,mdtype,&cnt);CHKERRMPI(ierr);} 9495c6c1daeSBarry Smith } 950ffc4695bSBarry Smith ierr = MPI_Bcast(&cnt,1,MPI_INT,0,comm);CHKERRMPI(ierr); 951ffc4695bSBarry Smith ierr = MPI_Bcast(data,cnt,mdtype,0,comm);CHKERRMPI(ierr); 95269e10bbaSLisandro Dalcin } 953ffc4695bSBarry Smith ierr = MPI_Type_get_extent(mdtype,&ul,&dsize);CHKERRMPI(ierr); 9545c6c1daeSBarry Smith vbinary->moff += dsize*cnt; 9559860990eSLisandro Dalcin if (count) *count = cnt; 9565c6c1daeSBarry Smith PetscFunctionReturn(0); 9575c6c1daeSBarry Smith } 9585c6c1daeSBarry Smith #endif 9595c6c1daeSBarry Smith 9605c6c1daeSBarry Smith /*@C 9615c6c1daeSBarry Smith PetscViewerBinaryRead - Reads from a binary file, all processors get the same result 9625c6c1daeSBarry Smith 963d083f849SBarry Smith Collective 9645c6c1daeSBarry Smith 9655c6c1daeSBarry Smith Input Parameters: 9665c6c1daeSBarry Smith + viewer - the binary viewer 967907376e6SBarry Smith . data - location of the data to be written 968060da220SMatthew G. Knepley . num - number of items of data to read 969907376e6SBarry Smith - dtype - type of data to read 9705c6c1daeSBarry Smith 971f8e4bde8SMatthew G. Knepley Output Parameters: 9725972f5f3SLisandro Dalcin . count - number of items of data actually read, or NULL. 973f8e4bde8SMatthew G. Knepley 9745c6c1daeSBarry Smith Level: beginner 9755c6c1daeSBarry Smith 9766a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), 9775c6c1daeSBarry Smith VecView(), MatView(), VecLoad(), MatLoad(), PetscViewerBinaryGetDescriptor(), 978d23dbae5SPatrick Sanan PetscViewerBinaryGetInfoPointer(), PetscFileMode, PetscViewer, PetscViewerBinaryRead() 9795c6c1daeSBarry Smith @*/ 980060da220SMatthew G. Knepley PetscErrorCode PetscViewerBinaryRead(PetscViewer viewer,void *data,PetscInt num,PetscInt *count,PetscDataType dtype) 9815c6c1daeSBarry Smith { 9825c6c1daeSBarry Smith PetscErrorCode ierr; 98322a8f86cSLisandro Dalcin PetscViewer_Binary *vbinary; 9845c6c1daeSBarry Smith 98522a8f86cSLisandro Dalcin PetscFunctionBegin; 98622a8f86cSLisandro Dalcin PetscValidHeaderSpecificType(viewer,PETSC_VIEWER_CLASSID,1,PETSCVIEWERBINARY); 98722a8f86cSLisandro Dalcin PetscValidLogicalCollectiveInt(viewer,num,3); 988c98fd787SBarry Smith ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr); 98922a8f86cSLisandro Dalcin vbinary = (PetscViewer_Binary*)viewer->data; 9905c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO) 991bc196f7cSDave May if (vbinary->usempiio) { 992060da220SMatthew G. Knepley ierr = PetscViewerBinaryWriteReadMPIIO(viewer,data,num,count,dtype,PETSC_FALSE);CHKERRQ(ierr); 9935c6c1daeSBarry Smith } else { 9945c6c1daeSBarry Smith #endif 9959860990eSLisandro Dalcin ierr = PetscBinarySynchronizedRead(PetscObjectComm((PetscObject)viewer),vbinary->fdes,data,num,count,dtype);CHKERRQ(ierr); 9965c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO) 9975c6c1daeSBarry Smith } 9985c6c1daeSBarry Smith #endif 9995c6c1daeSBarry Smith PetscFunctionReturn(0); 10005c6c1daeSBarry Smith } 10015c6c1daeSBarry Smith 10025c6c1daeSBarry Smith /*@C 10035c6c1daeSBarry Smith PetscViewerBinaryWrite - writes to a binary file, only from the first process 10045c6c1daeSBarry Smith 1005d083f849SBarry Smith Collective 10065c6c1daeSBarry Smith 10075c6c1daeSBarry Smith Input Parameters: 10085c6c1daeSBarry Smith + viewer - the binary viewer 10095c6c1daeSBarry Smith . data - location of data 10105824c9d0SPatrick Sanan . count - number of items of data to write 1011f253e43cSLisandro Dalcin - dtype - type of data to write 10125c6c1daeSBarry Smith 10135c6c1daeSBarry Smith Level: beginner 10145c6c1daeSBarry Smith 10156a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), 10165c6c1daeSBarry Smith VecView(), MatView(), VecLoad(), MatLoad(), PetscViewerBinaryGetDescriptor(), PetscDataType 1017d23dbae5SPatrick Sanan PetscViewerBinaryGetInfoPointer(), PetscFileMode, PetscViewer, PetscViewerBinaryRead() 10185c6c1daeSBarry Smith @*/ 1019f253e43cSLisandro Dalcin PetscErrorCode PetscViewerBinaryWrite(PetscViewer viewer,const void *data,PetscInt count,PetscDataType dtype) 10205c6c1daeSBarry Smith { 10215c6c1daeSBarry Smith PetscErrorCode ierr; 102222a8f86cSLisandro Dalcin PetscViewer_Binary *vbinary; 10235c6c1daeSBarry Smith 10245c6c1daeSBarry Smith PetscFunctionBegin; 102522a8f86cSLisandro Dalcin PetscValidHeaderSpecificType(viewer,PETSC_VIEWER_CLASSID,1,PETSCVIEWERBINARY); 102622a8f86cSLisandro Dalcin PetscValidLogicalCollectiveInt(viewer,count,3); 1027c98fd787SBarry Smith ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr); 102822a8f86cSLisandro Dalcin vbinary = (PetscViewer_Binary*)viewer->data; 10295c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO) 1030bc196f7cSDave May if (vbinary->usempiio) { 1031f253e43cSLisandro Dalcin ierr = PetscViewerBinaryWriteReadMPIIO(viewer,(void*)data,count,NULL,dtype,PETSC_TRUE);CHKERRQ(ierr); 10325c6c1daeSBarry Smith } else { 10335c6c1daeSBarry Smith #endif 1034f253e43cSLisandro Dalcin ierr = PetscBinarySynchronizedWrite(PetscObjectComm((PetscObject)viewer),vbinary->fdes,data,count,dtype);CHKERRQ(ierr); 10355c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO) 10365c6c1daeSBarry Smith } 10375c6c1daeSBarry Smith #endif 10385c6c1daeSBarry Smith PetscFunctionReturn(0); 10395c6c1daeSBarry Smith } 10405c6c1daeSBarry Smith 10415972f5f3SLisandro Dalcin static PetscErrorCode PetscViewerBinaryWriteReadAll(PetscViewer viewer,PetscBool write,void *data,PetscInt count,PetscInt start,PetscInt total,PetscDataType dtype) 10425972f5f3SLisandro Dalcin { 10435972f5f3SLisandro Dalcin MPI_Comm comm = PetscObjectComm((PetscObject)viewer); 10445972f5f3SLisandro Dalcin PetscMPIInt size,rank; 10455972f5f3SLisandro Dalcin MPI_Datatype mdtype; 1046ec4bef21SJose E. Roman PETSC_UNUSED MPI_Aint lb; 1047ec4bef21SJose E. Roman MPI_Aint dsize; 10485972f5f3SLisandro Dalcin PetscBool useMPIIO; 10495972f5f3SLisandro Dalcin PetscErrorCode ierr; 10505972f5f3SLisandro Dalcin 10515972f5f3SLisandro Dalcin PetscFunctionBegin; 10525972f5f3SLisandro Dalcin PetscValidHeaderSpecificType(viewer,PETSC_VIEWER_CLASSID,1,PETSCVIEWERBINARY); 1053064a246eSJacob Faibussowitsch PetscValidLogicalCollectiveBool(viewer,((start>=0)||(start==PETSC_DETERMINE)),5); 1054064a246eSJacob Faibussowitsch PetscValidLogicalCollectiveBool(viewer,((total>=0)||(total==PETSC_DETERMINE)),6); 1055064a246eSJacob Faibussowitsch PetscValidLogicalCollectiveInt(viewer,total,6); 10565972f5f3SLisandro Dalcin ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr); 10575972f5f3SLisandro Dalcin 10585972f5f3SLisandro Dalcin ierr = PetscDataTypeToMPIDataType(dtype,&mdtype);CHKERRQ(ierr); 1059ffc4695bSBarry Smith ierr = MPI_Type_get_extent(mdtype,&lb,&dsize);CHKERRMPI(ierr); 1060ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 1061ffc4695bSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 10625972f5f3SLisandro Dalcin 10635972f5f3SLisandro Dalcin ierr = PetscViewerBinaryGetUseMPIIO(viewer,&useMPIIO);CHKERRQ(ierr); 10645972f5f3SLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 10655972f5f3SLisandro Dalcin if (useMPIIO) { 10665972f5f3SLisandro Dalcin MPI_File mfdes; 10675972f5f3SLisandro Dalcin MPI_Offset off; 10685972f5f3SLisandro Dalcin PetscMPIInt cnt; 10695972f5f3SLisandro Dalcin 10705972f5f3SLisandro Dalcin if (start == PETSC_DETERMINE) { 1071ffc4695bSBarry Smith ierr = MPI_Scan(&count,&start,1,MPIU_INT,MPI_SUM,comm);CHKERRMPI(ierr); 10725972f5f3SLisandro Dalcin start -= count; 10735972f5f3SLisandro Dalcin } 10745972f5f3SLisandro Dalcin if (total == PETSC_DETERMINE) { 10755972f5f3SLisandro Dalcin total = start + count; 1076ffc4695bSBarry Smith ierr = MPI_Bcast(&total,1,MPIU_INT,size-1,comm);CHKERRMPI(ierr); 10775972f5f3SLisandro Dalcin } 10785972f5f3SLisandro Dalcin ierr = PetscMPIIntCast(count,&cnt);CHKERRQ(ierr); 10795972f5f3SLisandro Dalcin ierr = PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes);CHKERRQ(ierr); 10805972f5f3SLisandro Dalcin ierr = PetscViewerBinaryGetMPIIOOffset(viewer,&off);CHKERRQ(ierr); 10815972f5f3SLisandro Dalcin off += (MPI_Offset)(start*dsize); 10825972f5f3SLisandro Dalcin if (write) { 10835972f5f3SLisandro Dalcin ierr = MPIU_File_write_at_all(mfdes,off,data,cnt,mdtype,MPI_STATUS_IGNORE);CHKERRQ(ierr); 10845972f5f3SLisandro Dalcin } else { 10855972f5f3SLisandro Dalcin ierr = MPIU_File_read_at_all(mfdes,off,data,cnt,mdtype,MPI_STATUS_IGNORE);CHKERRQ(ierr); 10865972f5f3SLisandro Dalcin } 10875972f5f3SLisandro Dalcin off = (MPI_Offset)(total*dsize); 10885972f5f3SLisandro Dalcin ierr = PetscViewerBinaryAddMPIIOOffset(viewer,off);CHKERRQ(ierr); 10895972f5f3SLisandro Dalcin PetscFunctionReturn(0); 10905972f5f3SLisandro Dalcin } 10915972f5f3SLisandro Dalcin #endif 10925972f5f3SLisandro Dalcin { 10935972f5f3SLisandro Dalcin int fdes; 10945972f5f3SLisandro Dalcin char *workbuf = NULL; 1095dd400576SPatrick Sanan PetscInt tcount = rank == 0 ? 0 : count,maxcount=0,message_count,flowcontrolcount; 10965972f5f3SLisandro Dalcin PetscMPIInt tag,cnt,maxcnt,scnt=0,rcnt=0,j; 10975972f5f3SLisandro Dalcin MPI_Status status; 10985972f5f3SLisandro Dalcin 10995972f5f3SLisandro Dalcin ierr = PetscCommGetNewTag(comm,&tag);CHKERRQ(ierr); 1100ffc4695bSBarry Smith ierr = MPI_Reduce(&tcount,&maxcount,1,MPIU_INT,MPI_MAX,0,comm);CHKERRMPI(ierr); 11015972f5f3SLisandro Dalcin ierr = PetscMPIIntCast(maxcount,&maxcnt);CHKERRQ(ierr); 11025972f5f3SLisandro Dalcin ierr = PetscMPIIntCast(count,&cnt);CHKERRQ(ierr); 11035972f5f3SLisandro Dalcin 11045972f5f3SLisandro Dalcin ierr = PetscViewerBinaryGetDescriptor(viewer,&fdes);CHKERRQ(ierr); 11055972f5f3SLisandro Dalcin ierr = PetscViewerFlowControlStart(viewer,&message_count,&flowcontrolcount);CHKERRQ(ierr); 1106dd400576SPatrick Sanan if (rank == 0) { 11075972f5f3SLisandro Dalcin ierr = PetscMalloc(maxcnt*dsize,&workbuf);CHKERRQ(ierr); 11085972f5f3SLisandro Dalcin if (write) { 11095972f5f3SLisandro Dalcin ierr = PetscBinaryWrite(fdes,data,cnt,dtype);CHKERRQ(ierr); 11105972f5f3SLisandro Dalcin } else { 11115972f5f3SLisandro Dalcin ierr = PetscBinaryRead(fdes,data,cnt,NULL,dtype);CHKERRQ(ierr); 11125972f5f3SLisandro Dalcin } 11135972f5f3SLisandro Dalcin for (j=1; j<size; j++) { 11149dddd249SSatish Balay ierr = PetscViewerFlowControlStepMain(viewer,j,&message_count,flowcontrolcount);CHKERRQ(ierr); 11155972f5f3SLisandro Dalcin if (write) { 1116ffc4695bSBarry Smith ierr = MPI_Recv(workbuf,maxcnt,mdtype,j,tag,comm,&status);CHKERRMPI(ierr); 1117ffc4695bSBarry Smith ierr = MPI_Get_count(&status,mdtype,&rcnt);CHKERRMPI(ierr); 11185972f5f3SLisandro Dalcin ierr = PetscBinaryWrite(fdes,workbuf,rcnt,dtype);CHKERRQ(ierr); 11195972f5f3SLisandro Dalcin } else { 1120ffc4695bSBarry Smith ierr = MPI_Recv(&scnt,1,MPI_INT,j,tag,comm,MPI_STATUS_IGNORE);CHKERRMPI(ierr); 11215972f5f3SLisandro Dalcin ierr = PetscBinaryRead(fdes,workbuf,scnt,NULL,dtype);CHKERRQ(ierr); 1122ffc4695bSBarry Smith ierr = MPI_Send(workbuf,scnt,mdtype,j,tag,comm);CHKERRMPI(ierr); 11235972f5f3SLisandro Dalcin } 11245972f5f3SLisandro Dalcin } 11255972f5f3SLisandro Dalcin ierr = PetscFree(workbuf);CHKERRQ(ierr); 11269dddd249SSatish Balay ierr = PetscViewerFlowControlEndMain(viewer,&message_count);CHKERRQ(ierr); 11275972f5f3SLisandro Dalcin } else { 11285972f5f3SLisandro Dalcin ierr = PetscViewerFlowControlStepWorker(viewer,rank,&message_count);CHKERRQ(ierr); 11295972f5f3SLisandro Dalcin if (write) { 1130ffc4695bSBarry Smith ierr = MPI_Send(data,cnt,mdtype,0,tag,comm);CHKERRMPI(ierr); 11315972f5f3SLisandro Dalcin } else { 1132ffc4695bSBarry Smith ierr = MPI_Send(&cnt,1,MPI_INT,0,tag,comm);CHKERRMPI(ierr); 1133ffc4695bSBarry Smith ierr = MPI_Recv(data,cnt,mdtype,0,tag,comm,MPI_STATUS_IGNORE);CHKERRMPI(ierr); 11345972f5f3SLisandro Dalcin } 11355972f5f3SLisandro Dalcin ierr = PetscViewerFlowControlEndWorker(viewer,&message_count);CHKERRQ(ierr); 11365972f5f3SLisandro Dalcin } 11375972f5f3SLisandro Dalcin } 11385972f5f3SLisandro Dalcin PetscFunctionReturn(0); 11395972f5f3SLisandro Dalcin } 11405972f5f3SLisandro Dalcin 11415972f5f3SLisandro Dalcin /*@C 11425972f5f3SLisandro Dalcin PetscViewerBinaryReadAll - reads from a binary file from all processes 11435972f5f3SLisandro Dalcin 11445972f5f3SLisandro Dalcin Collective 11455972f5f3SLisandro Dalcin 11465972f5f3SLisandro Dalcin Input Parameters: 11475972f5f3SLisandro Dalcin + viewer - the binary viewer 11485972f5f3SLisandro Dalcin . data - location of data 11495972f5f3SLisandro Dalcin . count - local number of items of data to read 11505972f5f3SLisandro Dalcin . start - local start, can be PETSC_DETERMINE 11515972f5f3SLisandro Dalcin . total - global number of items of data to read, can be PETSC_DETERMINE 11525972f5f3SLisandro Dalcin - dtype - type of data to read 11535972f5f3SLisandro Dalcin 11545972f5f3SLisandro Dalcin Level: advanced 11555972f5f3SLisandro Dalcin 1156d23dbae5SPatrick Sanan .seealso: PetscViewerBinaryOpen(), PetscViewerBinarySetUseMPIIO(), PetscViewerBinaryRead(), PetscViewerBinaryWriteAll() 11575972f5f3SLisandro Dalcin @*/ 11585972f5f3SLisandro Dalcin PetscErrorCode PetscViewerBinaryReadAll(PetscViewer viewer,void *data,PetscInt count,PetscInt start,PetscInt total,PetscDataType dtype) 11595972f5f3SLisandro Dalcin { 11605972f5f3SLisandro Dalcin PetscErrorCode ierr; 11615972f5f3SLisandro Dalcin PetscFunctionBegin; 11625972f5f3SLisandro Dalcin ierr = PetscViewerBinaryWriteReadAll(viewer,PETSC_FALSE,data,count,start,total,dtype);CHKERRQ(ierr); 11635972f5f3SLisandro Dalcin PetscFunctionReturn(0); 11645972f5f3SLisandro Dalcin } 11655972f5f3SLisandro Dalcin 11665972f5f3SLisandro Dalcin /*@C 11675972f5f3SLisandro Dalcin PetscViewerBinaryWriteAll - writes to a binary file from all processes 11685972f5f3SLisandro Dalcin 11695972f5f3SLisandro Dalcin Collective 11705972f5f3SLisandro Dalcin 11715972f5f3SLisandro Dalcin Input Parameters: 11725972f5f3SLisandro Dalcin + viewer - the binary viewer 11735972f5f3SLisandro Dalcin . data - location of data 11745972f5f3SLisandro Dalcin . count - local number of items of data to write 11755972f5f3SLisandro Dalcin . start - local start, can be PETSC_DETERMINE 11765972f5f3SLisandro Dalcin . total - global number of items of data to write, can be PETSC_DETERMINE 11775972f5f3SLisandro Dalcin - dtype - type of data to write 11785972f5f3SLisandro Dalcin 11795972f5f3SLisandro Dalcin Level: advanced 11805972f5f3SLisandro Dalcin 1181d23dbae5SPatrick Sanan .seealso: PetscViewerBinaryOpen(), PetscViewerBinarySetUseMPIIO(), PetscViewerBinaryWriteAll(), PetscViewerBinaryReadAll() 11825972f5f3SLisandro Dalcin @*/ 11835972f5f3SLisandro Dalcin PetscErrorCode PetscViewerBinaryWriteAll(PetscViewer viewer,const void *data,PetscInt count,PetscInt start,PetscInt total,PetscDataType dtype) 11845972f5f3SLisandro Dalcin { 11855972f5f3SLisandro Dalcin PetscErrorCode ierr; 11865972f5f3SLisandro Dalcin PetscFunctionBegin; 11875972f5f3SLisandro Dalcin ierr = PetscViewerBinaryWriteReadAll(viewer,PETSC_TRUE,(void*)data,count,start,total,dtype);CHKERRQ(ierr); 11885972f5f3SLisandro Dalcin PetscFunctionReturn(0); 11895972f5f3SLisandro Dalcin } 11905972f5f3SLisandro Dalcin 11915c6c1daeSBarry Smith /*@C 11925c6c1daeSBarry Smith PetscViewerBinaryWriteStringArray - writes to a binary file, only from the first process an array of strings 11935c6c1daeSBarry Smith 1194d083f849SBarry Smith Collective 11955c6c1daeSBarry Smith 11965c6c1daeSBarry Smith Input Parameters: 11975c6c1daeSBarry Smith + viewer - the binary viewer 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 12056a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), 12065c6c1daeSBarry Smith VecView(), MatView(), VecLoad(), MatLoad(), PetscViewerBinaryGetDescriptor(), 1207d23dbae5SPatrick Sanan PetscViewerBinaryGetInfoPointer(), PetscFileMode, PetscViewer, PetscViewerBinaryRead() 12085c6c1daeSBarry Smith @*/ 120978fbdcc8SBarry Smith PetscErrorCode PetscViewerBinaryWriteStringArray(PetscViewer viewer,const char * const *data) 12105c6c1daeSBarry Smith { 12115c6c1daeSBarry Smith PetscErrorCode ierr; 12125c6c1daeSBarry Smith PetscInt i,n = 0,*sizes; 1213cc843e7aSLisandro Dalcin size_t len; 12145c6c1daeSBarry Smith 121522a8f86cSLisandro Dalcin PetscFunctionBegin; 1216c98fd787SBarry Smith ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr); 12175c6c1daeSBarry Smith /* count number of strings */ 12185c6c1daeSBarry Smith while (data[n++]); 12195c6c1daeSBarry Smith n--; 1220854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&sizes);CHKERRQ(ierr); 12215c6c1daeSBarry Smith sizes[0] = n; 12225c6c1daeSBarry Smith for (i=0; i<n; i++) { 1223cc843e7aSLisandro Dalcin ierr = PetscStrlen(data[i],&len);CHKERRQ(ierr); 1224cc843e7aSLisandro Dalcin sizes[i+1] = (PetscInt)len + 1; /* size includes space for the null terminator */ 12255c6c1daeSBarry Smith } 1226f253e43cSLisandro Dalcin ierr = PetscViewerBinaryWrite(viewer,sizes,n+1,PETSC_INT);CHKERRQ(ierr); 12275c6c1daeSBarry Smith for (i=0; i<n; i++) { 1228f253e43cSLisandro Dalcin ierr = PetscViewerBinaryWrite(viewer,(void*)data[i],sizes[i+1],PETSC_CHAR);CHKERRQ(ierr); 12295c6c1daeSBarry Smith } 12305c6c1daeSBarry Smith ierr = PetscFree(sizes);CHKERRQ(ierr); 12315c6c1daeSBarry Smith PetscFunctionReturn(0); 12325c6c1daeSBarry Smith } 12335c6c1daeSBarry Smith 12345c6c1daeSBarry Smith /*@C 12355c6c1daeSBarry Smith PetscViewerBinaryReadStringArray - reads a binary file an array of strings 12365c6c1daeSBarry Smith 1237d083f849SBarry Smith Collective 12385c6c1daeSBarry Smith 12395c6c1daeSBarry Smith Input Parameter: 12405c6c1daeSBarry Smith . viewer - the binary viewer 12415c6c1daeSBarry Smith 12425c6c1daeSBarry Smith Output Parameter: 12435c6c1daeSBarry Smith . data - location of the array of strings 12445c6c1daeSBarry Smith 12455c6c1daeSBarry Smith Level: intermediate 12465c6c1daeSBarry Smith 124795452b02SPatrick Sanan Notes: 124895452b02SPatrick Sanan array of strings is null terminated 12495c6c1daeSBarry Smith 12506a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), 12515c6c1daeSBarry Smith VecView(), MatView(), VecLoad(), MatLoad(), PetscViewerBinaryGetDescriptor(), 1252d23dbae5SPatrick Sanan PetscViewerBinaryGetInfoPointer(), PetscFileMode, PetscViewer, PetscViewerBinaryRead() 12535c6c1daeSBarry Smith @*/ 12545c6c1daeSBarry Smith PetscErrorCode PetscViewerBinaryReadStringArray(PetscViewer viewer,char ***data) 12555c6c1daeSBarry Smith { 12565c6c1daeSBarry Smith PetscErrorCode ierr; 1257060da220SMatthew G. Knepley PetscInt i,n,*sizes,N = 0; 12585c6c1daeSBarry Smith 125922a8f86cSLisandro Dalcin PetscFunctionBegin; 1260c98fd787SBarry Smith ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr); 12615c6c1daeSBarry Smith /* count number of strings */ 1262060da220SMatthew G. Knepley ierr = PetscViewerBinaryRead(viewer,&n,1,NULL,PETSC_INT);CHKERRQ(ierr); 1263785e854fSJed Brown ierr = PetscMalloc1(n,&sizes);CHKERRQ(ierr); 1264060da220SMatthew G. Knepley ierr = PetscViewerBinaryRead(viewer,sizes,n,NULL,PETSC_INT);CHKERRQ(ierr); 1265a297a907SKarl Rupp for (i=0; i<n; i++) N += sizes[i]; 1266854ce69bSBarry Smith ierr = PetscMalloc((n+1)*sizeof(char*) + N*sizeof(char),data);CHKERRQ(ierr); 12675c6c1daeSBarry Smith (*data)[0] = (char*)((*data) + n + 1); 1268a297a907SKarl Rupp for (i=1; i<n; i++) (*data)[i] = (*data)[i-1] + sizes[i-1]; 1269060da220SMatthew G. Knepley ierr = PetscViewerBinaryRead(viewer,(*data)[0],N,NULL,PETSC_CHAR);CHKERRQ(ierr); 127002c9f0b5SLisandro Dalcin (*data)[n] = NULL; 12715c6c1daeSBarry Smith ierr = PetscFree(sizes);CHKERRQ(ierr); 12725c6c1daeSBarry Smith PetscFunctionReturn(0); 12735c6c1daeSBarry Smith } 12745c6c1daeSBarry Smith 1275cc843e7aSLisandro Dalcin /*@C 1276cc843e7aSLisandro Dalcin PetscViewerFileSetMode - Sets the open mode of file 1277cc843e7aSLisandro Dalcin 1278cc843e7aSLisandro Dalcin Logically Collective on PetscViewer 1279cc843e7aSLisandro Dalcin 1280cc843e7aSLisandro Dalcin Input Parameters: 1281cc843e7aSLisandro Dalcin + viewer - the PetscViewer; must be a a PETSCVIEWERBINARY, PETSCVIEWERMATLAB, PETSCVIEWERHDF5, or PETSCVIEWERASCII PetscViewer 1282cc843e7aSLisandro Dalcin - mode - open mode of file 1283cc843e7aSLisandro Dalcin $ FILE_MODE_WRITE - create new file for output 1284cc843e7aSLisandro Dalcin $ FILE_MODE_READ - open existing file for input 1285cc843e7aSLisandro Dalcin $ FILE_MODE_APPEND - open existing file for output 1286cc843e7aSLisandro Dalcin 1287cc843e7aSLisandro Dalcin Level: advanced 1288cc843e7aSLisandro Dalcin 1289cc843e7aSLisandro Dalcin .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 1290cc843e7aSLisandro Dalcin 1291cc843e7aSLisandro Dalcin @*/ 1292cc843e7aSLisandro Dalcin PetscErrorCode PetscViewerFileSetMode(PetscViewer viewer,PetscFileMode mode) 1293cc843e7aSLisandro Dalcin { 1294cc843e7aSLisandro Dalcin PetscErrorCode ierr; 1295cc843e7aSLisandro Dalcin 1296cc843e7aSLisandro Dalcin PetscFunctionBegin; 1297cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1298cc843e7aSLisandro Dalcin PetscValidLogicalCollectiveEnum(viewer,mode,2); 12997e4fd573SVaclav Hapla if (mode == FILE_MODE_UNDEFINED) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Cannot set FILE_MODE_UNDEFINED"); 1300*98921bdaSJacob Faibussowitsch else if (mode < FILE_MODE_UNDEFINED || mode > FILE_MODE_APPEND_UPDATE) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_OUTOFRANGE,"Invalid file mode %d",(int)mode); 1301cc843e7aSLisandro Dalcin ierr = PetscTryMethod(viewer,"PetscViewerFileSetMode_C",(PetscViewer,PetscFileMode),(viewer,mode));CHKERRQ(ierr); 1302cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 1303cc843e7aSLisandro Dalcin } 1304cc843e7aSLisandro Dalcin 1305cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileSetMode_Binary(PetscViewer viewer,PetscFileMode mode) 1306cc843e7aSLisandro Dalcin { 1307cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 1308cc843e7aSLisandro Dalcin 1309cc843e7aSLisandro Dalcin PetscFunctionBegin; 1310*98921bdaSJacob Faibussowitsch if (viewer->setupcalled && vbinary->filemode != mode) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER,"Cannot change mode to %s after setup",PetscFileModes[mode]); 1311cc843e7aSLisandro Dalcin vbinary->filemode = mode; 1312cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 1313cc843e7aSLisandro Dalcin } 1314cc843e7aSLisandro Dalcin 1315cc843e7aSLisandro Dalcin /*@C 1316cc843e7aSLisandro Dalcin PetscViewerFileGetMode - Gets the open mode of file 1317cc843e7aSLisandro Dalcin 1318cc843e7aSLisandro Dalcin Not Collective 1319cc843e7aSLisandro Dalcin 1320cc843e7aSLisandro Dalcin Input Parameter: 1321cc843e7aSLisandro Dalcin . viewer - the PetscViewer; must be a PETSCVIEWERBINARY, PETSCVIEWERMATLAB, PETSCVIEWERHDF5, or PETSCVIEWERASCII PetscViewer 1322cc843e7aSLisandro Dalcin 1323cc843e7aSLisandro Dalcin Output Parameter: 1324cc843e7aSLisandro Dalcin . mode - open mode of file 1325cc843e7aSLisandro Dalcin $ FILE_MODE_WRITE - create new file for binary output 1326cc843e7aSLisandro Dalcin $ FILE_MODE_READ - open existing file for binary input 1327cc843e7aSLisandro Dalcin $ FILE_MODE_APPEND - open existing file for binary output 1328cc843e7aSLisandro Dalcin 1329cc843e7aSLisandro Dalcin Level: advanced 1330cc843e7aSLisandro Dalcin 1331cc843e7aSLisandro Dalcin .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 1332cc843e7aSLisandro Dalcin 1333cc843e7aSLisandro Dalcin @*/ 1334cc843e7aSLisandro Dalcin PetscErrorCode PetscViewerFileGetMode(PetscViewer viewer,PetscFileMode *mode) 1335cc843e7aSLisandro Dalcin { 1336cc843e7aSLisandro Dalcin PetscErrorCode ierr; 1337cc843e7aSLisandro Dalcin 1338cc843e7aSLisandro Dalcin PetscFunctionBegin; 1339cc843e7aSLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1340cc843e7aSLisandro Dalcin PetscValidPointer(mode,2); 1341cc843e7aSLisandro Dalcin ierr = PetscUseMethod(viewer,"PetscViewerFileGetMode_C",(PetscViewer,PetscFileMode*),(viewer,mode));CHKERRQ(ierr); 1342cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 1343cc843e7aSLisandro Dalcin } 1344cc843e7aSLisandro Dalcin 1345cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileGetMode_Binary(PetscViewer viewer,PetscFileMode *mode) 1346cc843e7aSLisandro Dalcin { 1347cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 1348cc843e7aSLisandro Dalcin 1349cc843e7aSLisandro Dalcin PetscFunctionBegin; 1350cc843e7aSLisandro Dalcin *mode = vbinary->filemode; 1351cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 1352cc843e7aSLisandro Dalcin } 1353cc843e7aSLisandro Dalcin 1354cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileSetName_Binary(PetscViewer viewer,const char name[]) 1355cc843e7aSLisandro Dalcin { 1356cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 1357cc843e7aSLisandro Dalcin PetscErrorCode ierr; 1358cc843e7aSLisandro Dalcin 1359cc843e7aSLisandro Dalcin PetscFunctionBegin; 1360cc843e7aSLisandro Dalcin if (viewer->setupcalled && vbinary->filename) { 1361cc843e7aSLisandro Dalcin /* gzip can be run after the file with the previous filename has been closed */ 1362cc843e7aSLisandro Dalcin ierr = PetscFree(vbinary->ogzfilename);CHKERRQ(ierr); 1363cc843e7aSLisandro Dalcin ierr = PetscStrallocpy(vbinary->filename,&vbinary->ogzfilename);CHKERRQ(ierr); 1364cc843e7aSLisandro Dalcin } 1365cc843e7aSLisandro Dalcin ierr = PetscFree(vbinary->filename);CHKERRQ(ierr); 1366cc843e7aSLisandro Dalcin ierr = PetscStrallocpy(name,&vbinary->filename);CHKERRQ(ierr); 1367cc843e7aSLisandro Dalcin viewer->setupcalled = PETSC_FALSE; 1368cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 1369cc843e7aSLisandro Dalcin } 1370cc843e7aSLisandro Dalcin 137181f0254dSBarry Smith static PetscErrorCode PetscViewerFileGetName_Binary(PetscViewer viewer,const char **name) 13725c6c1daeSBarry Smith { 13735c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 13745c6c1daeSBarry Smith 13755c6c1daeSBarry Smith PetscFunctionBegin; 13765c6c1daeSBarry Smith *name = vbinary->filename; 13775c6c1daeSBarry Smith PetscFunctionReturn(0); 13785c6c1daeSBarry Smith } 13795c6c1daeSBarry Smith 13805c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO) 1381e0385b85SDave May static PetscErrorCode PetscViewerFileSetUp_BinaryMPIIO(PetscViewer viewer) 13825c6c1daeSBarry Smith { 13835c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 1384e8a65b7dSLisandro Dalcin int amode; 1385cc843e7aSLisandro Dalcin PetscErrorCode ierr; 13865c6c1daeSBarry Smith 13875c6c1daeSBarry Smith PetscFunctionBegin; 1388a297a907SKarl Rupp vbinary->storecompressed = PETSC_FALSE; 13895c6c1daeSBarry Smith 1390cc843e7aSLisandro Dalcin vbinary->moff = 0; 1391cc843e7aSLisandro Dalcin switch (vbinary->filemode) { 1392e8a65b7dSLisandro Dalcin case FILE_MODE_READ: amode = MPI_MODE_RDONLY; break; 1393e8a65b7dSLisandro Dalcin case FILE_MODE_WRITE: amode = MPI_MODE_WRONLY | MPI_MODE_CREATE; break; 139422a8f86cSLisandro Dalcin case FILE_MODE_APPEND: amode = MPI_MODE_WRONLY | MPI_MODE_CREATE | MPI_MODE_APPEND; break; 13957e4fd573SVaclav Hapla case FILE_MODE_UNDEFINED: SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerSetUp()"); 1396*98921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Unsupported file mode %s",PetscFileModes[vbinary->filemode]); 13975c6c1daeSBarry Smith } 1398ffc4695bSBarry Smith ierr = MPI_File_open(PetscObjectComm((PetscObject)viewer),vbinary->filename,amode,MPI_INFO_NULL,&vbinary->mfdes);CHKERRMPI(ierr); 139922a8f86cSLisandro Dalcin /* 140022a8f86cSLisandro Dalcin The MPI standard does not have MPI_MODE_TRUNCATE. We emulate this behavior by setting the file size to zero. 140122a8f86cSLisandro Dalcin */ 1402ffc4695bSBarry Smith if (vbinary->filemode == FILE_MODE_WRITE) {ierr = MPI_File_set_size(vbinary->mfdes,0);CHKERRMPI(ierr);} 140322a8f86cSLisandro Dalcin /* 140422a8f86cSLisandro Dalcin Initially, all processes view the file as a linear byte stream. Therefore, for files opened with MPI_MODE_APPEND, 140522a8f86cSLisandro Dalcin MPI_File_get_position[_shared](fh, &offset) returns the absolute byte position at the end of file. 140622a8f86cSLisandro Dalcin Otherwise, we would need to call MPI_File_get_byte_offset(fh, offset, &byte_offset) to convert 140722a8f86cSLisandro Dalcin the offset in etype units to an absolute byte position. 140822a8f86cSLisandro Dalcin */ 1409ffc4695bSBarry Smith if (vbinary->filemode == FILE_MODE_APPEND) {ierr = MPI_File_get_position(vbinary->mfdes,&vbinary->moff);CHKERRMPI(ierr);} 1410cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 1411cc843e7aSLisandro Dalcin } 1412cc843e7aSLisandro Dalcin #endif 14135c6c1daeSBarry Smith 1414cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileSetUp_BinarySTDIO(PetscViewer viewer) 1415cc843e7aSLisandro Dalcin { 1416cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 1417cc843e7aSLisandro Dalcin const char *fname; 1418cc843e7aSLisandro Dalcin char bname[PETSC_MAX_PATH_LEN],*gz; 1419cc843e7aSLisandro Dalcin PetscBool found; 1420cc843e7aSLisandro Dalcin PetscMPIInt rank; 1421cc843e7aSLisandro Dalcin PetscErrorCode ierr; 14225c6c1daeSBarry Smith 1423cc843e7aSLisandro Dalcin PetscFunctionBegin; 1424ffc4695bSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank);CHKERRMPI(ierr); 14255c6c1daeSBarry Smith 1426cc843e7aSLisandro Dalcin /* if file name ends in .gz strip that off and note user wants file compressed */ 1427cc843e7aSLisandro Dalcin vbinary->storecompressed = PETSC_FALSE; 1428cc843e7aSLisandro Dalcin if (vbinary->filemode == FILE_MODE_WRITE) { 1429cc843e7aSLisandro Dalcin ierr = PetscStrstr(vbinary->filename,".gz",&gz);CHKERRQ(ierr); 1430cc843e7aSLisandro Dalcin if (gz && gz[3] == 0) {*gz = 0; vbinary->storecompressed = PETSC_TRUE;} 1431cc843e7aSLisandro Dalcin } 1432cc843e7aSLisandro Dalcin #if !defined(PETSC_HAVE_POPEN) 1433cc843e7aSLisandro Dalcin if (vbinary->storecompressed) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP_SYS,"Cannot run gzip on this machine"); 1434cc843e7aSLisandro Dalcin #endif 1435cc843e7aSLisandro Dalcin 1436cc843e7aSLisandro Dalcin fname = vbinary->filename; 1437cc843e7aSLisandro Dalcin if (vbinary->filemode == FILE_MODE_READ) { /* possibly get the file from remote site or compressed file */ 1438cc843e7aSLisandro Dalcin ierr = PetscFileRetrieve(PetscObjectComm((PetscObject)viewer),fname,bname,PETSC_MAX_PATH_LEN,&found);CHKERRQ(ierr); 1439*98921bdaSJacob Faibussowitsch if (!found) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_OPEN,"Cannot locate file: %s",fname); 1440cc843e7aSLisandro Dalcin fname = bname; 14415c6c1daeSBarry Smith } 14425c6c1daeSBarry Smith 1443cc843e7aSLisandro Dalcin vbinary->fdes = -1; 1444dd400576SPatrick Sanan if (rank == 0) { /* only first processor opens file*/ 1445cc843e7aSLisandro Dalcin PetscFileMode mode = vbinary->filemode; 1446cc843e7aSLisandro Dalcin if (mode == FILE_MODE_APPEND) { 1447cc843e7aSLisandro Dalcin /* check if asked to append to a non-existing file */ 1448cc843e7aSLisandro Dalcin ierr = PetscTestFile(fname,'\0',&found);CHKERRQ(ierr); 1449cc843e7aSLisandro Dalcin if (!found) mode = FILE_MODE_WRITE; 1450cc843e7aSLisandro Dalcin } 1451cc843e7aSLisandro Dalcin ierr = PetscBinaryOpen(fname,mode,&vbinary->fdes);CHKERRQ(ierr); 1452cc843e7aSLisandro Dalcin } 1453cc843e7aSLisandro Dalcin PetscFunctionReturn(0); 1454cc843e7aSLisandro Dalcin } 1455cc843e7aSLisandro Dalcin 1456cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerFileSetUp_BinaryInfo(PetscViewer viewer) 1457cc843e7aSLisandro Dalcin { 1458cc843e7aSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 1459cc843e7aSLisandro Dalcin PetscMPIInt rank; 1460cc843e7aSLisandro Dalcin PetscBool found; 1461cc843e7aSLisandro Dalcin PetscErrorCode ierr; 1462cc843e7aSLisandro Dalcin 1463cc843e7aSLisandro Dalcin PetscFunctionBegin; 1464cc843e7aSLisandro Dalcin vbinary->fdes_info = NULL; 1465ffc4695bSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank);CHKERRMPI(ierr); 1466dd400576SPatrick Sanan if (!vbinary->skipinfo && (vbinary->filemode == FILE_MODE_READ || rank == 0)) { 1467cc843e7aSLisandro Dalcin char infoname[PETSC_MAX_PATH_LEN],iname[PETSC_MAX_PATH_LEN],*gz; 1468cc843e7aSLisandro Dalcin 1469cc843e7aSLisandro Dalcin ierr = PetscStrncpy(infoname,vbinary->filename,sizeof(infoname));CHKERRQ(ierr); 1470cc843e7aSLisandro Dalcin /* remove .gz if it ends file name */ 1471cc843e7aSLisandro Dalcin ierr = PetscStrstr(infoname,".gz",&gz);CHKERRQ(ierr); 1472cc843e7aSLisandro Dalcin if (gz && gz[3] == 0) *gz = 0; 1473cc843e7aSLisandro Dalcin 1474a126751eSBarry Smith ierr = PetscStrlcat(infoname,".info",sizeof(infoname));CHKERRQ(ierr); 1475cc843e7aSLisandro Dalcin if (vbinary->filemode == FILE_MODE_READ) { 14765c6c1daeSBarry Smith ierr = PetscFixFilename(infoname,iname);CHKERRQ(ierr); 1477ce94432eSBarry Smith ierr = PetscFileRetrieve(PetscObjectComm((PetscObject)viewer),iname,infoname,PETSC_MAX_PATH_LEN,&found);CHKERRQ(ierr); 1478cc843e7aSLisandro Dalcin if (found) {ierr = PetscOptionsInsertFile(PetscObjectComm((PetscObject)viewer),((PetscObject)viewer)->options,infoname,PETSC_FALSE);CHKERRQ(ierr);} 1479dd400576SPatrick Sanan } else if (rank == 0) { /* write or append */ 1480cc843e7aSLisandro Dalcin const char *omode = (vbinary->filemode == FILE_MODE_APPEND) ? "a" : "w"; 1481cc843e7aSLisandro Dalcin vbinary->fdes_info = fopen(infoname,omode); 1482*98921bdaSJacob Faibussowitsch if (!vbinary->fdes_info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open .info file %s for writing",infoname); 14835c6c1daeSBarry Smith } 14845c6c1daeSBarry Smith } 14855c6c1daeSBarry Smith PetscFunctionReturn(0); 14865c6c1daeSBarry Smith } 14875c6c1daeSBarry Smith 1488cc843e7aSLisandro Dalcin static PetscErrorCode PetscViewerSetUp_Binary(PetscViewer viewer) 14895c6c1daeSBarry Smith { 14905c6c1daeSBarry Smith PetscViewer_Binary *vbinary = (PetscViewer_Binary*)viewer->data; 1491cc843e7aSLisandro Dalcin PetscBool usempiio; 1492cc843e7aSLisandro Dalcin PetscErrorCode ierr; 1493cc843e7aSLisandro Dalcin 14945c6c1daeSBarry Smith PetscFunctionBegin; 1495cc843e7aSLisandro Dalcin if (!vbinary->setfromoptionscalled) {ierr = PetscViewerSetFromOptions(viewer);CHKERRQ(ierr);} 1496cc843e7aSLisandro Dalcin if (!vbinary->filename) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call PetscViewerFileSetName()"); 1497cc843e7aSLisandro Dalcin if (vbinary->filemode == (PetscFileMode)-1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call PetscViewerFileSetMode()"); 1498cc843e7aSLisandro Dalcin ierr = PetscViewerFileClose_Binary(viewer);CHKERRQ(ierr); 1499cc843e7aSLisandro Dalcin 1500cc843e7aSLisandro Dalcin ierr = PetscViewerBinaryGetUseMPIIO(viewer,&usempiio);CHKERRQ(ierr); 1501cc843e7aSLisandro Dalcin if (usempiio) { 1502cc843e7aSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 1503cc843e7aSLisandro Dalcin ierr = PetscViewerFileSetUp_BinaryMPIIO(viewer);CHKERRQ(ierr); 1504cc843e7aSLisandro Dalcin #endif 1505cc843e7aSLisandro Dalcin } else { 1506cc843e7aSLisandro Dalcin ierr = PetscViewerFileSetUp_BinarySTDIO(viewer);CHKERRQ(ierr); 1507cc843e7aSLisandro Dalcin } 1508cc843e7aSLisandro Dalcin ierr = PetscViewerFileSetUp_BinaryInfo(viewer);CHKERRQ(ierr); 1509cc843e7aSLisandro Dalcin 1510cc843e7aSLisandro Dalcin ierr = PetscLogObjectState((PetscObject)viewer,"File: %s",vbinary->filename);CHKERRQ(ierr); 15115c6c1daeSBarry Smith PetscFunctionReturn(0); 15125c6c1daeSBarry Smith } 15135c6c1daeSBarry Smith 151481f0254dSBarry Smith static PetscErrorCode PetscViewerView_Binary(PetscViewer v,PetscViewer viewer) 15152bf49c77SBarry Smith { 1516cb6ad94fSLisandro Dalcin PetscViewer_Binary *vbinary = (PetscViewer_Binary*)v->data; 1517cb6ad94fSLisandro Dalcin const char *fname = vbinary->filename ? vbinary->filename : "not yet set"; 1518cc843e7aSLisandro Dalcin const char *fmode = vbinary->filemode != (PetscFileMode) -1 ? PetscFileModes[vbinary->filemode] : "not yet set"; 1519cb6ad94fSLisandro Dalcin PetscBool usempiio; 1520cc843e7aSLisandro Dalcin PetscErrorCode ierr; 15212bf49c77SBarry Smith 15222bf49c77SBarry Smith PetscFunctionBegin; 1523cb6ad94fSLisandro Dalcin ierr = PetscViewerBinaryGetUseMPIIO(v,&usempiio);CHKERRQ(ierr); 1524cb6ad94fSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer,"Filename: %s\n",fname);CHKERRQ(ierr); 1525cb6ad94fSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer,"Mode: %s (%s)\n",fmode,usempiio ? "mpiio" : "stdio");CHKERRQ(ierr); 15262bf49c77SBarry Smith PetscFunctionReturn(0); 15272bf49c77SBarry Smith } 15282bf49c77SBarry Smith 152922a8f86cSLisandro Dalcin static PetscErrorCode PetscViewerSetFromOptions_Binary(PetscOptionItems *PetscOptionsObject,PetscViewer viewer) 153003a1d59fSDave May { 153122a8f86cSLisandro Dalcin PetscViewer_Binary *binary = (PetscViewer_Binary*)viewer->data; 1532d22fd6bcSDave May char defaultname[PETSC_MAX_PATH_LEN]; 153303a1d59fSDave May PetscBool flg; 1534cc843e7aSLisandro Dalcin PetscErrorCode ierr; 1535e0385b85SDave May 153603a1d59fSDave May PetscFunctionBegin; 153722a8f86cSLisandro Dalcin if (viewer->setupcalled) PetscFunctionReturn(0); 153803a1d59fSDave May ierr = PetscOptionsHead(PetscOptionsObject,"Binary PetscViewer Options");CHKERRQ(ierr); 1539d22fd6bcSDave May ierr = PetscSNPrintf(defaultname,PETSC_MAX_PATH_LEN-1,"binaryoutput");CHKERRQ(ierr); 1540589a23caSBarry Smith ierr = PetscOptionsString("-viewer_binary_filename","Specify filename","PetscViewerFileSetName",defaultname,defaultname,sizeof(defaultname),&flg);CHKERRQ(ierr); 154122a8f86cSLisandro Dalcin if (flg) { ierr = PetscViewerFileSetName_Binary(viewer,defaultname);CHKERRQ(ierr); } 154222a8f86cSLisandro Dalcin ierr = PetscOptionsBool("-viewer_binary_skip_info","Skip writing/reading .info file","PetscViewerBinarySetSkipInfo",binary->skipinfo,&binary->skipinfo,NULL);CHKERRQ(ierr); 1543cc843e7aSLisandro Dalcin ierr = PetscOptionsBool("-viewer_binary_skip_options","Skip parsing Vec/Mat load options","PetscViewerBinarySetSkipOptions",binary->skipoptions,&binary->skipoptions,NULL);CHKERRQ(ierr); 154422a8f86cSLisandro Dalcin ierr = PetscOptionsBool("-viewer_binary_skip_header","Skip writing/reading header information","PetscViewerBinarySetSkipHeader",binary->skipheader,&binary->skipheader,NULL);CHKERRQ(ierr); 154503a1d59fSDave May #if defined(PETSC_HAVE_MPIIO) 154622a8f86cSLisandro Dalcin ierr = PetscOptionsBool("-viewer_binary_mpiio","Use MPI-IO functionality to write/read binary file","PetscViewerBinarySetUseMPIIO",binary->usempiio,&binary->usempiio,NULL);CHKERRQ(ierr); 15475972f5f3SLisandro Dalcin #else 15485972f5f3SLisandro Dalcin ierr = PetscOptionsBool("-viewer_binary_mpiio","Use MPI-IO functionality to write/read binary file (NOT AVAILABLE)","PetscViewerBinarySetUseMPIIO",PETSC_FALSE,NULL,NULL);CHKERRQ(ierr); 154903a1d59fSDave May #endif 155003a1d59fSDave May ierr = PetscOptionsTail();CHKERRQ(ierr); 1551bc196f7cSDave May binary->setfromoptionscalled = PETSC_TRUE; 155203a1d59fSDave May PetscFunctionReturn(0); 155303a1d59fSDave May } 155403a1d59fSDave May 15558556b5ebSBarry Smith /*MC 15568556b5ebSBarry Smith PETSCVIEWERBINARY - A viewer that saves to binary files 15578556b5ebSBarry Smith 15588556b5ebSBarry Smith .seealso: PetscViewerBinaryOpen(), PETSC_VIEWER_STDOUT_(),PETSC_VIEWER_STDOUT_SELF, PETSC_VIEWER_STDOUT_WORLD, PetscViewerCreate(), PetscViewerASCIIOpen(), 15598556b5ebSBarry Smith PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB, PETSCVIEWERDRAW, 156067918a83SBarry Smith PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType(), 156167918a83SBarry Smith PetscViewerBinaryGetUseMPIIO(), PetscViewerBinarySetUseMPIIO() 15628556b5ebSBarry Smith 15631b266c99SBarry Smith Level: beginner 15641b266c99SBarry Smith 15658556b5ebSBarry Smith M*/ 15668556b5ebSBarry Smith 15678cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_Binary(PetscViewer v) 15685c6c1daeSBarry Smith { 15695c6c1daeSBarry Smith PetscErrorCode ierr; 15705c6c1daeSBarry Smith PetscViewer_Binary *vbinary; 15715c6c1daeSBarry Smith 15725c6c1daeSBarry Smith PetscFunctionBegin; 1573b00a9115SJed Brown ierr = PetscNewLog(v,&vbinary);CHKERRQ(ierr); 15745c6c1daeSBarry Smith v->data = (void*)vbinary; 1575cc843e7aSLisandro Dalcin 157603a1d59fSDave May v->ops->setfromoptions = PetscViewerSetFromOptions_Binary; 15775c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_Binary; 15782bf49c77SBarry Smith v->ops->view = PetscViewerView_Binary; 1579c98fd787SBarry Smith v->ops->setup = PetscViewerSetUp_Binary; 1580cc843e7aSLisandro Dalcin v->ops->flush = NULL; /* Should we support Flush() ? */ 1581cc843e7aSLisandro Dalcin v->ops->getsubviewer = PetscViewerGetSubViewer_Binary; 1582cc843e7aSLisandro Dalcin v->ops->restoresubviewer = PetscViewerRestoreSubViewer_Binary; 1583cc843e7aSLisandro Dalcin v->ops->read = PetscViewerBinaryRead; 1584cc843e7aSLisandro Dalcin 1585cc843e7aSLisandro Dalcin vbinary->fdes = -1; 1586e8a65b7dSLisandro Dalcin #if defined(PETSC_HAVE_MPIIO) 1587cc843e7aSLisandro Dalcin vbinary->usempiio = PETSC_FALSE; 1588e8a65b7dSLisandro Dalcin vbinary->mfdes = MPI_FILE_NULL; 1589e8a65b7dSLisandro Dalcin vbinary->mfsub = MPI_FILE_NULL; 1590e8a65b7dSLisandro Dalcin #endif 1591cc843e7aSLisandro Dalcin vbinary->filename = NULL; 15927e4fd573SVaclav Hapla vbinary->filemode = FILE_MODE_UNDEFINED; 159302c9f0b5SLisandro Dalcin vbinary->fdes_info = NULL; 15945c6c1daeSBarry Smith vbinary->skipinfo = PETSC_FALSE; 15955c6c1daeSBarry Smith vbinary->skipoptions = PETSC_TRUE; 15965c6c1daeSBarry Smith vbinary->skipheader = PETSC_FALSE; 15975c6c1daeSBarry Smith vbinary->storecompressed = PETSC_FALSE; 1598f90597f1SStefano Zampini vbinary->ogzfilename = NULL; 15995c6c1daeSBarry Smith vbinary->flowcontrol = 256; /* seems a good number for Cray XT-5 */ 16005c6c1daeSBarry Smith 1601cc843e7aSLisandro Dalcin vbinary->setfromoptionscalled = PETSC_FALSE; 1602cc843e7aSLisandro Dalcin 1603bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetFlowControl_C",PetscViewerBinaryGetFlowControl_Binary);CHKERRQ(ierr); 1604bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetFlowControl_C",PetscViewerBinarySetFlowControl_Binary);CHKERRQ(ierr); 1605bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetSkipHeader_C",PetscViewerBinaryGetSkipHeader_Binary);CHKERRQ(ierr); 1606cc843e7aSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetSkipHeader_C",PetscViewerBinarySetSkipHeader_Binary);CHKERRQ(ierr); 1607807ea322SDave May ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetSkipOptions_C",PetscViewerBinaryGetSkipOptions_Binary);CHKERRQ(ierr); 1608807ea322SDave May ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetSkipOptions_C",PetscViewerBinarySetSkipOptions_Binary);CHKERRQ(ierr); 1609807ea322SDave May ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetSkipInfo_C",PetscViewerBinaryGetSkipInfo_Binary);CHKERRQ(ierr); 1610807ea322SDave May ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetSkipInfo_C",PetscViewerBinarySetSkipInfo_Binary);CHKERRQ(ierr); 1611bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetInfoPointer_C",PetscViewerBinaryGetInfoPointer_Binary);CHKERRQ(ierr); 1612bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_Binary);CHKERRQ(ierr); 1613cc843e7aSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_Binary);CHKERRQ(ierr); 1614cc843e7aSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_Binary);CHKERRQ(ierr); 1615cc843e7aSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_Binary);CHKERRQ(ierr); 16165c6c1daeSBarry Smith #if defined(PETSC_HAVE_MPIIO) 1617bc196f7cSDave May ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinaryGetUseMPIIO_C",PetscViewerBinaryGetUseMPIIO_Binary);CHKERRQ(ierr); 1618bc196f7cSDave May ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerBinarySetUseMPIIO_C",PetscViewerBinarySetUseMPIIO_Binary);CHKERRQ(ierr); 16195c6c1daeSBarry Smith #endif 16205c6c1daeSBarry Smith PetscFunctionReturn(0); 16215c6c1daeSBarry Smith } 16225c6c1daeSBarry Smith 16235c6c1daeSBarry Smith /* ---------------------------------------------------------------------*/ 16245c6c1daeSBarry Smith /* 16255c6c1daeSBarry Smith The variable Petsc_Viewer_Binary_keyval is used to indicate an MPI attribute that 16265c6c1daeSBarry Smith is attached to a communicator, in this case the attribute is a PetscViewer. 16275c6c1daeSBarry Smith */ 1628d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_Binary_keyval = MPI_KEYVAL_INVALID; 16295c6c1daeSBarry Smith 16305c6c1daeSBarry Smith /*@C 16315c6c1daeSBarry Smith PETSC_VIEWER_BINARY_ - Creates a binary PetscViewer shared by all processors 16325c6c1daeSBarry Smith in a communicator. 16335c6c1daeSBarry Smith 1634d083f849SBarry Smith Collective 16355c6c1daeSBarry Smith 16365c6c1daeSBarry Smith Input Parameter: 16375c6c1daeSBarry Smith . comm - the MPI communicator to share the binary PetscViewer 16385c6c1daeSBarry Smith 16395c6c1daeSBarry Smith Level: intermediate 16405c6c1daeSBarry Smith 16415c6c1daeSBarry Smith Options Database Keys: 164210699b91SBarry Smith + -viewer_binary_filename <name> - filename in which to store the binary data, defaults to binaryoutput 164310699b91SBarry Smith . -viewer_binary_skip_info - true means do not create .info file for this viewer 164410699b91SBarry Smith . -viewer_binary_skip_options - true means do not use the options database for this viewer 164510699b91SBarry Smith . -viewer_binary_skip_header - true means do not store the usual header information in the binary file 164610699b91SBarry Smith - -viewer_binary_mpiio - true means use the file via MPI-IO, maybe faster for large files and many MPI ranks 16475c6c1daeSBarry Smith 16485c6c1daeSBarry Smith Environmental variables: 164910699b91SBarry Smith - PETSC_VIEWER_BINARY_FILENAME - filename in which to store the binary data, defaults to binaryoutput 16505c6c1daeSBarry Smith 16515c6c1daeSBarry Smith Notes: 16525c6c1daeSBarry Smith Unlike almost all other PETSc routines, PETSC_VIEWER_BINARY_ does not return 16535c6c1daeSBarry Smith an error code. The binary PetscViewer is usually used in the form 16545c6c1daeSBarry Smith $ XXXView(XXX object,PETSC_VIEWER_BINARY_(comm)); 16555c6c1daeSBarry Smith 16565c6c1daeSBarry Smith .seealso: PETSC_VIEWER_BINARY_WORLD, PETSC_VIEWER_BINARY_SELF, PetscViewerBinaryOpen(), PetscViewerCreate(), 16575c6c1daeSBarry Smith PetscViewerDestroy() 16585c6c1daeSBarry Smith @*/ 16595c6c1daeSBarry Smith PetscViewer PETSC_VIEWER_BINARY_(MPI_Comm comm) 16605c6c1daeSBarry Smith { 16615c6c1daeSBarry Smith PetscErrorCode ierr; 16625c6c1daeSBarry Smith PetscBool flg; 16635c6c1daeSBarry Smith PetscViewer viewer; 16645c6c1daeSBarry Smith char fname[PETSC_MAX_PATH_LEN]; 16655c6c1daeSBarry Smith MPI_Comm ncomm; 16665c6c1daeSBarry Smith 16675c6c1daeSBarry Smith PetscFunctionBegin; 166802c9f0b5SLisandro 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);} 16695c6c1daeSBarry Smith if (Petsc_Viewer_Binary_keyval == MPI_KEYVAL_INVALID) { 167002c9f0b5SLisandro Dalcin ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_Binary_keyval,NULL); 167102c9f0b5SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 16725c6c1daeSBarry Smith } 167347435625SJed Brown ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_Binary_keyval,(void**)&viewer,(int*)&flg); 167402c9f0b5SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 16755c6c1daeSBarry Smith if (!flg) { /* PetscViewer not yet created */ 16765c6c1daeSBarry Smith ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_BINARY_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 16772cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 16785c6c1daeSBarry Smith if (!flg) { 16795c6c1daeSBarry Smith ierr = PetscStrcpy(fname,"binaryoutput"); 16802cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 16815c6c1daeSBarry Smith } 16825c6c1daeSBarry Smith ierr = PetscViewerBinaryOpen(ncomm,fname,FILE_MODE_WRITE,&viewer); 16832cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 16845c6c1daeSBarry Smith ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 16852cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 168647435625SJed Brown ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_Binary_keyval,(void*)viewer); 168702c9f0b5SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 16885c6c1daeSBarry Smith } 16895c6c1daeSBarry Smith ierr = PetscCommDestroy(&ncomm); 16902cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_BINARY_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 16915c6c1daeSBarry Smith PetscFunctionReturn(viewer); 16925c6c1daeSBarry Smith } 1693