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