xref: /petsc/src/vec/is/section/interface/hdf5/sectionhdf5.c (revision 2e1d0745e031e9eda42a183866b5a94ccb2b4e44)
1*2e1d0745SJose E. Roman #include <petsc/private/sectionimpl.h> /*I "petscsection.h" I*/
2*2e1d0745SJose E. Roman #include <petscsf.h>
3*2e1d0745SJose E. Roman #include <petscis.h>
4*2e1d0745SJose E. Roman #include <petscviewerhdf5.h>
5*2e1d0745SJose E. Roman #include <petsclayouthdf5.h>
6*2e1d0745SJose E. Roman 
7*2e1d0745SJose E. Roman static PetscErrorCode PetscSectionView_HDF5_SingleField(PetscSection s, PetscViewer viewer)
8*2e1d0745SJose E. Roman {
9*2e1d0745SJose E. Roman   MPI_Comm  comm;
10*2e1d0745SJose E. Roman   PetscInt  pStart, pEnd, p, n;
11*2e1d0745SJose E. Roman   PetscBool hasConstraints, includesConstraints;
12*2e1d0745SJose E. Roman   IS        dofIS, offIS, cdofIS, coffIS, cindIS;
13*2e1d0745SJose E. Roman   PetscInt *dofs, *offs, *cdofs, *coffs, *cinds, dof, cdof, m, moff, i;
14*2e1d0745SJose E. Roman 
15*2e1d0745SJose E. Roman   PetscFunctionBegin;
16*2e1d0745SJose E. Roman   PetscCall(PetscObjectGetComm((PetscObject)s, &comm));
17*2e1d0745SJose E. Roman   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
18*2e1d0745SJose E. Roman   hasConstraints = (s->bc) ? PETSC_TRUE : PETSC_FALSE;
19*2e1d0745SJose E. Roman   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &hasConstraints, 1, MPIU_BOOL, MPI_LOR, comm));
20*2e1d0745SJose E. Roman   for (p = pStart, n = 0, m = 0; p < pEnd; ++p) {
21*2e1d0745SJose E. Roman     PetscCall(PetscSectionGetDof(s, p, &dof));
22*2e1d0745SJose E. Roman     if (dof >= 0) {
23*2e1d0745SJose E. Roman       if (hasConstraints) {
24*2e1d0745SJose E. Roman         PetscCall(PetscSectionGetConstraintDof(s, p, &cdof));
25*2e1d0745SJose E. Roman         m += cdof;
26*2e1d0745SJose E. Roman       }
27*2e1d0745SJose E. Roman       n++;
28*2e1d0745SJose E. Roman     }
29*2e1d0745SJose E. Roman   }
30*2e1d0745SJose E. Roman   PetscCall(PetscMalloc1(n, &dofs));
31*2e1d0745SJose E. Roman   PetscCall(PetscMalloc1(n, &offs));
32*2e1d0745SJose E. Roman   if (hasConstraints) {
33*2e1d0745SJose E. Roman     PetscCall(PetscMalloc1(n, &cdofs));
34*2e1d0745SJose E. Roman     PetscCall(PetscMalloc1(n, &coffs));
35*2e1d0745SJose E. Roman     PetscCall(PetscMalloc1(m, &cinds));
36*2e1d0745SJose E. Roman   }
37*2e1d0745SJose E. Roman   for (p = pStart, n = 0, m = 0; p < pEnd; ++p) {
38*2e1d0745SJose E. Roman     PetscCall(PetscSectionGetDof(s, p, &dof));
39*2e1d0745SJose E. Roman     if (dof >= 0) {
40*2e1d0745SJose E. Roman       dofs[n] = dof;
41*2e1d0745SJose E. Roman       PetscCall(PetscSectionGetOffset(s, p, &offs[n]));
42*2e1d0745SJose E. Roman       if (hasConstraints) {
43*2e1d0745SJose E. Roman         const PetscInt *cpinds;
44*2e1d0745SJose E. Roman 
45*2e1d0745SJose E. Roman         PetscCall(PetscSectionGetConstraintDof(s, p, &cdof));
46*2e1d0745SJose E. Roman         PetscCall(PetscSectionGetConstraintIndices(s, p, &cpinds));
47*2e1d0745SJose E. Roman         cdofs[n] = cdof;
48*2e1d0745SJose E. Roman         coffs[n] = m;
49*2e1d0745SJose E. Roman         for (i = 0; i < cdof; ++i) cinds[m++] = cpinds[i];
50*2e1d0745SJose E. Roman       }
51*2e1d0745SJose E. Roman       n++;
52*2e1d0745SJose E. Roman     }
53*2e1d0745SJose E. Roman   }
54*2e1d0745SJose E. Roman   if (hasConstraints) {
55*2e1d0745SJose E. Roman     PetscCallMPI(MPI_Scan(&m, &moff, 1, MPIU_INT, MPI_SUM, comm));
56*2e1d0745SJose E. Roman     moff -= m;
57*2e1d0745SJose E. Roman     for (p = 0; p < n; ++p) coffs[p] += moff;
58*2e1d0745SJose E. Roman   }
59*2e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "hasConstraints", PETSC_BOOL, (void *)&hasConstraints));
60*2e1d0745SJose E. Roman   PetscCall(PetscSectionGetIncludesConstraints(s, &includesConstraints));
61*2e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "includesConstraints", PETSC_BOOL, (void *)&includesConstraints));
62*2e1d0745SJose E. Roman   PetscCall(ISCreateGeneral(comm, n, dofs, PETSC_OWN_POINTER, &dofIS));
63*2e1d0745SJose E. Roman   PetscCall(PetscObjectSetName((PetscObject)dofIS, "atlasDof"));
64*2e1d0745SJose E. Roman   PetscCall(ISView(dofIS, viewer));
65*2e1d0745SJose E. Roman   PetscCall(ISDestroy(&dofIS));
66*2e1d0745SJose E. Roman   PetscCall(ISCreateGeneral(comm, n, offs, PETSC_OWN_POINTER, &offIS));
67*2e1d0745SJose E. Roman   PetscCall(PetscObjectSetName((PetscObject)offIS, "atlasOff"));
68*2e1d0745SJose E. Roman   PetscCall(ISView(offIS, viewer));
69*2e1d0745SJose E. Roman   PetscCall(ISDestroy(&offIS));
70*2e1d0745SJose E. Roman   if (hasConstraints) {
71*2e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5PushGroup(viewer, "bc"));
72*2e1d0745SJose E. Roman     PetscCall(ISCreateGeneral(comm, n, cdofs, PETSC_OWN_POINTER, &cdofIS));
73*2e1d0745SJose E. Roman     PetscCall(PetscObjectSetName((PetscObject)cdofIS, "atlasDof"));
74*2e1d0745SJose E. Roman     PetscCall(ISView(cdofIS, viewer));
75*2e1d0745SJose E. Roman     PetscCall(ISDestroy(&cdofIS));
76*2e1d0745SJose E. Roman     PetscCall(ISCreateGeneral(comm, n, coffs, PETSC_OWN_POINTER, &coffIS));
77*2e1d0745SJose E. Roman     PetscCall(PetscObjectSetName((PetscObject)coffIS, "atlasOff"));
78*2e1d0745SJose E. Roman     PetscCall(ISView(coffIS, viewer));
79*2e1d0745SJose E. Roman     PetscCall(ISDestroy(&coffIS));
80*2e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5PopGroup(viewer));
81*2e1d0745SJose E. Roman     PetscCall(ISCreateGeneral(comm, m, cinds, PETSC_OWN_POINTER, &cindIS));
82*2e1d0745SJose E. Roman     PetscCall(PetscObjectSetName((PetscObject)cindIS, "bcIndices"));
83*2e1d0745SJose E. Roman     PetscCall(ISView(cindIS, viewer));
84*2e1d0745SJose E. Roman     PetscCall(ISDestroy(&cindIS));
85*2e1d0745SJose E. Roman   }
86*2e1d0745SJose E. Roman   PetscFunctionReturn(PETSC_SUCCESS);
87*2e1d0745SJose E. Roman }
88*2e1d0745SJose E. Roman 
89*2e1d0745SJose E. Roman PetscErrorCode PetscSectionView_HDF5_Internal(PetscSection s, PetscViewer viewer)
90*2e1d0745SJose E. Roman {
91*2e1d0745SJose E. Roman   PetscInt numFields, f;
92*2e1d0745SJose E. Roman 
93*2e1d0745SJose E. Roman   PetscFunctionBegin;
94*2e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5PushGroup(viewer, "section"));
95*2e1d0745SJose E. Roman   PetscCall(PetscSectionGetNumFields(s, &numFields));
96*2e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "numFields", PETSC_INT, (void *)&numFields));
97*2e1d0745SJose E. Roman   PetscCall(PetscSectionView_HDF5_SingleField(s, viewer));
98*2e1d0745SJose E. Roman   for (f = 0; f < numFields; ++f) {
99*2e1d0745SJose E. Roman     char        fname[PETSC_MAX_PATH_LEN];
100*2e1d0745SJose E. Roman     const char *fieldName;
101*2e1d0745SJose E. Roman     PetscInt    fieldComponents, c;
102*2e1d0745SJose E. Roman 
103*2e1d0745SJose E. Roman     PetscCall(PetscSNPrintf(fname, sizeof(fname), "field%" PetscInt_FMT, f));
104*2e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5PushGroup(viewer, fname));
105*2e1d0745SJose E. Roman     PetscCall(PetscSectionGetFieldName(s, f, &fieldName));
106*2e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "fieldName", PETSC_STRING, fieldName));
107*2e1d0745SJose E. Roman     PetscCall(PetscSectionGetFieldComponents(s, f, &fieldComponents));
108*2e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "fieldComponents", PETSC_INT, (void *)&fieldComponents));
109*2e1d0745SJose E. Roman     for (c = 0; c < fieldComponents; ++c) {
110*2e1d0745SJose E. Roman       char        cname[PETSC_MAX_PATH_LEN];
111*2e1d0745SJose E. Roman       const char *componentName;
112*2e1d0745SJose E. Roman 
113*2e1d0745SJose E. Roman       PetscCall(PetscSNPrintf(cname, sizeof(cname), "component%" PetscInt_FMT, c));
114*2e1d0745SJose E. Roman       PetscCall(PetscViewerHDF5PushGroup(viewer, cname));
115*2e1d0745SJose E. Roman       PetscCall(PetscSectionGetComponentName(s, f, c, &componentName));
116*2e1d0745SJose E. Roman       PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "componentName", PETSC_STRING, componentName));
117*2e1d0745SJose E. Roman       PetscCall(PetscViewerHDF5PopGroup(viewer));
118*2e1d0745SJose E. Roman     }
119*2e1d0745SJose E. Roman     PetscCall(PetscSectionView_HDF5_SingleField(s->field[f], viewer));
120*2e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5PopGroup(viewer));
121*2e1d0745SJose E. Roman   }
122*2e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5PopGroup(viewer));
123*2e1d0745SJose E. Roman   PetscFunctionReturn(PETSC_SUCCESS);
124*2e1d0745SJose E. Roman }
125*2e1d0745SJose E. Roman 
126*2e1d0745SJose E. Roman static PetscErrorCode PetscSectionLoad_HDF5_SingleField_SetConstraintIndices(PetscSection s, IS cindIS, IS coffIS)
127*2e1d0745SJose E. Roman {
128*2e1d0745SJose E. Roman   MPI_Comm        comm;
129*2e1d0745SJose E. Roman   PetscInt        pStart, pEnd, p, M, m, i, cdof;
130*2e1d0745SJose E. Roman   const PetscInt *data;
131*2e1d0745SJose E. Roman   PetscInt       *cinds;
132*2e1d0745SJose E. Roman   const PetscInt *coffs;
133*2e1d0745SJose E. Roman   PetscInt       *coffsets;
134*2e1d0745SJose E. Roman   PetscSF         sf;
135*2e1d0745SJose E. Roman   PetscLayout     layout;
136*2e1d0745SJose E. Roman 
137*2e1d0745SJose E. Roman   PetscFunctionBegin;
138*2e1d0745SJose E. Roman   PetscCall(PetscObjectGetComm((PetscObject)s, &comm));
139*2e1d0745SJose E. Roman   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
140*2e1d0745SJose E. Roman   PetscCall(ISGetSize(cindIS, &M));
141*2e1d0745SJose E. Roman   PetscCall(ISGetLocalSize(cindIS, &m));
142*2e1d0745SJose E. Roman   PetscCall(PetscMalloc1(m, &coffsets));
143*2e1d0745SJose E. Roman   PetscCall(ISGetIndices(coffIS, &coffs));
144*2e1d0745SJose E. Roman   for (p = pStart, m = 0; p < pEnd; ++p) {
145*2e1d0745SJose E. Roman     PetscCall(PetscSectionGetConstraintDof(s, p, &cdof));
146*2e1d0745SJose E. Roman     for (i = 0; i < cdof; ++i) coffsets[m++] = coffs[p - pStart] + i;
147*2e1d0745SJose E. Roman   }
148*2e1d0745SJose E. Roman   PetscCall(ISRestoreIndices(coffIS, &coffs));
149*2e1d0745SJose E. Roman   PetscCall(PetscSFCreate(comm, &sf));
150*2e1d0745SJose E. Roman   PetscCall(PetscLayoutCreate(comm, &layout));
151*2e1d0745SJose E. Roman   PetscCall(PetscLayoutSetSize(layout, M));
152*2e1d0745SJose E. Roman   PetscCall(PetscLayoutSetLocalSize(layout, m));
153*2e1d0745SJose E. Roman   PetscCall(PetscLayoutSetBlockSize(layout, 1));
154*2e1d0745SJose E. Roman   PetscCall(PetscLayoutSetUp(layout));
155*2e1d0745SJose E. Roman   PetscCall(PetscSFSetGraphLayout(sf, layout, m, NULL, PETSC_OWN_POINTER, coffsets));
156*2e1d0745SJose E. Roman   PetscCall(PetscLayoutDestroy(&layout));
157*2e1d0745SJose E. Roman   PetscCall(PetscFree(coffsets));
158*2e1d0745SJose E. Roman   PetscCall(PetscMalloc1(m, &cinds));
159*2e1d0745SJose E. Roman   PetscCall(ISGetIndices(cindIS, &data));
160*2e1d0745SJose E. Roman   PetscCall(PetscSFBcastBegin(sf, MPIU_INT, data, cinds, MPI_REPLACE));
161*2e1d0745SJose E. Roman   PetscCall(PetscSFBcastEnd(sf, MPIU_INT, data, cinds, MPI_REPLACE));
162*2e1d0745SJose E. Roman   PetscCall(ISRestoreIndices(cindIS, &data));
163*2e1d0745SJose E. Roman   PetscCall(PetscSFDestroy(&sf));
164*2e1d0745SJose E. Roman   PetscCall(PetscSectionSetUpBC(s));
165*2e1d0745SJose E. Roman   for (p = pStart, m = 0; p < pEnd; ++p) {
166*2e1d0745SJose E. Roman     PetscCall(PetscSectionGetConstraintDof(s, p, &cdof));
167*2e1d0745SJose E. Roman     PetscCall(PetscSectionSetConstraintIndices(s, p, PetscSafePointerPlusOffset(cinds, m)));
168*2e1d0745SJose E. Roman     m += cdof;
169*2e1d0745SJose E. Roman   }
170*2e1d0745SJose E. Roman   PetscCall(PetscFree(cinds));
171*2e1d0745SJose E. Roman   PetscFunctionReturn(PETSC_SUCCESS);
172*2e1d0745SJose E. Roman }
173*2e1d0745SJose E. Roman 
174*2e1d0745SJose E. Roman static PetscErrorCode PetscSectionLoad_HDF5_SingleField(PetscSection s, PetscViewer viewer)
175*2e1d0745SJose E. Roman {
176*2e1d0745SJose E. Roman   MPI_Comm comm;
177*2e1d0745SJose E. Roman   PetscInt pStart, pEnd, p, N, n, M, m;
178*2e1d0745SJose E. Roman #if defined(PETSC_USE_DEBUG)
179*2e1d0745SJose E. Roman   PetscInt N1, M1;
180*2e1d0745SJose E. Roman #endif
181*2e1d0745SJose E. Roman   PetscBool       hasConstraints, includesConstraints;
182*2e1d0745SJose E. Roman   IS              dofIS, offIS, cdofIS, coffIS, cindIS;
183*2e1d0745SJose E. Roman   const PetscInt *dofs, *offs, *cdofs;
184*2e1d0745SJose E. Roman   PetscLayout     map;
185*2e1d0745SJose E. Roman 
186*2e1d0745SJose E. Roman   PetscFunctionBegin;
187*2e1d0745SJose E. Roman   PetscCall(PetscObjectGetComm((PetscObject)s, &comm));
188*2e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "includesConstraints", PETSC_BOOL, NULL, (void *)&includesConstraints));
189*2e1d0745SJose E. Roman   PetscCall(PetscSectionSetIncludesConstraints(s, includesConstraints));
190*2e1d0745SJose E. Roman   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
191*2e1d0745SJose E. Roman   n = pEnd - pStart;
192*2e1d0745SJose E. Roman #if defined(PETSC_USE_DEBUG)
193*2e1d0745SJose E. Roman   PetscCallMPI(MPIU_Allreduce(&n, &N1, 1, MPIU_INT, MPI_SUM, comm));
194*2e1d0745SJose E. Roman #endif
195*2e1d0745SJose E. Roman   PetscCall(ISCreate(comm, &dofIS));
196*2e1d0745SJose E. Roman   PetscCall(PetscObjectSetName((PetscObject)dofIS, "atlasDof"));
197*2e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5ReadSizes(viewer, "atlasDof", NULL, &N));
198*2e1d0745SJose E. Roman #if defined(PETSC_USE_DEBUG)
199*2e1d0745SJose 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);
200*2e1d0745SJose E. Roman #endif
201*2e1d0745SJose E. Roman   PetscCall(ISGetLayout(dofIS, &map));
202*2e1d0745SJose E. Roman   PetscCall(PetscLayoutSetSize(map, N));
203*2e1d0745SJose E. Roman   PetscCall(PetscLayoutSetLocalSize(map, n));
204*2e1d0745SJose E. Roman   PetscCall(ISLoad(dofIS, viewer));
205*2e1d0745SJose E. Roman   PetscCall(ISCreate(comm, &offIS));
206*2e1d0745SJose E. Roman   PetscCall(PetscObjectSetName((PetscObject)offIS, "atlasOff"));
207*2e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5ReadSizes(viewer, "atlasOff", NULL, &N));
208*2e1d0745SJose E. Roman #if defined(PETSC_USE_DEBUG)
209*2e1d0745SJose 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);
210*2e1d0745SJose E. Roman #endif
211*2e1d0745SJose E. Roman   PetscCall(ISGetLayout(offIS, &map));
212*2e1d0745SJose E. Roman   PetscCall(PetscLayoutSetSize(map, N));
213*2e1d0745SJose E. Roman   PetscCall(PetscLayoutSetLocalSize(map, n));
214*2e1d0745SJose E. Roman   PetscCall(ISLoad(offIS, viewer));
215*2e1d0745SJose E. Roman   PetscCall(ISGetIndices(dofIS, &dofs));
216*2e1d0745SJose E. Roman   PetscCall(ISGetIndices(offIS, &offs));
217*2e1d0745SJose E. Roman   for (p = pStart, n = 0; p < pEnd; ++p, ++n) {
218*2e1d0745SJose E. Roman     PetscCall(PetscSectionSetDof(s, p, dofs[n]));
219*2e1d0745SJose E. Roman     PetscCall(PetscSectionSetOffset(s, p, offs[n]));
220*2e1d0745SJose E. Roman   }
221*2e1d0745SJose E. Roman   PetscCall(ISRestoreIndices(dofIS, &dofs));
222*2e1d0745SJose E. Roman   PetscCall(ISRestoreIndices(offIS, &offs));
223*2e1d0745SJose E. Roman   PetscCall(ISDestroy(&dofIS));
224*2e1d0745SJose E. Roman   PetscCall(ISDestroy(&offIS));
225*2e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "hasConstraints", PETSC_BOOL, NULL, (void *)&hasConstraints));
226*2e1d0745SJose E. Roman   if (hasConstraints) {
227*2e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5PushGroup(viewer, "bc"));
228*2e1d0745SJose E. Roman     PetscCall(ISCreate(comm, &cdofIS));
229*2e1d0745SJose E. Roman     PetscCall(PetscObjectSetName((PetscObject)cdofIS, "atlasDof"));
230*2e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5ReadSizes(viewer, "atlasDof", NULL, &N));
231*2e1d0745SJose E. Roman #if defined(PETSC_USE_DEBUG)
232*2e1d0745SJose 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);
233*2e1d0745SJose E. Roman #endif
234*2e1d0745SJose E. Roman     PetscCall(ISGetLayout(cdofIS, &map));
235*2e1d0745SJose E. Roman     PetscCall(PetscLayoutSetSize(map, N));
236*2e1d0745SJose E. Roman     PetscCall(PetscLayoutSetLocalSize(map, n));
237*2e1d0745SJose E. Roman     PetscCall(ISLoad(cdofIS, viewer));
238*2e1d0745SJose E. Roman     PetscCall(ISGetIndices(cdofIS, &cdofs));
239*2e1d0745SJose E. Roman     for (p = pStart, n = 0; p < pEnd; ++p, ++n) PetscCall(PetscSectionSetConstraintDof(s, p, cdofs[n]));
240*2e1d0745SJose E. Roman     PetscCall(ISRestoreIndices(cdofIS, &cdofs));
241*2e1d0745SJose E. Roman     PetscCall(ISDestroy(&cdofIS));
242*2e1d0745SJose E. Roman     PetscCall(ISCreate(comm, &coffIS));
243*2e1d0745SJose E. Roman     PetscCall(PetscObjectSetName((PetscObject)coffIS, "atlasOff"));
244*2e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5ReadSizes(viewer, "atlasOff", NULL, &N));
245*2e1d0745SJose E. Roman #if defined(PETSC_USE_DEBUG)
246*2e1d0745SJose 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);
247*2e1d0745SJose E. Roman #endif
248*2e1d0745SJose E. Roman     PetscCall(ISGetLayout(coffIS, &map));
249*2e1d0745SJose E. Roman     PetscCall(PetscLayoutSetSize(map, N));
250*2e1d0745SJose E. Roman     PetscCall(PetscLayoutSetLocalSize(map, n));
251*2e1d0745SJose E. Roman     PetscCall(ISLoad(coffIS, viewer));
252*2e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5PopGroup(viewer));
253*2e1d0745SJose E. Roman     PetscCall(ISCreate(comm, &cindIS));
254*2e1d0745SJose E. Roman     PetscCall(PetscObjectSetName((PetscObject)cindIS, "bcIndices"));
255*2e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5ReadSizes(viewer, "bcIndices", NULL, &M));
256*2e1d0745SJose E. Roman     if (!s->bc) m = 0;
257*2e1d0745SJose E. Roman     else PetscCall(PetscSectionGetStorageSize(s->bc, &m));
258*2e1d0745SJose E. Roman #if defined(PETSC_USE_DEBUG)
259*2e1d0745SJose E. Roman     PetscCallMPI(MPIU_Allreduce(&m, &M1, 1, MPIU_INT, MPI_SUM, comm));
260*2e1d0745SJose 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);
261*2e1d0745SJose E. Roman #endif
262*2e1d0745SJose E. Roman     PetscCall(ISGetLayout(cindIS, &map));
263*2e1d0745SJose E. Roman     PetscCall(PetscLayoutSetSize(map, M));
264*2e1d0745SJose E. Roman     PetscCall(PetscLayoutSetLocalSize(map, m));
265*2e1d0745SJose E. Roman     PetscCall(ISLoad(cindIS, viewer));
266*2e1d0745SJose E. Roman     PetscCall(PetscSectionLoad_HDF5_SingleField_SetConstraintIndices(s, cindIS, coffIS));
267*2e1d0745SJose E. Roman     PetscCall(ISDestroy(&coffIS));
268*2e1d0745SJose E. Roman     PetscCall(ISDestroy(&cindIS));
269*2e1d0745SJose E. Roman   }
270*2e1d0745SJose E. Roman   PetscFunctionReturn(PETSC_SUCCESS);
271*2e1d0745SJose E. Roman }
272*2e1d0745SJose E. Roman 
273*2e1d0745SJose E. Roman PetscErrorCode PetscSectionLoad_HDF5_Internal(PetscSection s, PetscViewer viewer)
274*2e1d0745SJose E. Roman {
275*2e1d0745SJose E. Roman   MPI_Comm comm;
276*2e1d0745SJose E. Roman   PetscInt N, n, numFields, f;
277*2e1d0745SJose E. Roman 
278*2e1d0745SJose E. Roman   PetscFunctionBegin;
279*2e1d0745SJose E. Roman   PetscCall(PetscObjectGetComm((PetscObject)s, &comm));
280*2e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5PushGroup(viewer, "section"));
281*2e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "numFields", PETSC_INT, NULL, (void *)&numFields));
282*2e1d0745SJose E. Roman   if (s->pStart < 0 && s->pEnd < 0) n = PETSC_DECIDE;
283*2e1d0745SJose E. Roman   else {
284*2e1d0745SJose E. Roman     PetscCheck(s->pStart == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "s->pStart must be 0 (got %" PetscInt_FMT ")", s->pStart);
285*2e1d0745SJose E. Roman     PetscCheck(s->pEnd >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "s->pEnd must be >= 0, (got %" PetscInt_FMT ")", s->pEnd);
286*2e1d0745SJose E. Roman     n = s->pEnd;
287*2e1d0745SJose E. Roman   }
288*2e1d0745SJose E. Roman   if (numFields > 0) PetscCall(PetscSectionSetNumFields(s, numFields));
289*2e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5ReadSizes(viewer, "atlasDof", NULL, &N));
290*2e1d0745SJose E. Roman   if (n == PETSC_DECIDE) PetscCall(PetscSplitOwnership(comm, &n, &N));
291*2e1d0745SJose E. Roman   PetscCall(PetscSectionSetChart(s, 0, n));
292*2e1d0745SJose E. Roman   PetscCall(PetscSectionLoad_HDF5_SingleField(s, viewer));
293*2e1d0745SJose E. Roman   for (f = 0; f < numFields; ++f) {
294*2e1d0745SJose E. Roman     char     fname[PETSC_MAX_PATH_LEN];
295*2e1d0745SJose E. Roman     char    *fieldName;
296*2e1d0745SJose E. Roman     PetscInt fieldComponents, c;
297*2e1d0745SJose E. Roman 
298*2e1d0745SJose E. Roman     PetscCall(PetscSNPrintf(fname, sizeof(fname), "field%" PetscInt_FMT, f));
299*2e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5PushGroup(viewer, fname));
300*2e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "fieldName", PETSC_STRING, NULL, &fieldName));
301*2e1d0745SJose E. Roman     PetscCall(PetscSectionSetFieldName(s, f, fieldName));
302*2e1d0745SJose E. Roman     PetscCall(PetscFree(fieldName));
303*2e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "fieldComponents", PETSC_INT, NULL, (void *)&fieldComponents));
304*2e1d0745SJose E. Roman     PetscCall(PetscSectionSetFieldComponents(s, f, fieldComponents));
305*2e1d0745SJose E. Roman     for (c = 0; c < fieldComponents; ++c) {
306*2e1d0745SJose E. Roman       char  cname[PETSC_MAX_PATH_LEN];
307*2e1d0745SJose E. Roman       char *componentName;
308*2e1d0745SJose E. Roman 
309*2e1d0745SJose E. Roman       PetscCall(PetscSNPrintf(cname, sizeof(cname), "component%" PetscInt_FMT, c));
310*2e1d0745SJose E. Roman       PetscCall(PetscViewerHDF5PushGroup(viewer, cname));
311*2e1d0745SJose E. Roman       PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "componentName", PETSC_STRING, NULL, &componentName));
312*2e1d0745SJose E. Roman       PetscCall(PetscSectionSetComponentName(s, f, c, componentName));
313*2e1d0745SJose E. Roman       PetscCall(PetscFree(componentName));
314*2e1d0745SJose E. Roman       PetscCall(PetscViewerHDF5PopGroup(viewer));
315*2e1d0745SJose E. Roman     }
316*2e1d0745SJose E. Roman     PetscCall(PetscSectionLoad_HDF5_SingleField(s->field[f], viewer));
317*2e1d0745SJose E. Roman     PetscCall(PetscViewerHDF5PopGroup(viewer));
318*2e1d0745SJose E. Roman   }
319*2e1d0745SJose E. Roman   PetscCall(PetscViewerHDF5PopGroup(viewer));
320*2e1d0745SJose E. Roman   PetscFunctionReturn(PETSC_SUCCESS);
321*2e1d0745SJose E. Roman }
322