1235f7792SMatthew G. Knepley #include <petscis.h> /*I "petscis.h" I*/ 2af0996ceSBarry Smith #include <petsc/private/isimpl.h> 390e28553SVaclav Hapla #include <petsc/private/viewerimpl.h> 420e823e8SBarry Smith #include <petsclayouthdf5.h> 5235f7792SMatthew G. Knepley 6*9371c9d4SSatish Balay PetscErrorCode ISView_Binary(IS is, PetscViewer viewer) { 7f92d9ac3SLisandro Dalcin PetscBool skipHeader; 8f92d9ac3SLisandro Dalcin PetscLayout map; 9f92d9ac3SLisandro Dalcin PetscInt tr[2], n, s, N; 10f92d9ac3SLisandro Dalcin const PetscInt *iarray; 11f92d9ac3SLisandro Dalcin 12f92d9ac3SLisandro Dalcin PetscFunctionBegin; 139566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer)); 149566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &skipHeader)); 15f92d9ac3SLisandro Dalcin 169566063dSJacob Faibussowitsch PetscCall(ISGetLayout(is, &map)); 179566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(map, &n)); 189566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(map, &s, NULL)); 199566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(map, &N)); 20f92d9ac3SLisandro Dalcin 21f92d9ac3SLisandro Dalcin /* write IS header */ 22*9371c9d4SSatish Balay tr[0] = IS_FILE_CLASSID; 23*9371c9d4SSatish Balay tr[1] = N; 249566063dSJacob Faibussowitsch if (!skipHeader) PetscCall(PetscViewerBinaryWrite(viewer, tr, 2, PETSC_INT)); 25f92d9ac3SLisandro Dalcin 26f92d9ac3SLisandro Dalcin /* write IS indices */ 279566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &iarray)); 289566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWriteAll(viewer, iarray, n, s, N, PETSC_INT)); 299566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &iarray)); 30f92d9ac3SLisandro Dalcin PetscFunctionReturn(0); 31f92d9ac3SLisandro Dalcin } 32f92d9ac3SLisandro Dalcin 33235f7792SMatthew G. Knepley #if defined(PETSC_HAVE_HDF5) 34235f7792SMatthew G. Knepley /* 35235f7792SMatthew G. Knepley This should handle properly the cases where PetscInt is 32 or 64 and hsize_t is 32 or 64. These means properly casting with 36235f7792SMatthew G. Knepley checks back and forth between the two types of variables. 37235f7792SMatthew G. Knepley */ 38*9371c9d4SSatish Balay PetscErrorCode ISLoad_HDF5(IS is, PetscViewer viewer) { 39235f7792SMatthew G. Knepley hid_t inttype; /* int type (H5T_NATIVE_INT or H5T_NATIVE_LLONG) */ 4090e28553SVaclav Hapla PetscInt *ind; 41235f7792SMatthew G. Knepley const char *isname; 42235f7792SMatthew G. Knepley 43235f7792SMatthew G. Knepley PetscFunctionBegin; 44fa7c0e95SVaclav Hapla PetscCheck(((PetscObject)is)->name, PetscObjectComm((PetscObject)is), PETSC_ERR_SUP, "IS name must be given using PetscObjectSetName() before ISLoad() since HDF5 can store multiple objects in a single file"); 45235f7792SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 46235f7792SMatthew G. Knepley inttype = H5T_NATIVE_LLONG; 47235f7792SMatthew G. Knepley #else 48235f7792SMatthew G. Knepley inttype = H5T_NATIVE_INT; 49235f7792SMatthew G. Knepley #endif 509566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)is, &isname)); 519566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Load(viewer, isname, is->map, inttype, (void **)&ind)); 529566063dSJacob Faibussowitsch PetscCall(ISGeneralSetIndices(is, is->map->n, ind, PETSC_OWN_POINTER)); 53235f7792SMatthew G. Knepley PetscFunctionReturn(0); 54235f7792SMatthew G. Knepley } 55235f7792SMatthew G. Knepley #endif 56235f7792SMatthew G. Knepley 57*9371c9d4SSatish Balay PetscErrorCode ISLoad_Binary(IS is, PetscViewer viewer) { 58f92d9ac3SLisandro Dalcin PetscBool isgeneral, skipHeader; 59f92d9ac3SLisandro Dalcin PetscInt tr[2], rows, N, n, s, *idx; 60f92d9ac3SLisandro Dalcin PetscLayout map; 61ede126feSBarry Smith 62ede126feSBarry Smith PetscFunctionBegin; 639566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)is, ISGENERAL, &isgeneral)); 64fa7c0e95SVaclav Hapla PetscCheck(isgeneral, PetscObjectComm((PetscObject)is), PETSC_ERR_ARG_INCOMP, "IS must be of type ISGENERAL to load into it"); 659566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer)); 669566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &skipHeader)); 67ede126feSBarry Smith 689566063dSJacob Faibussowitsch PetscCall(ISGetLayout(is, &map)); 699566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(map, &N)); 70f92d9ac3SLisandro Dalcin 71f92d9ac3SLisandro Dalcin /* read IS header */ 72f92d9ac3SLisandro Dalcin if (!skipHeader) { 739566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, tr, 2, NULL, PETSC_INT)); 74fa7c0e95SVaclav Hapla PetscCheck(tr[0] == IS_FILE_CLASSID, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Not an IS next in file"); 75fa7c0e95SVaclav Hapla PetscCheck(tr[1] >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "IS size (%" PetscInt_FMT ") in file is negative", tr[1]); 76c9cc58a2SBarry Smith PetscCheck(N < 0 || N == tr[1], PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "IS in file different size (%" PetscInt_FMT ") than input IS (%" PetscInt_FMT ")", tr[1], N); 77f92d9ac3SLisandro Dalcin rows = tr[1]; 78ede126feSBarry Smith } else { 79fa7c0e95SVaclav Hapla PetscCheck(N >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "IS binary file header was skipped, thus the user must specify the global size of input IS"); 80f92d9ac3SLisandro Dalcin rows = N; 81ede126feSBarry Smith } 82f92d9ac3SLisandro Dalcin 83f92d9ac3SLisandro Dalcin /* set IS size if not already set */ 849566063dSJacob Faibussowitsch if (N < 0) PetscCall(PetscLayoutSetSize(map, rows)); 859566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(map)); 86f92d9ac3SLisandro Dalcin 87f92d9ac3SLisandro Dalcin /* get IS sizes and check global size */ 889566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(map, &N)); 899566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(map, &n)); 909566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(map, &s, NULL)); 91fa7c0e95SVaclav Hapla PetscCheck(N == rows, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "IS in file different size (%" PetscInt_FMT ") than input IS (%" PetscInt_FMT ")", rows, N); 92f92d9ac3SLisandro Dalcin 93f92d9ac3SLisandro Dalcin /* read IS indices */ 949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &idx)); 959566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryReadAll(viewer, idx, n, s, N, PETSC_INT)); 969566063dSJacob Faibussowitsch PetscCall(ISGeneralSetIndices(is, n, idx, PETSC_OWN_POINTER)); 97ede126feSBarry Smith PetscFunctionReturn(0); 98ede126feSBarry Smith } 99ede126feSBarry Smith 100*9371c9d4SSatish Balay PetscErrorCode ISLoad_Default(IS is, PetscViewer viewer) { 101f92d9ac3SLisandro Dalcin PetscBool isbinary, ishdf5; 102235f7792SMatthew G. Knepley 103235f7792SMatthew G. Knepley PetscFunctionBegin; 1049566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 1059566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5)); 106235f7792SMatthew G. Knepley if (isbinary) { 1079566063dSJacob Faibussowitsch PetscCall(ISLoad_Binary(is, viewer)); 108235f7792SMatthew G. Knepley } else if (ishdf5) { 109f92d9ac3SLisandro Dalcin #if defined(PETSC_HAVE_HDF5) 1109566063dSJacob Faibussowitsch PetscCall(ISLoad_HDF5(is, viewer)); 111235f7792SMatthew G. Knepley #endif 112235f7792SMatthew G. Knepley } 113235f7792SMatthew G. Knepley PetscFunctionReturn(0); 114235f7792SMatthew G. Knepley } 115