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 6d71ae5a4SJacob Faibussowitsch PetscErrorCode ISView_Binary(IS is, PetscViewer viewer) 7d71ae5a4SJacob Faibussowitsch { 8f92d9ac3SLisandro Dalcin PetscBool skipHeader; 9f92d9ac3SLisandro Dalcin PetscLayout map; 10f92d9ac3SLisandro Dalcin PetscInt tr[2], n, s, N; 11f92d9ac3SLisandro Dalcin const PetscInt *iarray; 12f92d9ac3SLisandro Dalcin 13f92d9ac3SLisandro Dalcin PetscFunctionBegin; 149566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer)); 159566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &skipHeader)); 16f92d9ac3SLisandro Dalcin 179566063dSJacob Faibussowitsch PetscCall(ISGetLayout(is, &map)); 189566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(map, &n)); 199566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(map, &s, NULL)); 209566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(map, &N)); 21f92d9ac3SLisandro Dalcin 22f92d9ac3SLisandro Dalcin /* write IS header */ 239371c9d4SSatish Balay tr[0] = IS_FILE_CLASSID; 249371c9d4SSatish Balay tr[1] = N; 259566063dSJacob Faibussowitsch if (!skipHeader) PetscCall(PetscViewerBinaryWrite(viewer, tr, 2, PETSC_INT)); 26f92d9ac3SLisandro Dalcin 27f92d9ac3SLisandro Dalcin /* write IS indices */ 289566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &iarray)); 299566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWriteAll(viewer, iarray, n, s, N, PETSC_INT)); 309566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &iarray)); 313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 32f92d9ac3SLisandro Dalcin } 33f92d9ac3SLisandro Dalcin 34235f7792SMatthew G. Knepley #if defined(PETSC_HAVE_HDF5) 35235f7792SMatthew G. Knepley /* 36*822ddb90SBarry Smith This should handle properly the cases where PetscInt is 32 or 64 and hsize_t is 32 or 64. That is properly casting with 37235f7792SMatthew G. Knepley checks back and forth between the two types of variables. 38235f7792SMatthew G. Knepley */ 3966976f2fSJacob Faibussowitsch static PetscErrorCode ISLoad_HDF5(IS is, PetscViewer viewer) 40d71ae5a4SJacob Faibussowitsch { 4190e28553SVaclav Hapla PetscInt *ind; 42235f7792SMatthew G. Knepley const char *isname; 43235f7792SMatthew G. Knepley 44235f7792SMatthew G. Knepley PetscFunctionBegin; 45fa7c0e95SVaclav 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"); 469566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)is, &isname)); 471dc74096SMartin Diehl #if defined(PETSC_USE_64BIT_INDICES) 481dc74096SMartin Diehl PetscCall(PetscViewerHDF5Load(viewer, isname, is->map, H5T_NATIVE_LLONG, (void **)&ind)); 491dc74096SMartin Diehl #else 501dc74096SMartin Diehl PetscCall(PetscViewerHDF5Load(viewer, isname, is->map, H5T_NATIVE_INT, (void **)&ind)); 511dc74096SMartin Diehl #endif 529566063dSJacob Faibussowitsch PetscCall(ISGeneralSetIndices(is, is->map->n, ind, PETSC_OWN_POINTER)); 5321c42226SMatthew G. Knepley PetscCall(PetscInfo(is, "Read IS object with name %s of size %" PetscInt_FMT ":%" PetscInt_FMT "\n", isname, is->map->n, is->map->N)); 543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 55235f7792SMatthew G. Knepley } 56235f7792SMatthew G. Knepley #endif 57235f7792SMatthew G. Knepley 58da8c939bSJacob Faibussowitsch static PetscErrorCode ISLoad_Binary(IS is, PetscViewer viewer) 59d71ae5a4SJacob Faibussowitsch { 60f92d9ac3SLisandro Dalcin PetscBool isgeneral, skipHeader; 61f92d9ac3SLisandro Dalcin PetscInt tr[2], rows, N, n, s, *idx; 62f92d9ac3SLisandro Dalcin PetscLayout map; 63ede126feSBarry Smith 64ede126feSBarry Smith PetscFunctionBegin; 659566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)is, ISGENERAL, &isgeneral)); 66fa7c0e95SVaclav Hapla PetscCheck(isgeneral, PetscObjectComm((PetscObject)is), PETSC_ERR_ARG_INCOMP, "IS must be of type ISGENERAL to load into it"); 679566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer)); 689566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &skipHeader)); 69ede126feSBarry Smith 709566063dSJacob Faibussowitsch PetscCall(ISGetLayout(is, &map)); 719566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(map, &N)); 72f92d9ac3SLisandro Dalcin 73f92d9ac3SLisandro Dalcin /* read IS header */ 74f92d9ac3SLisandro Dalcin if (!skipHeader) { 759566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, tr, 2, NULL, PETSC_INT)); 76fa7c0e95SVaclav Hapla PetscCheck(tr[0] == IS_FILE_CLASSID, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Not an IS next in file"); 77fa7c0e95SVaclav Hapla PetscCheck(tr[1] >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "IS size (%" PetscInt_FMT ") in file is negative", tr[1]); 78c9cc58a2SBarry 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); 79f92d9ac3SLisandro Dalcin rows = tr[1]; 80ede126feSBarry Smith } else { 81fa7c0e95SVaclav 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"); 82f92d9ac3SLisandro Dalcin rows = N; 83ede126feSBarry Smith } 84f92d9ac3SLisandro Dalcin 85f92d9ac3SLisandro Dalcin /* set IS size if not already set */ 869566063dSJacob Faibussowitsch if (N < 0) PetscCall(PetscLayoutSetSize(map, rows)); 879566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(map)); 88f92d9ac3SLisandro Dalcin 89f92d9ac3SLisandro Dalcin /* get IS sizes and check global size */ 909566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(map, &N)); 919566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(map, &n)); 929566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(map, &s, NULL)); 93fa7c0e95SVaclav 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); 94f92d9ac3SLisandro Dalcin 95f92d9ac3SLisandro Dalcin /* read IS indices */ 969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &idx)); 979566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryReadAll(viewer, idx, n, s, N, PETSC_INT)); 989566063dSJacob Faibussowitsch PetscCall(ISGeneralSetIndices(is, n, idx, PETSC_OWN_POINTER)); 993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 100ede126feSBarry Smith } 101ede126feSBarry Smith 102d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLoad_Default(IS is, PetscViewer viewer) 103d71ae5a4SJacob Faibussowitsch { 104f92d9ac3SLisandro Dalcin PetscBool isbinary, ishdf5; 105235f7792SMatthew G. Knepley 106235f7792SMatthew G. Knepley PetscFunctionBegin; 1079566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 1089566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5)); 109235f7792SMatthew G. Knepley if (isbinary) { 1109566063dSJacob Faibussowitsch PetscCall(ISLoad_Binary(is, viewer)); 111235f7792SMatthew G. Knepley } else if (ishdf5) { 112f92d9ac3SLisandro Dalcin #if defined(PETSC_HAVE_HDF5) 1139566063dSJacob Faibussowitsch PetscCall(ISLoad_HDF5(is, viewer)); 114235f7792SMatthew G. Knepley #endif 115235f7792SMatthew G. Knepley } 1163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 117235f7792SMatthew G. Knepley } 118