xref: /petsc/src/vec/is/section/interface/hdf5/sectionhdf5.c (revision 8112c1cbf372cb53bf7c5aca994a84a6a303db4d)
12e1d0745SJose E. Roman #include <petsc/private/sectionimpl.h> /*I "petscsection.h" I*/
22e1d0745SJose E. Roman #include <petscsf.h>
32e1d0745SJose E. Roman #include <petscis.h>
42e1d0745SJose E. Roman #include <petscviewerhdf5.h>
52e1d0745SJose E. Roman #include <petsclayouthdf5.h>
62e1d0745SJose E. Roman 
PetscSectionView_HDF5_SingleField(PetscSection s,PetscViewer viewer)72e1d0745SJose E. Roman static PetscErrorCode PetscSectionView_HDF5_SingleField(PetscSection s, PetscViewer viewer)
82e1d0745SJose E. Roman {
92e1d0745SJose E. Roman   MPI_Comm  comm;
102e1d0745SJose E. Roman   PetscInt  pStart, pEnd, p, n;
112e1d0745SJose E. Roman   PetscBool hasConstraints, includesConstraints;
122e1d0745SJose E. Roman   IS        dofIS, offIS, cdofIS, coffIS, cindIS;
132e1d0745SJose E. Roman   PetscInt *dofs, *offs, *cdofs, *coffs, *cinds, dof, cdof, m, moff, i;
142e1d0745SJose E. Roman 
152e1d0745SJose E. Roman   PetscFunctionBegin;
162e1d0745SJose E. Roman   PetscCall(PetscObjectGetComm((PetscObject)s, &comm));
172e1d0745SJose E. Roman   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
182e1d0745SJose E. Roman   hasConstraints = (s->bc) ? PETSC_TRUE : PETSC_FALSE;
19*5440e5dcSBarry Smith   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &hasConstraints, 1, MPI_C_BOOL, MPI_LOR, comm));
202e1d0745SJose E. Roman   for (p = pStart, n = 0, m = 0; p < pEnd; ++p) {
212e1d0745SJose E. Roman     PetscCall(PetscSectionGetDof(s, p, &dof));
222e1d0745SJose E. Roman     if (dof >= 0) {
232e1d0745SJose E. Roman       if (hasConstraints) {
242e1d0745SJose E. Roman         PetscCall(PetscSectionGetConstraintDof(s, p, &cdof));
252e1d0745SJose E. Roman         m += cdof;
262e1d0745SJose E. Roman       }
272e1d0745SJose E. Roman       n++;
282e1d0745SJose E. Roman     }
292e1d0745SJose E. Roman   }
302e1d0745SJose E. Roman   PetscCall(PetscMalloc1(n, &dofs));
312e1d0745SJose E. Roman   PetscCall(PetscMalloc1(n, &offs));
322e1d0745SJose E. Roman   if (hasConstraints) {
332e1d0745SJose E. Roman     PetscCall(PetscMalloc1(n, &cdofs));
342e1d0745SJose E. Roman     PetscCall(PetscMalloc1(n, &coffs));
352e1d0745SJose E. Roman     PetscCall(PetscMalloc1(m, &cinds));
362e1d0745SJose E. Roman   }
372e1d0745SJose E. Roman   for (p = pStart, n = 0, m = 0; p < pEnd; ++p) {
382e1d0745SJose E. Roman     PetscCall(PetscSectionGetDof(s, p, &dof));
392e1d0745SJose E. Roman     if (dof >= 0) {
402e1d0745SJose E. Roman       dofs[n] = dof;
412e1d0745SJose E. Roman       PetscCall(PetscSectionGetOffset(s, p, &offs[n]));
422e1d0745SJose E. Roman       if (hasConstraints) {
432e1d0745SJose E. Roman         const PetscInt *cpinds;
442e1d0745SJose E. Roman 
452e1d0745SJose E. Roman         PetscCall(PetscSectionGetConstraintDof(s, p, &cdof));
462e1d0745SJose E. Roman         PetscCall(PetscSectionGetConstraintIndices(s, p, &cpinds));
472e1d0745SJose E. Roman         cdofs[n] = cdof;
482e1d0745SJose E. Roman         coffs[n] = m;
492e1d0745SJose E. Roman         for (i = 0; i < cdof; ++i) cinds[m++] = cpinds[i];
502e1d0745SJose E. Roman       }
512e1d0745SJose E. Roman       n++;
522e1d0745SJose E. Roman     }
532e1d0745SJose E. Roman   }
542e1d0745SJose E. Roman   if (hasConstraints) {
552e1d0745SJose E. Roman     PetscCallMPI(MPI_Scan(&m, &moff, 1, MPIU_INT, MPI_SUM, comm));
562e1d0745SJose E. Roman     moff -= m;
572e1d0745SJose E. Roman     for (p = 0; p < n; ++p) coffs[p] += moff;
582e1d0745SJose E. Roman   }
592e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "hasConstraints", PETSC_BOOL, (void *)&hasConstraints));
602e1d0745SJose E. Roman   PetscCall(PetscSectionGetIncludesConstraints(s, &includesConstraints));
612e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "includesConstraints", PETSC_BOOL, (void *)&includesConstraints));
622e1d0745SJose E. Roman   PetscCall(ISCreateGeneral(comm, n, dofs, PETSC_OWN_POINTER, &dofIS));
632e1d0745SJose E. Roman   PetscCall(PetscObjectSetName((PetscObject)dofIS, "atlasDof"));
642e1d0745SJose E. Roman   PetscCall(ISView(dofIS, viewer));
652e1d0745SJose E. Roman   PetscCall(ISDestroy(&dofIS));
662e1d0745SJose E. Roman   PetscCall(ISCreateGeneral(comm, n, offs, PETSC_OWN_POINTER, &offIS));
672e1d0745SJose E. Roman   PetscCall(PetscObjectSetName((PetscObject)offIS, "atlasOff"));
682e1d0745SJose E. Roman   PetscCall(ISView(offIS, viewer));
692e1d0745SJose E. Roman   PetscCall(ISDestroy(&offIS));
702e1d0745SJose E. Roman   if (hasConstraints) {
712e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5PushGroup(viewer, "bc"));
722e1d0745SJose E. Roman     PetscCall(ISCreateGeneral(comm, n, cdofs, PETSC_OWN_POINTER, &cdofIS));
732e1d0745SJose E. Roman     PetscCall(PetscObjectSetName((PetscObject)cdofIS, "atlasDof"));
742e1d0745SJose E. Roman     PetscCall(ISView(cdofIS, viewer));
752e1d0745SJose E. Roman     PetscCall(ISDestroy(&cdofIS));
762e1d0745SJose E. Roman     PetscCall(ISCreateGeneral(comm, n, coffs, PETSC_OWN_POINTER, &coffIS));
772e1d0745SJose E. Roman     PetscCall(PetscObjectSetName((PetscObject)coffIS, "atlasOff"));
782e1d0745SJose E. Roman     PetscCall(ISView(coffIS, viewer));
792e1d0745SJose E. Roman     PetscCall(ISDestroy(&coffIS));
802e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5PopGroup(viewer));
812e1d0745SJose E. Roman     PetscCall(ISCreateGeneral(comm, m, cinds, PETSC_OWN_POINTER, &cindIS));
822e1d0745SJose E. Roman     PetscCall(PetscObjectSetName((PetscObject)cindIS, "bcIndices"));
832e1d0745SJose E. Roman     PetscCall(ISView(cindIS, viewer));
842e1d0745SJose E. Roman     PetscCall(ISDestroy(&cindIS));
852e1d0745SJose E. Roman   }
862e1d0745SJose E. Roman   PetscFunctionReturn(PETSC_SUCCESS);
872e1d0745SJose E. Roman }
882e1d0745SJose E. Roman 
PetscSectionView_HDF5_Internal(PetscSection s,PetscViewer viewer)892e1d0745SJose E. Roman PetscErrorCode PetscSectionView_HDF5_Internal(PetscSection s, PetscViewer viewer)
902e1d0745SJose E. Roman {
912e1d0745SJose E. Roman   PetscInt numFields, f;
922e1d0745SJose E. Roman 
932e1d0745SJose E. Roman   PetscFunctionBegin;
942e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5PushGroup(viewer, "section"));
952e1d0745SJose E. Roman   PetscCall(PetscSectionGetNumFields(s, &numFields));
962e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "numFields", PETSC_INT, (void *)&numFields));
972e1d0745SJose E. Roman   PetscCall(PetscSectionView_HDF5_SingleField(s, viewer));
982e1d0745SJose E. Roman   for (f = 0; f < numFields; ++f) {
992e1d0745SJose E. Roman     char        fname[PETSC_MAX_PATH_LEN];
1002e1d0745SJose E. Roman     const char *fieldName;
1012e1d0745SJose E. Roman     PetscInt    fieldComponents, c;
1022e1d0745SJose E. Roman 
1032e1d0745SJose E. Roman     PetscCall(PetscSNPrintf(fname, sizeof(fname), "field%" PetscInt_FMT, f));
1042e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5PushGroup(viewer, fname));
1052e1d0745SJose E. Roman     PetscCall(PetscSectionGetFieldName(s, f, &fieldName));
1062e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "fieldName", PETSC_STRING, fieldName));
1072e1d0745SJose E. Roman     PetscCall(PetscSectionGetFieldComponents(s, f, &fieldComponents));
1082e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "fieldComponents", PETSC_INT, (void *)&fieldComponents));
1092e1d0745SJose E. Roman     for (c = 0; c < fieldComponents; ++c) {
1102e1d0745SJose E. Roman       char        cname[PETSC_MAX_PATH_LEN];
1112e1d0745SJose E. Roman       const char *componentName;
1122e1d0745SJose E. Roman 
1132e1d0745SJose E. Roman       PetscCall(PetscSNPrintf(cname, sizeof(cname), "component%" PetscInt_FMT, c));
1142e1d0745SJose E. Roman       PetscCall(PetscViewerHDF5PushGroup(viewer, cname));
1152e1d0745SJose E. Roman       PetscCall(PetscSectionGetComponentName(s, f, c, &componentName));
1162e1d0745SJose E. Roman       PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "componentName", PETSC_STRING, componentName));
1172e1d0745SJose E. Roman       PetscCall(PetscViewerHDF5PopGroup(viewer));
1182e1d0745SJose E. Roman     }
1192e1d0745SJose E. Roman     PetscCall(PetscSectionView_HDF5_SingleField(s->field[f], viewer));
1202e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5PopGroup(viewer));
1212e1d0745SJose E. Roman   }
1222e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5PopGroup(viewer));
1232e1d0745SJose E. Roman   PetscFunctionReturn(PETSC_SUCCESS);
1242e1d0745SJose E. Roman }
1252e1d0745SJose E. Roman 
PetscSectionLoad_HDF5_SingleField_SetConstraintIndices(PetscSection s,IS cindIS,IS coffIS)1262e1d0745SJose E. Roman static PetscErrorCode PetscSectionLoad_HDF5_SingleField_SetConstraintIndices(PetscSection s, IS cindIS, IS coffIS)
1272e1d0745SJose E. Roman {
1282e1d0745SJose E. Roman   MPI_Comm        comm;
1292e1d0745SJose E. Roman   PetscInt        pStart, pEnd, p, M, m, i, cdof;
1302e1d0745SJose E. Roman   const PetscInt *data;
1312e1d0745SJose E. Roman   PetscInt       *cinds;
1322e1d0745SJose E. Roman   const PetscInt *coffs;
1332e1d0745SJose E. Roman   PetscInt       *coffsets;
1342e1d0745SJose E. Roman   PetscSF         sf;
1352e1d0745SJose E. Roman   PetscLayout     layout;
1362e1d0745SJose E. Roman 
1372e1d0745SJose E. Roman   PetscFunctionBegin;
1382e1d0745SJose E. Roman   PetscCall(PetscObjectGetComm((PetscObject)s, &comm));
1392e1d0745SJose E. Roman   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
1402e1d0745SJose E. Roman   PetscCall(ISGetSize(cindIS, &M));
1412e1d0745SJose E. Roman   PetscCall(ISGetLocalSize(cindIS, &m));
1422e1d0745SJose E. Roman   PetscCall(PetscMalloc1(m, &coffsets));
1432e1d0745SJose E. Roman   PetscCall(ISGetIndices(coffIS, &coffs));
1442e1d0745SJose E. Roman   for (p = pStart, m = 0; p < pEnd; ++p) {
1452e1d0745SJose E. Roman     PetscCall(PetscSectionGetConstraintDof(s, p, &cdof));
1462e1d0745SJose E. Roman     for (i = 0; i < cdof; ++i) coffsets[m++] = coffs[p - pStart] + i;
1472e1d0745SJose E. Roman   }
1482e1d0745SJose E. Roman   PetscCall(ISRestoreIndices(coffIS, &coffs));
1492e1d0745SJose E. Roman   PetscCall(PetscSFCreate(comm, &sf));
1502e1d0745SJose E. Roman   PetscCall(PetscLayoutCreate(comm, &layout));
1512e1d0745SJose E. Roman   PetscCall(PetscLayoutSetSize(layout, M));
1522e1d0745SJose E. Roman   PetscCall(PetscLayoutSetLocalSize(layout, m));
1532e1d0745SJose E. Roman   PetscCall(PetscLayoutSetBlockSize(layout, 1));
1542e1d0745SJose E. Roman   PetscCall(PetscLayoutSetUp(layout));
1552e1d0745SJose E. Roman   PetscCall(PetscSFSetGraphLayout(sf, layout, m, NULL, PETSC_OWN_POINTER, coffsets));
1562e1d0745SJose E. Roman   PetscCall(PetscLayoutDestroy(&layout));
1572e1d0745SJose E. Roman   PetscCall(PetscFree(coffsets));
1582e1d0745SJose E. Roman   PetscCall(PetscMalloc1(m, &cinds));
1592e1d0745SJose E. Roman   PetscCall(ISGetIndices(cindIS, &data));
1602e1d0745SJose E. Roman   PetscCall(PetscSFBcastBegin(sf, MPIU_INT, data, cinds, MPI_REPLACE));
1612e1d0745SJose E. Roman   PetscCall(PetscSFBcastEnd(sf, MPIU_INT, data, cinds, MPI_REPLACE));
1622e1d0745SJose E. Roman   PetscCall(ISRestoreIndices(cindIS, &data));
1632e1d0745SJose E. Roman   PetscCall(PetscSFDestroy(&sf));
1642e1d0745SJose E. Roman   PetscCall(PetscSectionSetUpBC(s));
1652e1d0745SJose E. Roman   for (p = pStart, m = 0; p < pEnd; ++p) {
1662e1d0745SJose E. Roman     PetscCall(PetscSectionGetConstraintDof(s, p, &cdof));
1672e1d0745SJose E. Roman     PetscCall(PetscSectionSetConstraintIndices(s, p, PetscSafePointerPlusOffset(cinds, m)));
1682e1d0745SJose E. Roman     m += cdof;
1692e1d0745SJose E. Roman   }
1702e1d0745SJose E. Roman   PetscCall(PetscFree(cinds));
1712e1d0745SJose E. Roman   PetscFunctionReturn(PETSC_SUCCESS);
1722e1d0745SJose E. Roman }
1732e1d0745SJose E. Roman 
PetscSectionLoad_HDF5_SingleField(PetscSection s,PetscViewer viewer)1742e1d0745SJose E. Roman static PetscErrorCode PetscSectionLoad_HDF5_SingleField(PetscSection s, PetscViewer viewer)
1752e1d0745SJose E. Roman {
1762e1d0745SJose E. Roman   MPI_Comm comm;
1772e1d0745SJose E. Roman   PetscInt pStart, pEnd, p, N, n, M, m;
1782e1d0745SJose E. Roman #if defined(PETSC_USE_DEBUG)
1792e1d0745SJose E. Roman   PetscInt N1, M1;
1802e1d0745SJose E. Roman #endif
1812e1d0745SJose E. Roman   PetscBool       hasConstraints, includesConstraints;
1822e1d0745SJose E. Roman   IS              dofIS, offIS, cdofIS, coffIS, cindIS;
1832e1d0745SJose E. Roman   const PetscInt *dofs, *offs, *cdofs;
1842e1d0745SJose E. Roman   PetscLayout     map;
1852e1d0745SJose E. Roman 
1862e1d0745SJose E. Roman   PetscFunctionBegin;
1872e1d0745SJose E. Roman   PetscCall(PetscObjectGetComm((PetscObject)s, &comm));
1882e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "includesConstraints", PETSC_BOOL, NULL, (void *)&includesConstraints));
1892e1d0745SJose E. Roman   PetscCall(PetscSectionSetIncludesConstraints(s, includesConstraints));
1902e1d0745SJose E. Roman   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
1912e1d0745SJose E. Roman   n = pEnd - pStart;
1922e1d0745SJose E. Roman #if defined(PETSC_USE_DEBUG)
1932e1d0745SJose E. Roman   PetscCallMPI(MPIU_Allreduce(&n, &N1, 1, MPIU_INT, MPI_SUM, comm));
1942e1d0745SJose E. Roman #endif
1952e1d0745SJose E. Roman   PetscCall(ISCreate(comm, &dofIS));
1962e1d0745SJose E. Roman   PetscCall(PetscObjectSetName((PetscObject)dofIS, "atlasDof"));
1972e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5ReadSizes(viewer, "atlasDof", NULL, &N));
1982e1d0745SJose E. Roman #if defined(PETSC_USE_DEBUG)
1992e1d0745SJose E. Roman   PetscCheck(N1 == N, comm, PETSC_ERR_ARG_SIZ, "Unable to load s->atlasDof: sum of local sizes (%" PetscInt_FMT ") != global size (%" PetscInt_FMT "): local size on this process is %" PetscInt_FMT, N1, N, n);
2002e1d0745SJose E. Roman #endif
2012e1d0745SJose E. Roman   PetscCall(ISGetLayout(dofIS, &map));
2022e1d0745SJose E. Roman   PetscCall(PetscLayoutSetSize(map, N));
2032e1d0745SJose E. Roman   PetscCall(PetscLayoutSetLocalSize(map, n));
2042e1d0745SJose E. Roman   PetscCall(ISLoad(dofIS, viewer));
2052e1d0745SJose E. Roman   PetscCall(ISCreate(comm, &offIS));
2062e1d0745SJose E. Roman   PetscCall(PetscObjectSetName((PetscObject)offIS, "atlasOff"));
2072e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5ReadSizes(viewer, "atlasOff", NULL, &N));
2082e1d0745SJose E. Roman #if defined(PETSC_USE_DEBUG)
2092e1d0745SJose E. Roman   PetscCheck(N1 == N, comm, PETSC_ERR_ARG_SIZ, "Unable to load s->atlasOff: sum of local sizes (%" PetscInt_FMT ") != global size (%" PetscInt_FMT "): local size on this process is %" PetscInt_FMT, N1, N, n);
2102e1d0745SJose E. Roman #endif
2112e1d0745SJose E. Roman   PetscCall(ISGetLayout(offIS, &map));
2122e1d0745SJose E. Roman   PetscCall(PetscLayoutSetSize(map, N));
2132e1d0745SJose E. Roman   PetscCall(PetscLayoutSetLocalSize(map, n));
2142e1d0745SJose E. Roman   PetscCall(ISLoad(offIS, viewer));
2152e1d0745SJose E. Roman   PetscCall(ISGetIndices(dofIS, &dofs));
2162e1d0745SJose E. Roman   PetscCall(ISGetIndices(offIS, &offs));
2172e1d0745SJose E. Roman   for (p = pStart, n = 0; p < pEnd; ++p, ++n) {
2182e1d0745SJose E. Roman     PetscCall(PetscSectionSetDof(s, p, dofs[n]));
2192e1d0745SJose E. Roman     PetscCall(PetscSectionSetOffset(s, p, offs[n]));
2202e1d0745SJose E. Roman   }
2212e1d0745SJose E. Roman   PetscCall(ISRestoreIndices(dofIS, &dofs));
2222e1d0745SJose E. Roman   PetscCall(ISRestoreIndices(offIS, &offs));
2232e1d0745SJose E. Roman   PetscCall(ISDestroy(&dofIS));
2242e1d0745SJose E. Roman   PetscCall(ISDestroy(&offIS));
2252e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "hasConstraints", PETSC_BOOL, NULL, (void *)&hasConstraints));
2262e1d0745SJose E. Roman   if (hasConstraints) {
2272e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5PushGroup(viewer, "bc"));
2282e1d0745SJose E. Roman     PetscCall(ISCreate(comm, &cdofIS));
2292e1d0745SJose E. Roman     PetscCall(PetscObjectSetName((PetscObject)cdofIS, "atlasDof"));
2302e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5ReadSizes(viewer, "atlasDof", NULL, &N));
2312e1d0745SJose E. Roman #if defined(PETSC_USE_DEBUG)
2322e1d0745SJose E. Roman     PetscCheck(N1 == N, comm, PETSC_ERR_ARG_SIZ, "Unable to load s->bc->atlasDof: sum of local sizes (%" PetscInt_FMT ") != global size (%" PetscInt_FMT "): local size on this process is %" PetscInt_FMT, N1, N, n);
2332e1d0745SJose E. Roman #endif
2342e1d0745SJose E. Roman     PetscCall(ISGetLayout(cdofIS, &map));
2352e1d0745SJose E. Roman     PetscCall(PetscLayoutSetSize(map, N));
2362e1d0745SJose E. Roman     PetscCall(PetscLayoutSetLocalSize(map, n));
2372e1d0745SJose E. Roman     PetscCall(ISLoad(cdofIS, viewer));
2382e1d0745SJose E. Roman     PetscCall(ISGetIndices(cdofIS, &cdofs));
2392e1d0745SJose E. Roman     for (p = pStart, n = 0; p < pEnd; ++p, ++n) PetscCall(PetscSectionSetConstraintDof(s, p, cdofs[n]));
2402e1d0745SJose E. Roman     PetscCall(ISRestoreIndices(cdofIS, &cdofs));
2412e1d0745SJose E. Roman     PetscCall(ISDestroy(&cdofIS));
2422e1d0745SJose E. Roman     PetscCall(ISCreate(comm, &coffIS));
2432e1d0745SJose E. Roman     PetscCall(PetscObjectSetName((PetscObject)coffIS, "atlasOff"));
2442e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5ReadSizes(viewer, "atlasOff", NULL, &N));
2452e1d0745SJose E. Roman #if defined(PETSC_USE_DEBUG)
2462e1d0745SJose E. Roman     PetscCheck(N1 == N, comm, PETSC_ERR_ARG_SIZ, "Unable to load s->bc->atlasOff: sum of local sizes (%" PetscInt_FMT ") != global size (%" PetscInt_FMT "): local size on this process is %" PetscInt_FMT, N1, N, n);
2472e1d0745SJose E. Roman #endif
2482e1d0745SJose E. Roman     PetscCall(ISGetLayout(coffIS, &map));
2492e1d0745SJose E. Roman     PetscCall(PetscLayoutSetSize(map, N));
2502e1d0745SJose E. Roman     PetscCall(PetscLayoutSetLocalSize(map, n));
2512e1d0745SJose E. Roman     PetscCall(ISLoad(coffIS, viewer));
2522e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5PopGroup(viewer));
2532e1d0745SJose E. Roman     PetscCall(ISCreate(comm, &cindIS));
2542e1d0745SJose E. Roman     PetscCall(PetscObjectSetName((PetscObject)cindIS, "bcIndices"));
2552e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5ReadSizes(viewer, "bcIndices", NULL, &M));
2562e1d0745SJose E. Roman     if (!s->bc) m = 0;
2572e1d0745SJose E. Roman     else PetscCall(PetscSectionGetStorageSize(s->bc, &m));
2582e1d0745SJose E. Roman #if defined(PETSC_USE_DEBUG)
2592e1d0745SJose E. Roman     PetscCallMPI(MPIU_Allreduce(&m, &M1, 1, MPIU_INT, MPI_SUM, comm));
2602e1d0745SJose E. Roman     PetscCheck(M1 == M, comm, PETSC_ERR_ARG_SIZ, "Unable to load s->bcIndices: sum of local sizes (%" PetscInt_FMT ") != global size (%" PetscInt_FMT "): local size on this process is %" PetscInt_FMT, M1, M, m);
2612e1d0745SJose E. Roman #endif
2622e1d0745SJose E. Roman     PetscCall(ISGetLayout(cindIS, &map));
2632e1d0745SJose E. Roman     PetscCall(PetscLayoutSetSize(map, M));
2642e1d0745SJose E. Roman     PetscCall(PetscLayoutSetLocalSize(map, m));
2652e1d0745SJose E. Roman     PetscCall(ISLoad(cindIS, viewer));
2662e1d0745SJose E. Roman     PetscCall(PetscSectionLoad_HDF5_SingleField_SetConstraintIndices(s, cindIS, coffIS));
2672e1d0745SJose E. Roman     PetscCall(ISDestroy(&coffIS));
2682e1d0745SJose E. Roman     PetscCall(ISDestroy(&cindIS));
2692e1d0745SJose E. Roman   }
2702e1d0745SJose E. Roman   PetscFunctionReturn(PETSC_SUCCESS);
2712e1d0745SJose E. Roman }
2722e1d0745SJose E. Roman 
PetscSectionLoad_HDF5_Internal(PetscSection s,PetscViewer viewer)2732e1d0745SJose E. Roman PetscErrorCode PetscSectionLoad_HDF5_Internal(PetscSection s, PetscViewer viewer)
2742e1d0745SJose E. Roman {
2752e1d0745SJose E. Roman   MPI_Comm comm;
2762e1d0745SJose E. Roman   PetscInt N, n, numFields, f;
2772e1d0745SJose E. Roman 
2782e1d0745SJose E. Roman   PetscFunctionBegin;
2792e1d0745SJose E. Roman   PetscCall(PetscObjectGetComm((PetscObject)s, &comm));
2802e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5PushGroup(viewer, "section"));
2812e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "numFields", PETSC_INT, NULL, (void *)&numFields));
2822e1d0745SJose E. Roman   if (s->pStart < 0 && s->pEnd < 0) n = PETSC_DECIDE;
2832e1d0745SJose E. Roman   else {
2842e1d0745SJose E. Roman     PetscCheck(s->pStart == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "s->pStart must be 0 (got %" PetscInt_FMT ")", s->pStart);
2852e1d0745SJose E. Roman     PetscCheck(s->pEnd >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "s->pEnd must be >= 0, (got %" PetscInt_FMT ")", s->pEnd);
2862e1d0745SJose E. Roman     n = s->pEnd;
2872e1d0745SJose E. Roman   }
2882e1d0745SJose E. Roman   if (numFields > 0) PetscCall(PetscSectionSetNumFields(s, numFields));
2892e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5ReadSizes(viewer, "atlasDof", NULL, &N));
2902e1d0745SJose E. Roman   if (n == PETSC_DECIDE) PetscCall(PetscSplitOwnership(comm, &n, &N));
2912e1d0745SJose E. Roman   PetscCall(PetscSectionSetChart(s, 0, n));
2922e1d0745SJose E. Roman   PetscCall(PetscSectionLoad_HDF5_SingleField(s, viewer));
2932e1d0745SJose E. Roman   for (f = 0; f < numFields; ++f) {
2942e1d0745SJose E. Roman     char     fname[PETSC_MAX_PATH_LEN];
2952e1d0745SJose E. Roman     char    *fieldName;
2962e1d0745SJose E. Roman     PetscInt fieldComponents, c;
2972e1d0745SJose E. Roman 
2982e1d0745SJose E. Roman     PetscCall(PetscSNPrintf(fname, sizeof(fname), "field%" PetscInt_FMT, f));
2992e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5PushGroup(viewer, fname));
3002e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "fieldName", PETSC_STRING, NULL, &fieldName));
3012e1d0745SJose E. Roman     PetscCall(PetscSectionSetFieldName(s, f, fieldName));
3022e1d0745SJose E. Roman     PetscCall(PetscFree(fieldName));
3032e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "fieldComponents", PETSC_INT, NULL, (void *)&fieldComponents));
3042e1d0745SJose E. Roman     PetscCall(PetscSectionSetFieldComponents(s, f, fieldComponents));
3052e1d0745SJose E. Roman     for (c = 0; c < fieldComponents; ++c) {
3062e1d0745SJose E. Roman       char  cname[PETSC_MAX_PATH_LEN];
3072e1d0745SJose E. Roman       char *componentName;
3082e1d0745SJose E. Roman 
3092e1d0745SJose E. Roman       PetscCall(PetscSNPrintf(cname, sizeof(cname), "component%" PetscInt_FMT, c));
3102e1d0745SJose E. Roman       PetscCall(PetscViewerHDF5PushGroup(viewer, cname));
3112e1d0745SJose E. Roman       PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "componentName", PETSC_STRING, NULL, &componentName));
3122e1d0745SJose E. Roman       PetscCall(PetscSectionSetComponentName(s, f, c, componentName));
3132e1d0745SJose E. Roman       PetscCall(PetscFree(componentName));
3142e1d0745SJose E. Roman       PetscCall(PetscViewerHDF5PopGroup(viewer));
3152e1d0745SJose E. Roman     }
3162e1d0745SJose E. Roman     PetscCall(PetscSectionLoad_HDF5_SingleField(s->field[f], viewer));
3172e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5PopGroup(viewer));
3182e1d0745SJose E. Roman   }
3192e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5PopGroup(viewer));
3202e1d0745SJose E. Roman   PetscFunctionReturn(PETSC_SUCCESS);
3212e1d0745SJose E. Roman }
322