Actual source code: plexcgns2.c

  1: #define PETSCDM_DLL
  2: #include <petsc/private/dmpleximpl.h>
  3: #include <petsc/private/viewercgnsimpl.h>

  5: #include <pcgnslib.h>
  6: #include <cgns_io.h>

  8: #if !defined(CGNS_ENUMT)
  9:   #define CGNS_ENUMT(a) a
 10: #endif
 11: #if !defined(CGNS_ENUMV)
 12:   #define CGNS_ENUMV(a) a
 13: #endif
 14: // Permute plex closure ordering to CGNS
 15: static PetscErrorCode DMPlexCGNSGetPermutation_Internal(DMPolytopeType cell_type, PetscInt closure_size, CGNS_ENUMT(ElementType_t) * element_type, const int **perm)
 16: {
 17:   CGNS_ENUMT(ElementType_t) element_type_tmp;

 19:   // https://cgns.github.io/CGNS_docs_current/sids/conv.html#unstructgrid
 20:   static const int bar_2[2]   = {0, 1};
 21:   static const int bar_3[3]   = {1, 2, 0};
 22:   static const int bar_4[4]   = {2, 3, 0, 1};
 23:   static const int bar_5[5]   = {3, 4, 0, 1, 2};
 24:   static const int tri_3[3]   = {0, 1, 2};
 25:   static const int tri_6[6]   = {3, 4, 5, 0, 1, 2};
 26:   static const int tri_10[10] = {7, 8, 9, 1, 2, 3, 4, 5, 6, 0};
 27:   static const int quad_4[4]  = {0, 1, 2, 3};
 28:   static const int quad_9[9]  = {
 29:     5, 6, 7, 8, // vertices
 30:     1, 2, 3, 4, // edges
 31:     0,          // center
 32:   };
 33:   static const int quad_16[] = {
 34:     12, 13, 14, 15,               // vertices
 35:     4,  5,  6,  7,  8, 9, 10, 11, // edges
 36:     0,  1,  3,  2,                // centers
 37:   };
 38:   static const int quad_25[] = {
 39:     21, 22, 23, 24,                                 // vertices
 40:     9,  10, 11, 12, 13, 14, 15, 16, 17, 18, 19, 20, // edges
 41:     0,  1,  2,  5,  8,  7,  6,  3,  4,              // centers
 42:   };
 43:   static const int tetra_4[4]   = {0, 2, 1, 3};
 44:   static const int tetra_10[10] = {6, 8, 7, 9, 2, 1, 0, 3, 5, 4};
 45:   static const int tetra_20[20] = {
 46:     16, 18, 17, 19,         // vertices
 47:     9,  8,  7,  6,  5,  4,  // bottom edges
 48:     10, 11, 14, 15, 13, 12, // side edges
 49:     0,  2,  3,  1,          // faces
 50:   };
 51:   static const int hexa_8[8]   = {0, 3, 2, 1, 4, 5, 6, 7};
 52:   static const int hexa_27[27] = {
 53:     19, 22, 21, 20, 23, 24, 25, 26, // vertices
 54:     10, 9,  8,  7,                  // bottom edges
 55:     16, 15, 18, 17,                 // mid edges
 56:     11, 12, 13, 14,                 // top edges
 57:     1,  3,  5,  4,  6,  2,          // faces
 58:     0,                              // center
 59:   };
 60:   static const int hexa_64[64] = {
 61:     // debug with $PETSC_ARCH/tests/dm/impls/plex/tests/ex49 -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 1,1,1 -dm_coord_petscspace_degree 3
 62:     56, 59, 58, 57, 60, 61, 62, 63, // vertices
 63:     39, 38, 37, 36, 35, 34, 33, 32, // bottom edges
 64:     51, 50, 48, 49, 52, 53, 55, 54, // mid edges; Paraview needs edge 21-22 swapped with 23-24
 65:     40, 41, 42, 43, 44, 45, 46, 47, // top edges
 66:     8,  10, 11, 9,                  // z-minus face
 67:     16, 17, 19, 18,                 // y-minus face
 68:     24, 25, 27, 26,                 // x-plus face
 69:     20, 21, 23, 22,                 // y-plus face
 70:     30, 28, 29, 31,                 // x-minus face
 71:     12, 13, 15, 14,                 // z-plus face
 72:     0,  1,  3,  2,  4,  5,  7,  6,  // center
 73:   };

 75:   PetscFunctionBegin;
 76:   element_type_tmp = CGNS_ENUMV(ElementTypeNull);
 77:   *perm            = NULL;
 78:   switch (cell_type) {
 79:   case DM_POLYTOPE_SEGMENT:
 80:     switch (closure_size) {
 81:     case 2:
 82:       element_type_tmp = CGNS_ENUMV(BAR_2);
 83:       *perm            = bar_2;
 84:     case 3:
 85:       element_type_tmp = CGNS_ENUMV(BAR_3);
 86:       *perm            = bar_3;
 87:     case 4:
 88:       element_type_tmp = CGNS_ENUMV(BAR_4);
 89:       *perm            = bar_4;
 90:       break;
 91:     case 5:
 92:       element_type_tmp = CGNS_ENUMV(BAR_5);
 93:       *perm            = bar_5;
 94:       break;
 95:     default:
 96:       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
 97:     }
 98:     break;
 99:   case DM_POLYTOPE_TRIANGLE:
100:     switch (closure_size) {
101:     case 3:
102:       element_type_tmp = CGNS_ENUMV(TRI_3);
103:       *perm            = tri_3;
104:       break;
105:     case 6:
106:       element_type_tmp = CGNS_ENUMV(TRI_6);
107:       *perm            = tri_6;
108:       break;
109:     case 10:
110:       element_type_tmp = CGNS_ENUMV(TRI_10);
111:       *perm            = tri_10;
112:       break;
113:     default:
114:       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
115:     }
116:     break;
117:   case DM_POLYTOPE_QUADRILATERAL:
118:     switch (closure_size) {
119:     case 4:
120:       element_type_tmp = CGNS_ENUMV(QUAD_4);
121:       *perm            = quad_4;
122:       break;
123:     case 9:
124:       element_type_tmp = CGNS_ENUMV(QUAD_9);
125:       *perm            = quad_9;
126:       break;
127:     case 16:
128:       element_type_tmp = CGNS_ENUMV(QUAD_16);
129:       *perm            = quad_16;
130:       break;
131:     case 25:
132:       element_type_tmp = CGNS_ENUMV(QUAD_25);
133:       *perm            = quad_25;
134:       break;
135:     default:
136:       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
137:     }
138:     break;
139:   case DM_POLYTOPE_TETRAHEDRON:
140:     switch (closure_size) {
141:     case 4:
142:       element_type_tmp = CGNS_ENUMV(TETRA_4);
143:       *perm            = tetra_4;
144:       break;
145:     case 10:
146:       element_type_tmp = CGNS_ENUMV(TETRA_10);
147:       *perm            = tetra_10;
148:       break;
149:     case 20:
150:       element_type_tmp = CGNS_ENUMV(TETRA_20);
151:       *perm            = tetra_20;
152:       break;
153:     default:
154:       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
155:     }
156:     break;
157:   case DM_POLYTOPE_HEXAHEDRON:
158:     switch (closure_size) {
159:     case 8:
160:       element_type_tmp = CGNS_ENUMV(HEXA_8);
161:       *perm            = hexa_8;
162:       break;
163:     case 27:
164:       element_type_tmp = CGNS_ENUMV(HEXA_27);
165:       *perm            = hexa_27;
166:       break;
167:     case 64:
168:       element_type_tmp = CGNS_ENUMV(HEXA_64);
169:       *perm            = hexa_64;
170:       break;
171:     default:
172:       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
173:     }
174:     break;
175:   default:
176:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
177:   }
178:   if (element_type) *element_type = element_type_tmp;
179:   PetscFunctionReturn(PETSC_SUCCESS);
180: }

182: /*
183:   Input Parameters:
184: + cellType  - The CGNS-defined element type

186:   Output Parameters:
187: + dmcelltype  - The equivalent DMPolytopeType for the cellType
188: . numCorners - Number of corners of the polytope
189: - dim - The topological dimension of the polytope

191: CGNS elements defined in: https://cgns.github.io/CGNS_docs_current/sids/conv.html#unstructgrid
192: */
193: static inline PetscErrorCode CGNSElementTypeGetTopologyInfo(CGNS_ENUMT(ElementType_t) cellType, DMPolytopeType *dmcelltype, PetscInt *numCorners, PetscInt *dim)
194: {
195:   DMPolytopeType _dmcelltype;

197:   PetscFunctionBeginUser;
198:   switch (cellType) {
199:   case CGNS_ENUMV(BAR_2):
200:   case CGNS_ENUMV(BAR_3):
201:   case CGNS_ENUMV(BAR_4):
202:   case CGNS_ENUMV(BAR_5):
203:     _dmcelltype = DM_POLYTOPE_SEGMENT;
204:     break;
205:   case CGNS_ENUMV(TRI_3):
206:   case CGNS_ENUMV(TRI_6):
207:   case CGNS_ENUMV(TRI_9):
208:   case CGNS_ENUMV(TRI_10):
209:   case CGNS_ENUMV(TRI_12):
210:   case CGNS_ENUMV(TRI_15):
211:     _dmcelltype = DM_POLYTOPE_TRIANGLE;
212:     break;
213:   case CGNS_ENUMV(QUAD_4):
214:   case CGNS_ENUMV(QUAD_8):
215:   case CGNS_ENUMV(QUAD_9):
216:   case CGNS_ENUMV(QUAD_12):
217:   case CGNS_ENUMV(QUAD_16):
218:   case CGNS_ENUMV(QUAD_P4_16):
219:   case CGNS_ENUMV(QUAD_25):
220:     _dmcelltype = DM_POLYTOPE_QUADRILATERAL;
221:     break;
222:   case CGNS_ENUMV(TETRA_4):
223:   case CGNS_ENUMV(TETRA_10):
224:   case CGNS_ENUMV(TETRA_16):
225:   case CGNS_ENUMV(TETRA_20):
226:   case CGNS_ENUMV(TETRA_22):
227:   case CGNS_ENUMV(TETRA_34):
228:   case CGNS_ENUMV(TETRA_35):
229:     _dmcelltype = DM_POLYTOPE_TETRAHEDRON;
230:     break;
231:   case CGNS_ENUMV(PYRA_5):
232:   case CGNS_ENUMV(PYRA_13):
233:   case CGNS_ENUMV(PYRA_14):
234:   case CGNS_ENUMV(PYRA_21):
235:   case CGNS_ENUMV(PYRA_29):
236:   case CGNS_ENUMV(PYRA_P4_29):
237:   case CGNS_ENUMV(PYRA_30):
238:   case CGNS_ENUMV(PYRA_50):
239:   case CGNS_ENUMV(PYRA_55):
240:     _dmcelltype = DM_POLYTOPE_PYRAMID;
241:     break;
242:   case CGNS_ENUMV(PENTA_6):
243:   case CGNS_ENUMV(PENTA_15):
244:   case CGNS_ENUMV(PENTA_18):
245:   case CGNS_ENUMV(PENTA_24):
246:   case CGNS_ENUMV(PENTA_33):
247:   case CGNS_ENUMV(PENTA_38):
248:   case CGNS_ENUMV(PENTA_40):
249:   case CGNS_ENUMV(PENTA_66):
250:   case CGNS_ENUMV(PENTA_75):
251:     _dmcelltype = DM_POLYTOPE_TRI_PRISM;
252:     break;
253:   case CGNS_ENUMV(HEXA_8):
254:   case CGNS_ENUMV(HEXA_20):
255:   case CGNS_ENUMV(HEXA_27):
256:   case CGNS_ENUMV(HEXA_32):
257:   case CGNS_ENUMV(HEXA_44):
258:   case CGNS_ENUMV(HEXA_56):
259:   case CGNS_ENUMV(HEXA_64):
260:   case CGNS_ENUMV(HEXA_98):
261:   case CGNS_ENUMV(HEXA_125):
262:     _dmcelltype = DM_POLYTOPE_HEXAHEDRON;
263:     break;
264:   case CGNS_ENUMV(MIXED):
265:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid CGNS ElementType_t: MIXED");
266:   default:
267:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid CGNS ElementType_t: %d", (int)cellType);
268:   }

270:   if (dmcelltype) *dmcelltype = _dmcelltype;
271:   if (numCorners) *numCorners = DMPolytopeTypeGetNumVertices(_dmcelltype);
272:   if (dim) *dim = DMPolytopeTypeGetDim(_dmcelltype);
273:   PetscFunctionReturn(PETSC_SUCCESS);
274: }

276: /*
277:   Input Parameters:
278: + cellType  - The CGNS-defined cell type

280:   Output Parameters:
281: + numClosure - Number of nodes that define the function space on the cell
282: - pOrder - The polynomial order of the cell

284: CGNS elements defined in: https://cgns.github.io/CGNS_docs_current/sids/conv.html#unstructgrid

286: Note: we only support "full" elements, ie. not seredipity elements
287: */
288: static inline PetscErrorCode CGNSElementTypeGetDiscretizationInfo(CGNS_ENUMT(ElementType_t) cellType, PetscInt *numClosure, PetscInt *pOrder)
289: {
290:   PetscInt _numClosure, _pOrder;

292:   PetscFunctionBeginUser;
293:   switch (cellType) {
294:   case CGNS_ENUMV(BAR_2):
295:     _numClosure = 2;
296:     _pOrder     = 1;
297:     break;
298:   case CGNS_ENUMV(BAR_3):
299:     _numClosure = 3;
300:     _pOrder     = 2;
301:     break;
302:   case CGNS_ENUMV(BAR_4):
303:     _numClosure = 4;
304:     _pOrder     = 3;
305:     break;
306:   case CGNS_ENUMV(BAR_5):
307:     _numClosure = 5;
308:     _pOrder     = 4;
309:     break;
310:   case CGNS_ENUMV(TRI_3):
311:     _numClosure = 3;
312:     _pOrder     = 1;
313:     break;
314:   case CGNS_ENUMV(TRI_6):
315:     _numClosure = 6;
316:     _pOrder     = 2;
317:     break;
318:   case CGNS_ENUMV(TRI_10):
319:     _numClosure = 10;
320:     _pOrder     = 3;
321:     break;
322:   case CGNS_ENUMV(TRI_15):
323:     _numClosure = 15;
324:     _pOrder     = 4;
325:     break;
326:   case CGNS_ENUMV(QUAD_4):
327:     _numClosure = 4;
328:     _pOrder     = 1;
329:     break;
330:   case CGNS_ENUMV(QUAD_9):
331:     _numClosure = 9;
332:     _pOrder     = 2;
333:     break;
334:   case CGNS_ENUMV(QUAD_16):
335:     _numClosure = 16;
336:     _pOrder     = 3;
337:     break;
338:   case CGNS_ENUMV(QUAD_25):
339:     _numClosure = 25;
340:     _pOrder     = 4;
341:     break;
342:   case CGNS_ENUMV(TETRA_4):
343:     _numClosure = 4;
344:     _pOrder     = 1;
345:     break;
346:   case CGNS_ENUMV(TETRA_10):
347:     _numClosure = 10;
348:     _pOrder     = 2;
349:     break;
350:   case CGNS_ENUMV(TETRA_20):
351:     _numClosure = 20;
352:     _pOrder     = 3;
353:     break;
354:   case CGNS_ENUMV(TETRA_35):
355:     _numClosure = 35;
356:     _pOrder     = 4;
357:     break;
358:   case CGNS_ENUMV(PYRA_5):
359:     _numClosure = 5;
360:     _pOrder     = 1;
361:     break;
362:   case CGNS_ENUMV(PYRA_14):
363:     _numClosure = 14;
364:     _pOrder     = 2;
365:     break;
366:   case CGNS_ENUMV(PYRA_30):
367:     _numClosure = 30;
368:     _pOrder     = 3;
369:     break;
370:   case CGNS_ENUMV(PYRA_55):
371:     _numClosure = 55;
372:     _pOrder     = 4;
373:     break;
374:   case CGNS_ENUMV(PENTA_6):
375:     _numClosure = 6;
376:     _pOrder     = 1;
377:     break;
378:   case CGNS_ENUMV(PENTA_18):
379:     _numClosure = 18;
380:     _pOrder     = 2;
381:     break;
382:   case CGNS_ENUMV(PENTA_40):
383:     _numClosure = 40;
384:     _pOrder     = 3;
385:     break;
386:   case CGNS_ENUMV(PENTA_75):
387:     _numClosure = 75;
388:     _pOrder     = 4;
389:     break;
390:   case CGNS_ENUMV(HEXA_8):
391:     _numClosure = 8;
392:     _pOrder     = 1;
393:     break;
394:   case CGNS_ENUMV(HEXA_27):
395:     _numClosure = 27;
396:     _pOrder     = 2;
397:     break;
398:   case CGNS_ENUMV(HEXA_64):
399:     _numClosure = 64;
400:     _pOrder     = 3;
401:     break;
402:   case CGNS_ENUMV(HEXA_125):
403:     _numClosure = 125;
404:     _pOrder     = 4;
405:     break;
406:   case CGNS_ENUMV(MIXED):
407:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid CGNS ElementType_t: MIXED");
408:   default:
409:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported or Invalid cell type %d", (int)cellType);
410:   }
411:   if (numClosure) *numClosure = _numClosure;
412:   if (pOrder) *pOrder = _pOrder;
413:   PetscFunctionReturn(PETSC_SUCCESS);
414: }

416: static PetscErrorCode PetscCGNSDataType(PetscDataType pd, CGNS_ENUMT(DataType_t) * cd)
417: {
418:   PetscFunctionBegin;
419:   switch (pd) {
420:   case PETSC_FLOAT:
421:     *cd = CGNS_ENUMV(RealSingle);
422:     break;
423:   case PETSC_DOUBLE:
424:     *cd = CGNS_ENUMV(RealDouble);
425:     break;
426:   case PETSC_COMPLEX:
427:     *cd = CGNS_ENUMV(ComplexDouble);
428:     break;
429:   default:
430:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Data type %s", PetscDataTypes[pd]);
431:   }
432:   PetscFunctionReturn(PETSC_SUCCESS);
433: }

435: PetscErrorCode DMPlexCreateCGNSFromFile_Internal(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
436: {
437:   int       cgid                = -1;
438:   PetscBool use_parallel_viewer = PETSC_FALSE;

440:   PetscFunctionBegin;
441:   PetscAssertPointer(filename, 2);
442:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_plex_cgns_parallel", &use_parallel_viewer, NULL));

444:   if (use_parallel_viewer) {
445:     PetscCallCGNS(cgp_mpi_comm(comm));
446:     PetscCallCGNS(cgp_open(filename, CG_MODE_READ, &cgid));
447:     PetscCheck(cgid > 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "cgp_open(\"%s\",...) did not return a valid file ID", filename);
448:     PetscCall(DMPlexCreateCGNS(comm, cgid, interpolate, dm));
449:     PetscCallCGNS(cgp_close(cgid));
450:   } else {
451:     PetscCallCGNS(cg_open(filename, CG_MODE_READ, &cgid));
452:     PetscCheck(cgid > 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "cg_open(\"%s\",...) did not return a valid file ID", filename);
453:     PetscCall(DMPlexCreateCGNS(comm, cgid, interpolate, dm));
454:     PetscCallCGNS(cg_close(cgid));
455:   }
456:   PetscFunctionReturn(PETSC_SUCCESS);
457: }

459: PetscErrorCode DMPlexCreateCGNS_Internal_Serial(MPI_Comm comm, PetscInt cgid, PetscBool interpolate, DM *dm)
460: {
461:   PetscMPIInt  num_proc, rank;
462:   DM           cdm;
463:   DMLabel      label;
464:   PetscSection coordSection;
465:   Vec          coordinates;
466:   PetscScalar *coords;
467:   PetscInt    *cellStart, *vertStart, v;
468:   PetscInt     labelIdRange[2], labelId;
469:   /* Read from file */
470:   char      basename[CGIO_MAX_NAME_LENGTH + 1];
471:   char      buffer[CGIO_MAX_NAME_LENGTH + 1];
472:   int       dim = 0, physDim = 0, coordDim = 0, numVertices = 0, numCells = 0;
473:   int       nzones = 0;
474:   const int B      = 1; // Only support single base

476:   PetscFunctionBegin;
477:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
478:   PetscCallMPI(MPI_Comm_size(comm, &num_proc));
479:   PetscCall(DMCreate(comm, dm));
480:   PetscCall(DMSetType(*dm, DMPLEX));

482:   /* Open CGNS II file and read basic information on rank 0, then broadcast to all processors */
483:   if (rank == 0) {
484:     int nbases, z;

486:     PetscCallCGNS(cg_nbases(cgid, &nbases));
487:     PetscCheck(nbases <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single base, not %d", nbases);
488:     PetscCallCGNS(cg_base_read(cgid, B, basename, &dim, &physDim));
489:     PetscCallCGNS(cg_nzones(cgid, B, &nzones));
490:     PetscCall(PetscCalloc2(nzones + 1, &cellStart, nzones + 1, &vertStart));
491:     for (z = 1; z <= nzones; ++z) {
492:       cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */

494:       PetscCallCGNS(cg_zone_read(cgid, B, z, buffer, sizes));
495:       numVertices += sizes[0];
496:       numCells += sizes[1];
497:       cellStart[z] += sizes[1] + cellStart[z - 1];
498:       vertStart[z] += sizes[0] + vertStart[z - 1];
499:     }
500:     for (z = 1; z <= nzones; ++z) vertStart[z] += numCells;
501:     coordDim = dim;
502:   }
503:   PetscCallMPI(MPI_Bcast(basename, CGIO_MAX_NAME_LENGTH + 1, MPI_CHAR, 0, comm));
504:   PetscCallMPI(MPI_Bcast(&dim, 1, MPI_INT, 0, comm));
505:   PetscCallMPI(MPI_Bcast(&coordDim, 1, MPI_INT, 0, comm));
506:   PetscCallMPI(MPI_Bcast(&nzones, 1, MPI_INT, 0, comm));

508:   PetscCall(PetscObjectSetName((PetscObject)*dm, basename));
509:   PetscCall(DMSetDimension(*dm, dim));
510:   PetscCall(DMCreateLabel(*dm, "celltype"));
511:   PetscCall(DMPlexSetChart(*dm, 0, numCells + numVertices));

513:   /* Read zone information */
514:   if (rank == 0) {
515:     int z, c, c_loc;

517:     /* Read the cell set connectivity table and build mesh topology
518:        CGNS standard requires that cells in a zone be numbered sequentially and be pairwise disjoint. */
519:     /* First set sizes */
520:     for (z = 1, c = 0; z <= nzones; ++z) {
521:       CGNS_ENUMT(ZoneType_t) zonetype;
522:       int nsections;
523:       CGNS_ENUMT(ElementType_t) cellType;
524:       cgsize_t       start, end;
525:       int            nbndry, parentFlag;
526:       PetscInt       numCorners, pOrder;
527:       DMPolytopeType ctype;
528:       const int      S = 1; // Only support single section

530:       PetscCallCGNS(cg_zone_type(cgid, B, z, &zonetype));
531:       PetscCheck(zonetype != CGNS_ENUMV(Structured), PETSC_COMM_SELF, PETSC_ERR_LIB, "Can only handle Unstructured zones for CGNS");
532:       PetscCallCGNS(cg_nsections(cgid, B, z, &nsections));
533:       PetscCheck(nsections <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single section, not %d", nsections);
534:       PetscCallCGNS(cg_section_read(cgid, B, z, S, buffer, &cellType, &start, &end, &nbndry, &parentFlag));
535:       if (cellType == CGNS_ENUMV(MIXED)) {
536:         cgsize_t elementDataSize, *elements;
537:         PetscInt off;

539:         PetscCallCGNS(cg_ElementDataSize(cgid, B, z, S, &elementDataSize));
540:         PetscCall(PetscMalloc1(elementDataSize, &elements));
541:         PetscCallCGNS(cg_poly_elements_read(cgid, B, z, S, elements, NULL, NULL));
542:         for (c_loc = start, off = 0; c_loc <= end; ++c_loc, ++c) {
543:           PetscCall(CGNSElementTypeGetTopologyInfo(elements[off], &ctype, &numCorners, NULL));
544:           PetscCall(CGNSElementTypeGetDiscretizationInfo(elements[off], NULL, &pOrder));
545:           PetscCheck(pOrder == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Serial CGNS reader only supports first order elements, not %" PetscInt_FMT " order", pOrder);
546:           PetscCall(DMPlexSetConeSize(*dm, c, numCorners));
547:           PetscCall(DMPlexSetCellType(*dm, c, ctype));
548:           off += numCorners + 1;
549:         }
550:         PetscCall(PetscFree(elements));
551:       } else {
552:         PetscCall(CGNSElementTypeGetTopologyInfo(cellType, &ctype, &numCorners, NULL));
553:         PetscCall(CGNSElementTypeGetDiscretizationInfo(cellType, NULL, &pOrder));
554:         PetscCheck(pOrder == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Serial CGNS reader only supports first order elements, not %" PetscInt_FMT " order", pOrder);
555:         for (c_loc = start; c_loc <= end; ++c_loc, ++c) {
556:           PetscCall(DMPlexSetConeSize(*dm, c, numCorners));
557:           PetscCall(DMPlexSetCellType(*dm, c, ctype));
558:         }
559:       }
560:     }
561:     for (v = numCells; v < numCells + numVertices; ++v) PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT));
562:   }

564:   PetscCall(DMSetUp(*dm));

566:   PetscCall(DMCreateLabel(*dm, "zone"));
567:   if (rank == 0) {
568:     int z, c, c_loc, v_loc;

570:     PetscCall(DMGetLabel(*dm, "zone", &label));
571:     for (z = 1, c = 0; z <= nzones; ++z) {
572:       CGNS_ENUMT(ElementType_t) cellType;
573:       cgsize_t  elementDataSize, *elements, start, end;
574:       int       nbndry, parentFlag;
575:       PetscInt *cone, numc, numCorners, maxCorners = 27, pOrder;
576:       const int S = 1; // Only support single section

578:       PetscCallCGNS(cg_section_read(cgid, B, z, S, buffer, &cellType, &start, &end, &nbndry, &parentFlag));
579:       numc = end - start;
580:       PetscCallCGNS(cg_ElementDataSize(cgid, B, z, S, &elementDataSize));
581:       PetscCall(PetscMalloc2(elementDataSize, &elements, maxCorners, &cone));
582:       PetscCallCGNS(cg_poly_elements_read(cgid, B, z, S, elements, NULL, NULL));
583:       if (cellType == CGNS_ENUMV(MIXED)) {
584:         /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */
585:         for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) {
586:           PetscCall(CGNSElementTypeGetTopologyInfo(elements[v], NULL, &numCorners, NULL));
587:           PetscCall(CGNSElementTypeGetDiscretizationInfo(elements[v], NULL, &pOrder));
588:           PetscCheck(pOrder == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Serial CGNS reader only supports first order elements, not %" PetscInt_FMT " order", pOrder);
589:           ++v;
590:           for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) cone[v_loc] = elements[v] + numCells - 1;
591:           PetscCall(DMPlexReorderCell(*dm, c, cone));
592:           PetscCall(DMPlexSetCone(*dm, c, cone));
593:           PetscCall(DMLabelSetValue(label, c, z));
594:         }
595:       } else {
596:         PetscCall(CGNSElementTypeGetTopologyInfo(cellType, NULL, &numCorners, NULL));
597:         PetscCall(CGNSElementTypeGetDiscretizationInfo(cellType, NULL, &pOrder));
598:         PetscCheck(pOrder == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Serial CGNS reader only supports first order elements, not %" PetscInt_FMT " order", pOrder);
599:         /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */
600:         for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) {
601:           for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) cone[v_loc] = elements[v] + numCells - 1;
602:           PetscCall(DMPlexReorderCell(*dm, c, cone));
603:           PetscCall(DMPlexSetCone(*dm, c, cone));
604:           PetscCall(DMLabelSetValue(label, c, z));
605:         }
606:       }
607:       PetscCall(PetscFree2(elements, cone));
608:     }
609:   }

611:   PetscCall(DMPlexSymmetrize(*dm));
612:   PetscCall(DMPlexStratify(*dm));
613:   if (interpolate) {
614:     DM idm;

616:     PetscCall(DMPlexInterpolate(*dm, &idm));
617:     PetscCall(DMDestroy(dm));
618:     *dm = idm;
619:   }

621:   /* Read coordinates */
622:   PetscCall(DMSetCoordinateDim(*dm, coordDim));
623:   PetscCall(DMGetCoordinateDM(*dm, &cdm));
624:   PetscCall(DMGetLocalSection(cdm, &coordSection));
625:   PetscCall(PetscSectionSetNumFields(coordSection, 1));
626:   PetscCall(PetscSectionSetFieldComponents(coordSection, 0, coordDim));
627:   PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices));
628:   for (v = numCells; v < numCells + numVertices; ++v) {
629:     PetscCall(PetscSectionSetDof(coordSection, v, dim));
630:     PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, coordDim));
631:   }
632:   PetscCall(PetscSectionSetUp(coordSection));

634:   PetscCall(DMCreateLocalVector(cdm, &coordinates));
635:   PetscCall(VecGetArray(coordinates, &coords));
636:   if (rank == 0) {
637:     PetscInt off = 0;
638:     float   *x[3];
639:     int      z, d;

641:     PetscCall(PetscMalloc3(numVertices, &x[0], numVertices, &x[1], numVertices, &x[2]));
642:     for (z = 1; z <= nzones; ++z) {
643:       CGNS_ENUMT(DataType_t) datatype;
644:       cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */
645:       cgsize_t range_min[3] = {1, 1, 1};
646:       cgsize_t range_max[3] = {1, 1, 1};
647:       int      ngrids, ncoords;

649:       PetscCallCGNS(cg_zone_read(cgid, B, z, buffer, sizes));
650:       range_max[0] = sizes[0];
651:       PetscCallCGNS(cg_ngrids(cgid, B, z, &ngrids));
652:       PetscCheck(ngrids <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single grid, not %d", ngrids);
653:       PetscCallCGNS(cg_ncoords(cgid, B, z, &ncoords));
654:       PetscCheck(ncoords == coordDim, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a coordinate array for each dimension, not %d", ncoords);
655:       for (d = 0; d < coordDim; ++d) {
656:         PetscCallCGNS(cg_coord_info(cgid, B, z, 1 + d, &datatype, buffer));
657:         PetscCallCGNS(cg_coord_read(cgid, B, z, buffer, CGNS_ENUMV(RealSingle), range_min, range_max, x[d]));
658:       }
659:       if (coordDim >= 1) {
660:         for (v = 0; v < sizes[0]; ++v) coords[(v + off) * coordDim + 0] = x[0][v];
661:       }
662:       if (coordDim >= 2) {
663:         for (v = 0; v < sizes[0]; ++v) coords[(v + off) * coordDim + 1] = x[1][v];
664:       }
665:       if (coordDim >= 3) {
666:         for (v = 0; v < sizes[0]; ++v) coords[(v + off) * coordDim + 2] = x[2][v];
667:       }
668:       off += sizes[0];
669:     }
670:     PetscCall(PetscFree3(x[0], x[1], x[2]));
671:   }
672:   PetscCall(VecRestoreArray(coordinates, &coords));

674:   PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
675:   PetscCall(VecSetBlockSize(coordinates, coordDim));
676:   PetscCall(DMSetCoordinatesLocal(*dm, coordinates));
677:   PetscCall(VecDestroy(&coordinates));

679:   /* Read boundary conditions */
680:   PetscCall(DMGetNumLabels(*dm, &labelIdRange[0]));
681:   if (rank == 0) {
682:     CGNS_ENUMT(BCType_t) bctype;
683:     CGNS_ENUMT(DataType_t) datatype;
684:     CGNS_ENUMT(PointSetType_t) pointtype;
685:     cgsize_t  *points;
686:     PetscReal *normals;
687:     int        normal[3];
688:     char      *bcname = buffer;
689:     cgsize_t   npoints, nnormals;
690:     int        z, nbc, bc, c, ndatasets;

692:     for (z = 1; z <= nzones; ++z) {
693:       PetscCallCGNS(cg_nbocos(cgid, B, z, &nbc));
694:       for (bc = 1; bc <= nbc; ++bc) {
695:         PetscCallCGNS(cg_boco_info(cgid, B, z, bc, bcname, &bctype, &pointtype, &npoints, normal, &nnormals, &datatype, &ndatasets));
696:         PetscCall(DMCreateLabel(*dm, bcname));
697:         PetscCall(DMGetLabel(*dm, bcname, &label));
698:         PetscCall(PetscMalloc2(npoints, &points, nnormals, &normals));
699:         PetscCallCGNS(cg_boco_read(cgid, B, z, bc, points, (void *)normals));
700:         if (pointtype == CGNS_ENUMV(ElementRange)) {
701:           // Range of cells: assuming half-open interval
702:           for (c = points[0]; c < points[1]; ++c) PetscCall(DMLabelSetValue(label, c - cellStart[z - 1], 1));
703:         } else if (pointtype == CGNS_ENUMV(ElementList)) {
704:           // List of cells
705:           for (c = 0; c < npoints; ++c) PetscCall(DMLabelSetValue(label, points[c] - cellStart[z - 1], 1));
706:         } else if (pointtype == CGNS_ENUMV(PointRange)) {
707:           CGNS_ENUMT(GridLocation_t) gridloc;

709:           // List of points:
710:           PetscCallCGNS(cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end"));
711:           PetscCallCGNS(cg_gridlocation_read(&gridloc));
712:           // Range of points: assuming half-open interval
713:           for (c = points[0]; c < points[1]; ++c) {
714:             if (gridloc == CGNS_ENUMV(Vertex)) PetscCall(DMLabelSetValue(label, c - vertStart[z - 1], 1));
715:             else PetscCall(DMLabelSetValue(label, c - cellStart[z - 1], 1));
716:           }
717:         } else if (pointtype == CGNS_ENUMV(PointList)) {
718:           CGNS_ENUMT(GridLocation_t) gridloc;

720:           // List of points:
721:           PetscCallCGNS(cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end"));
722:           PetscCallCGNS(cg_gridlocation_read(&gridloc));
723:           for (c = 0; c < npoints; ++c) {
724:             if (gridloc == CGNS_ENUMV(Vertex)) PetscCall(DMLabelSetValue(label, points[c] - vertStart[z - 1], 1));
725:             else PetscCall(DMLabelSetValue(label, points[c] - cellStart[z - 1], 1));
726:           }
727:         } else SETERRQ(comm, PETSC_ERR_SUP, "Unsupported point set type %d", (int)pointtype);
728:         PetscCall(PetscFree2(points, normals));
729:       }
730:     }
731:     PetscCall(PetscFree2(cellStart, vertStart));
732:   }
733:   PetscCall(DMGetNumLabels(*dm, &labelIdRange[1]));
734:   PetscCallMPI(MPI_Bcast(labelIdRange, 2, MPIU_INT, 0, comm));

736:   /* Create BC labels at all processes */
737:   for (labelId = labelIdRange[0]; labelId < labelIdRange[1]; ++labelId) {
738:     char       *labelName = buffer;
739:     size_t      len       = sizeof(buffer);
740:     const char *locName;

742:     if (rank == 0) {
743:       PetscCall(DMGetLabelByNum(*dm, labelId, &label));
744:       PetscCall(PetscObjectGetName((PetscObject)label, &locName));
745:       PetscCall(PetscStrncpy(labelName, locName, len));
746:     }
747:     PetscCallMPI(MPI_Bcast(labelName, (PetscMPIInt)len, MPIU_INT, 0, comm));
748:     PetscCallMPI(DMCreateLabel(*dm, labelName));
749:   }
750:   PetscFunctionReturn(PETSC_SUCCESS);
751: }

753: PetscErrorCode DMPlexCreateCGNS_Internal_Parallel(MPI_Comm comm, PetscInt cgid, PetscBool interpolate, DM *dm)
754: {
755:   PetscMPIInt num_proc, rank;
756:   /* Read from file */
757:   char     basename[CGIO_MAX_NAME_LENGTH + 1];
758:   char     buffer[CGIO_MAX_NAME_LENGTH + 1];
759:   int      dim = 0, physDim = 0, coordDim = 0;
760:   PetscInt NVertices = 0, NCells = 0;
761:   int      nzones = 0, nbases;
762:   int      z      = 1; // Only supports single zone files
763:   int      B      = 1; // Only supports single base

765:   PetscFunctionBegin;
766:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
767:   PetscCallMPI(MPI_Comm_size(comm, &num_proc));
768:   PetscCall(DMCreate(comm, dm));
769:   PetscCall(DMSetType(*dm, DMPLEX));

771:   PetscCallCGNS(cg_nbases(cgid, &nbases));
772:   PetscCheck(nbases <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single base, not %d", nbases);
773:   //  From the CGNS web page                 cell_dim  phys_dim (embedding space in PETSC) CGNS defines as length of spatial vectors/components)
774:   PetscCallCGNS(cg_base_read(cgid, B, basename, &dim, &physDim));
775:   PetscCallCGNS(cg_nzones(cgid, B, &nzones));
776:   PetscCheck(nzones == 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "Parallel reader limited to one zone, not %d", nzones);
777:   {
778:     cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */

780:     PetscCallCGNS(cg_zone_read(cgid, B, z, buffer, sizes));
781:     NVertices = sizes[0];
782:     NCells    = sizes[1];
783:   }

785:   PetscCall(PetscObjectSetName((PetscObject)*dm, basename));
786:   PetscCall(DMSetDimension(*dm, dim));
787:   coordDim = dim;

789:   // This is going to be a headache for mixed-topology and multiple sections. We may have to restore reading the data twice (once before  the SetChart
790:   // call to get this right but continuing for now with single section, single topology, one zone.
791:   // establish element ranges for my rank
792:   PetscInt    mystarte, myende, mystartv, myendv, myownede, myownedv;
793:   PetscLayout elem_map, vtx_map;
794:   PetscCall(PetscLayoutCreateFromSizes(comm, PETSC_DECIDE, NCells, 1, &elem_map));
795:   PetscCall(PetscLayoutCreateFromSizes(comm, PETSC_DECIDE, NVertices, 1, &vtx_map));
796:   PetscCall(PetscLayoutGetRange(elem_map, &mystarte, &myende));
797:   PetscCall(PetscLayoutGetLocalSize(elem_map, &myownede));
798:   PetscCall(PetscLayoutGetRange(vtx_map, &mystartv, &myendv));
799:   PetscCall(PetscLayoutGetLocalSize(vtx_map, &myownedv));

801:   // -- Build Plex in parallel
802:   DMPolytopeType dm_cell_type = DM_POLYTOPE_UNKNOWN;
803:   PetscInt       pOrder = 1, numClosure = -1;
804:   cgsize_t      *elements;
805:   {
806:     int        nsections;
807:     PetscInt  *elementsQ1, numCorners = -1;
808:     const int *perm;
809:     cgsize_t   start, end; // Throwaway

811:     cg_nsections(cgid, B, z, &nsections);
812:     // Read element connectivity
813:     for (int index_sect = 1; index_sect <= nsections; index_sect++) {
814:       int      nbndry, parentFlag;
815:       PetscInt cell_dim;
816:       CGNS_ENUMT(ElementType_t) cellType;

818:       PetscCallCGNS(cg_section_read(cgid, B, z, index_sect, buffer, &cellType, &start, &end, &nbndry, &parentFlag));

820:       PetscCall(CGNSElementTypeGetTopologyInfo(cellType, &dm_cell_type, &numCorners, &cell_dim));
821:       // Skip over element that are not max dimension (ie. boundary elements)
822:       if (cell_dim != dim) continue;
823:       PetscCall(CGNSElementTypeGetDiscretizationInfo(cellType, &numClosure, &pOrder));
824:       PetscCall(PetscMalloc1(myownede * numClosure, &elements));
825:       PetscCallCGNS(cgp_elements_read_data(cgid, B, z, index_sect, mystarte + 1, myende, elements));
826:       for (PetscInt v = 0; v < myownede * numClosure; ++v) elements[v] -= 1; // 0 based
827:       break;
828:     }

830:     // Create corners-only connectivity
831:     PetscCall(PetscMalloc1(myownede * numCorners, &elementsQ1));
832:     PetscCall(DMPlexCGNSGetPermutation_Internal(dm_cell_type, numCorners, NULL, &perm));
833:     for (PetscInt e = 0; e < myownede; ++e) {
834:       for (PetscInt v = 0; v < numCorners; ++v) elementsQ1[e * numCorners + perm[v]] = elements[e * numClosure + v];
835:     }

837:     // Build cell-vertex Plex
838:     PetscCall(DMPlexBuildFromCellListParallel(*dm, myownede, myownedv, NVertices, numCorners, elementsQ1, NULL, NULL));
839:     PetscCall(DMViewFromOptions(*dm, NULL, "-corner_dm_view"));
840:     PetscCall(PetscFree(elementsQ1));
841:   }

843:   if (interpolate) {
844:     DM idm;
845:     PetscCall(DMPlexInterpolate(*dm, &idm));
846:     PetscCall(DMDestroy(dm));
847:     *dm = idm;
848:     PetscCall(DMViewFromOptions(*dm, NULL, "-interpolate_dm_view"));
849:   }

851:   // -- Create SF for naive nodal-data read to elements
852:   PetscSF plex_to_cgns_sf;
853:   {
854:     PetscInt     nleaves, num_comp;
855:     PetscInt    *leaf, num_leaves = 0;
856:     PetscInt     cStart, cEnd;
857:     const int   *perm;
858:     PetscSF      cgns_to_local_sf;
859:     PetscSection local_section;
860:     PetscFE      fe;

862:     // sfNatural requires PetscSection to handle DMDistribute, so we use PetscFE to define the section
863:     // Use number of components = 1 to work with just the nodes themselves
864:     PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, 1, dm_cell_type, pOrder, PETSC_DETERMINE, &fe));
865:     PetscCall(PetscObjectSetName((PetscObject)fe, "FE for sfNatural"));
866:     PetscCall(DMAddField(*dm, NULL, (PetscObject)fe));
867:     PetscCall(DMCreateDS(*dm));
868:     PetscCall(PetscFEDestroy(&fe));

870:     PetscCall(DMGetLocalSection(*dm, &local_section));
871:     PetscCall(PetscSectionViewFromOptions(local_section, NULL, "-fe_natural_section_view"));
872:     PetscCall(PetscSectionGetFieldComponents(local_section, 0, &num_comp));
873:     PetscCall(PetscSectionGetStorageSize(local_section, &nleaves));
874:     nleaves /= num_comp;
875:     PetscCall(PetscMalloc1(nleaves, &leaf));
876:     for (PetscInt i = 0; i < nleaves; i++) leaf[i] = -1;

878:     // Get the permutation from CGNS closure numbering to PLEX closure numbering
879:     PetscCall(DMPlexCGNSGetPermutation_Internal(dm_cell_type, numClosure, NULL, &perm));
880:     PetscCall(DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd));
881:     for (PetscInt cell = cStart; cell < cEnd; ++cell) {
882:       PetscInt num_closure_dof, *closure_idx = NULL;

884:       PetscCall(DMPlexGetClosureIndices(*dm, local_section, local_section, cell, PETSC_FALSE, &num_closure_dof, &closure_idx, NULL, NULL));
885:       PetscAssert(numClosure * num_comp == num_closure_dof, comm, PETSC_ERR_PLIB, "Closure dof size does not match polytope");
886:       for (PetscInt i = 0; i < numClosure; i++) {
887:         PetscInt li = closure_idx[perm[i] * num_comp] / num_comp;
888:         if (li < 0) continue;

890:         PetscInt cgns_idx = elements[cell * numClosure + i];
891:         if (leaf[li] == -1) {
892:           leaf[li] = cgns_idx;
893:           num_leaves++;
894:         } else PetscAssert(leaf[li] == cgns_idx, PETSC_COMM_SELF, PETSC_ERR_PLIB, "leaf does not match previously set");
895:       }
896:       PetscCall(DMPlexRestoreClosureIndices(*dm, local_section, local_section, cell, PETSC_FALSE, &num_closure_dof, &closure_idx, NULL, NULL));
897:     }
898:     PetscAssert(num_leaves == nleaves, PETSC_COMM_SELF, PETSC_ERR_PLIB, "leaf count in closure does not match nleaves");
899:     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)*dm), &cgns_to_local_sf));
900:     PetscCall(PetscSFSetGraphLayout(cgns_to_local_sf, vtx_map, nleaves, NULL, PETSC_USE_POINTER, leaf));
901:     PetscCall(PetscObjectSetName((PetscObject)cgns_to_local_sf, "CGNS to Plex SF"));
902:     PetscCall(PetscSFViewFromOptions(cgns_to_local_sf, NULL, "-CGNStoPlex_sf_view"));
903:     PetscCall(PetscFree(leaf));
904:     PetscCall(PetscFree(elements));

906:     { // Convert cgns_to_local to global_to_cgns
907:       PetscSF sectionsf, cgns_to_global_sf;

909:       PetscCall(DMGetSectionSF(*dm, &sectionsf));
910:       PetscCall(PetscSFComposeInverse(cgns_to_local_sf, sectionsf, &cgns_to_global_sf));
911:       PetscCall(PetscSFDestroy(&cgns_to_local_sf));
912:       PetscCall(PetscSFCreateInverseSF(cgns_to_global_sf, &plex_to_cgns_sf));
913:       PetscCall(PetscObjectSetName((PetscObject)plex_to_cgns_sf, "Global Plex to CGNS"));
914:       PetscCall(PetscSFDestroy(&cgns_to_global_sf));
915:     }
916:   }

918:   { // -- Set coordinates for DM
919:     PetscScalar *coords;
920:     float       *x[3];
921:     double      *xd[3];
922:     PetscBool    read_with_double;
923:     PetscFE      cfe;

925:     // Setup coordinate space first. Use pOrder here for isoparametric; revisit with CPEX-0045 High Order.
926:     PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, coordDim, dm_cell_type, pOrder, PETSC_DETERMINE, &cfe));
927:     PetscCall(DMSetCoordinateDisc(*dm, cfe, PETSC_FALSE));
928:     PetscCall(PetscFEDestroy(&cfe));

930:     { // Determine if coords are written in single or double precision
931:       CGNS_ENUMT(DataType_t) datatype;

933:       PetscCallCGNS(cg_coord_info(cgid, B, z, 1, &datatype, buffer));
934:       if (datatype == CGNS_ENUMV(RealDouble)) read_with_double = PETSC_TRUE;
935:       else read_with_double = PETSC_TRUE;
936:     }

938:     // Read coords from file and set into component-major ordering
939:     if (read_with_double) PetscCall(PetscMalloc3(myownedv, &xd[0], myownedv, &xd[1], myownedv, &xd[2]));
940:     else PetscCall(PetscMalloc3(myownedv, &x[0], myownedv, &x[1], myownedv, &x[2]));
941:     PetscCall(PetscMalloc1(myownedv * coordDim, &coords));
942:     {
943:       cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */
944:       cgsize_t range_min[3] = {mystartv + 1, 1, 1};
945:       cgsize_t range_max[3] = {myendv, 1, 1};
946:       int      ngrids, ncoords;

948:       PetscCallCGNS(cg_zone_read(cgid, B, z, buffer, sizes));
949:       PetscCallCGNS(cg_ngrids(cgid, B, z, &ngrids));
950:       PetscCheck(ngrids <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single grid, not %d", ngrids);
951:       PetscCallCGNS(cg_ncoords(cgid, B, z, &ncoords));
952:       PetscCheck(ncoords == coordDim, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a coordinate array for each dimension, not %d", ncoords);
953:       if (read_with_double) {
954:         for (int d = 0; d < coordDim; ++d) PetscCallCGNS(cgp_coord_read_data(cgid, B, z, (d + 1), range_min, range_max, xd[d]));
955:         if (coordDim >= 1) {
956:           for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 0] = xd[0][v];
957:         }
958:         if (coordDim >= 2) {
959:           for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 1] = xd[1][v];
960:         }
961:         if (coordDim >= 3) {
962:           for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 2] = xd[2][v];
963:         }
964:       } else {
965:         for (int d = 0; d < coordDim; ++d) PetscCallCGNS(cgp_coord_read_data(cgid, 1, z, (d + 1), range_min, range_max, x[d]));
966:         if (coordDim >= 1) {
967:           for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 0] = x[0][v];
968:         }
969:         if (coordDim >= 2) {
970:           for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 1] = x[1][v];
971:         }
972:         if (coordDim >= 3) {
973:           for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 2] = x[2][v];
974:         }
975:       }
976:     }
977:     if (read_with_double) PetscCall(PetscFree3(xd[0], xd[1], xd[2]));
978:     else PetscCall(PetscFree3(x[0], x[1], x[2]));

980:     { // Reduce CGNS-ordered coordinate nodes to Plex ordering, and set DM's coordinates
981:       Vec          coord_global;
982:       MPI_Datatype unit;
983:       PetscScalar *coord_global_array;
984:       DM           cdm;

986:       PetscCall(DMGetCoordinateDM(*dm, &cdm));
987:       PetscCall(DMCreateGlobalVector(cdm, &coord_global));
988:       PetscCall(VecGetArrayWrite(coord_global, &coord_global_array));
989:       PetscCallMPI(MPI_Type_contiguous(coordDim, MPIU_SCALAR, &unit));
990:       PetscCallMPI(MPI_Type_commit(&unit));
991:       PetscCall(PetscSFReduceBegin(plex_to_cgns_sf, unit, coords, coord_global_array, MPI_REPLACE));
992:       PetscCall(PetscSFReduceEnd(plex_to_cgns_sf, unit, coords, coord_global_array, MPI_REPLACE));
993:       PetscCall(VecRestoreArrayWrite(coord_global, &coord_global_array));
994:       PetscCallMPI(MPI_Type_free(&unit));
995:       PetscCall(DMSetCoordinates(*dm, coord_global));
996:       PetscCall(VecDestroy(&coord_global));
997:     }
998:     PetscCall(PetscFree(coords));
999:   }

1001:   // -- Set sfNatural for solution vectors in CGNS file
1002:   // NOTE: We set sfNatural to be the map between the original CGNS ordering of nodes and the Plex ordering of nodes.
1003:   PetscCall(PetscSFViewFromOptions(plex_to_cgns_sf, NULL, "-sfNatural_init_view"));
1004:   PetscCall(DMSetNaturalSF(*dm, plex_to_cgns_sf));
1005:   PetscCall(DMSetUseNatural(*dm, PETSC_TRUE));
1006:   PetscCall(PetscSFDestroy(&plex_to_cgns_sf));

1008:   PetscCall(PetscLayoutDestroy(&elem_map));
1009:   PetscCall(PetscLayoutDestroy(&vtx_map));
1010:   PetscFunctionReturn(PETSC_SUCCESS);
1011: }

1013: // node_l2g must be freed
1014: static PetscErrorCode DMPlexCreateNodeNumbering(DM dm, PetscInt *num_local_nodes, PetscInt *num_global_nodes, PetscInt *nStart, PetscInt *nEnd, const PetscInt **node_l2g)
1015: {
1016:   PetscSection    local_section;
1017:   PetscSF         point_sf;
1018:   PetscInt        pStart, pEnd, spStart, spEnd, *points, nleaves, ncomp, *nodes;
1019:   PetscMPIInt     comm_size;
1020:   const PetscInt *ilocal, field = 0;

1022:   PetscFunctionBegin;
1023:   *num_local_nodes  = -1;
1024:   *num_global_nodes = -1;
1025:   *nStart           = -1;
1026:   *nEnd             = -1;
1027:   *node_l2g         = NULL;

1029:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &comm_size));
1030:   PetscCall(DMGetLocalSection(dm, &local_section));
1031:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
1032:   PetscCall(PetscSectionGetChart(local_section, &spStart, &spEnd));
1033:   PetscCall(DMGetPointSF(dm, &point_sf));
1034:   if (comm_size == 1) nleaves = 0;
1035:   else PetscCall(PetscSFGetGraph(point_sf, NULL, &nleaves, &ilocal, NULL));
1036:   PetscCall(PetscSectionGetFieldComponents(local_section, field, &ncomp));

1038:   PetscInt local_node = 0, owned_node = 0, owned_start = 0;
1039:   for (PetscInt p = spStart, leaf = 0; p < spEnd; p++) {
1040:     PetscInt dof;
1041:     PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof));
1042:     PetscAssert(dof % ncomp == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Field dof %" PetscInt_FMT " must be divisible by components %" PetscInt_FMT, dof, ncomp);
1043:     local_node += dof / ncomp;
1044:     if (leaf < nleaves && p == ilocal[leaf]) { // skip points owned by a different process
1045:       leaf++;
1046:     } else {
1047:       owned_node += dof / ncomp;
1048:     }
1049:   }
1050:   PetscCallMPI(MPI_Exscan(&owned_node, &owned_start, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
1051:   PetscCall(PetscMalloc1(pEnd - pStart, &points));
1052:   owned_node = 0;
1053:   for (PetscInt p = spStart, leaf = 0; p < spEnd; p++) {
1054:     if (leaf < nleaves && p == ilocal[leaf]) { // skip points owned by a different process
1055:       points[p - pStart] = -1;
1056:       leaf++;
1057:       continue;
1058:     }
1059:     PetscInt dof, offset;
1060:     PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof));
1061:     PetscCall(PetscSectionGetFieldOffset(local_section, p, field, &offset));
1062:     PetscAssert(offset % ncomp == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Field offset %" PetscInt_FMT " must be divisible by components %" PetscInt_FMT, offset, ncomp);
1063:     points[p - pStart] = owned_start + owned_node;
1064:     owned_node += dof / ncomp;
1065:   }
1066:   if (comm_size > 1) {
1067:     PetscCall(PetscSFBcastBegin(point_sf, MPIU_INT, points, points, MPI_REPLACE));
1068:     PetscCall(PetscSFBcastEnd(point_sf, MPIU_INT, points, points, MPI_REPLACE));
1069:   }

1071:   // Set up global indices for each local node
1072:   PetscCall(PetscMalloc1(local_node, &nodes));
1073:   for (PetscInt p = spStart; p < spEnd; p++) {
1074:     PetscInt dof, offset;
1075:     PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof));
1076:     PetscCall(PetscSectionGetFieldOffset(local_section, p, field, &offset));
1077:     for (PetscInt n = 0; n < dof / ncomp; n++) nodes[offset / ncomp + n] = points[p - pStart] + n;
1078:   }
1079:   PetscCall(PetscFree(points));
1080:   *num_local_nodes = local_node;
1081:   *nStart          = owned_start;
1082:   *nEnd            = owned_start + owned_node;
1083:   PetscCallMPI(MPIU_Allreduce(&owned_node, num_global_nodes, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
1084:   *node_l2g = nodes;
1085:   PetscFunctionReturn(PETSC_SUCCESS);
1086: }

1088: PetscErrorCode DMView_PlexCGNS(DM dm, PetscViewer viewer)
1089: {
1090:   PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;
1091:   PetscInt          fvGhostStart;
1092:   PetscInt          topo_dim, coord_dim, num_global_elems;
1093:   PetscInt          cStart, cEnd, num_local_nodes, num_global_nodes, nStart, nEnd;
1094:   const PetscInt   *node_l2g;
1095:   Vec               coord;
1096:   DM                colloc_dm, cdm;
1097:   PetscMPIInt       size;
1098:   const char       *dm_name;
1099:   int               base, zone;
1100:   cgsize_t          isize[3];

1102:   PetscFunctionBegin;
1103:   if (!cgv->file_num) {
1104:     PetscInt time_step;
1105:     PetscCall(DMGetOutputSequenceNumber(dm, &time_step, NULL));
1106:     PetscCall(PetscViewerCGNSFileOpen_Internal(viewer, time_step));
1107:   }
1108:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
1109:   PetscCall(DMGetDimension(dm, &topo_dim));
1110:   PetscCall(DMGetCoordinateDim(dm, &coord_dim));
1111:   PetscCall(PetscObjectGetName((PetscObject)dm, &dm_name));
1112:   PetscCallCGNS(cg_base_write(cgv->file_num, dm_name, topo_dim, coord_dim, &base));
1113:   PetscCallCGNS(cg_goto(cgv->file_num, base, NULL));
1114:   PetscCallCGNS(cg_dataclass_write(CGNS_ENUMV(NormalizedByDimensional)));

1116:   {
1117:     PetscFE        fe, fe_coord;
1118:     PetscClassId   ds_id;
1119:     PetscDualSpace dual_space, dual_space_coord;
1120:     PetscInt       num_fields, field_order = -1, field_order_coord;
1121:     PetscBool      is_simplex;
1122:     PetscCall(DMGetNumFields(dm, &num_fields));
1123:     if (num_fields > 0) {
1124:       PetscCall(DMGetField(dm, 0, NULL, (PetscObject *)&fe));
1125:       PetscCall(PetscObjectGetClassId((PetscObject)fe, &ds_id));
1126:       if (ds_id != PETSCFE_CLASSID) {
1127:         fe = NULL;
1128:         if (ds_id == PETSCFV_CLASSID) field_order = -1; // use whatever is present for coords; field will be CellCenter
1129:         else field_order = 1;                           // assume vertex-based linear elements
1130:       }
1131:     } else fe = NULL;
1132:     if (fe) {
1133:       PetscCall(PetscFEGetDualSpace(fe, &dual_space));
1134:       PetscCall(PetscDualSpaceGetOrder(dual_space, &field_order));
1135:     }
1136:     PetscCall(DMGetCoordinateDM(dm, &cdm));
1137:     PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *)&fe_coord));
1138:     {
1139:       PetscClassId id;
1140:       PetscCall(PetscObjectGetClassId((PetscObject)fe_coord, &id));
1141:       if (id != PETSCFE_CLASSID) fe_coord = NULL;
1142:     }
1143:     if (fe_coord) {
1144:       PetscCall(PetscFEGetDualSpace(fe_coord, &dual_space_coord));
1145:       PetscCall(PetscDualSpaceGetOrder(dual_space_coord, &field_order_coord));
1146:     } else field_order_coord = 1;
1147:     if (field_order > 0 && field_order != field_order_coord) {
1148:       PetscInt quadrature_order = field_order;
1149:       PetscCall(DMClone(dm, &colloc_dm));
1150:       { // Inform the new colloc_dm that it is a coordinate DM so isoperiodic affine corrections can be applied
1151:         const PetscSF *face_sfs;
1152:         PetscInt       num_face_sfs;
1153:         PetscCall(DMPlexGetIsoperiodicFaceSF(dm, &num_face_sfs, &face_sfs));
1154:         PetscCall(DMPlexSetIsoperiodicFaceSF(colloc_dm, num_face_sfs, (PetscSF *)face_sfs));
1155:         if (face_sfs) colloc_dm->periodic.setup = DMPeriodicCoordinateSetUp_Internal;
1156:       }
1157:       PetscCall(DMPlexIsSimplex(dm, &is_simplex));
1158:       PetscCall(PetscFECreateLagrange(PetscObjectComm((PetscObject)dm), topo_dim, coord_dim, is_simplex, field_order, quadrature_order, &fe));
1159:       PetscCall(DMSetCoordinateDisc(colloc_dm, fe, PETSC_TRUE));
1160:       PetscCall(PetscFEDestroy(&fe));
1161:     } else {
1162:       PetscCall(PetscObjectReference((PetscObject)dm));
1163:       colloc_dm = dm;
1164:     }
1165:   }
1166:   PetscCall(DMGetCoordinateDM(colloc_dm, &cdm));
1167:   PetscCall(DMPlexCreateNodeNumbering(cdm, &num_local_nodes, &num_global_nodes, &nStart, &nEnd, &node_l2g));
1168:   PetscCall(DMGetCoordinatesLocal(colloc_dm, &coord));
1169:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1170:   PetscCall(DMPlexGetCellTypeStratum(dm, DM_POLYTOPE_FV_GHOST, &fvGhostStart, NULL));
1171:   if (fvGhostStart >= 0) cEnd = fvGhostStart;
1172:   num_global_elems = cEnd - cStart;
1173:   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &num_global_elems, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
1174:   isize[0] = num_global_nodes;
1175:   isize[1] = num_global_elems;
1176:   isize[2] = 0;
1177:   PetscCallCGNS(cg_zone_write(cgv->file_num, base, "Zone", isize, CGNS_ENUMV(Unstructured), &zone));

1179:   {
1180:     const PetscScalar *X;
1181:     PetscScalar       *x;
1182:     int                coord_ids[3];

1184:     PetscCall(VecGetArrayRead(coord, &X));
1185:     for (PetscInt d = 0; d < coord_dim; d++) {
1186:       const double exponents[] = {0, 1, 0, 0, 0};
1187:       char         coord_name[64];
1188:       PetscCall(PetscSNPrintf(coord_name, sizeof coord_name, "Coordinate%c", 'X' + (int)d));
1189:       PetscCallCGNS(cgp_coord_write(cgv->file_num, base, zone, CGNS_ENUMV(RealDouble), coord_name, &coord_ids[d]));
1190:       PetscCallCGNS(cg_goto(cgv->file_num, base, "Zone_t", zone, "GridCoordinates", 0, coord_name, 0, NULL));
1191:       PetscCallCGNS(cg_exponents_write(CGNS_ENUMV(RealDouble), exponents));
1192:     }

1194:     int        section;
1195:     cgsize_t   e_owned, e_global, e_start, *conn = NULL;
1196:     const int *perm;
1197:     CGNS_ENUMT(ElementType_t) element_type = CGNS_ENUMV(ElementTypeNull);
1198:     {
1199:       PetscCall(PetscMalloc1(nEnd - nStart, &x));
1200:       for (PetscInt d = 0; d < coord_dim; d++) {
1201:         for (PetscInt n = 0; n < num_local_nodes; n++) {
1202:           PetscInt gn = node_l2g[n];
1203:           if (gn < nStart || nEnd <= gn) continue;
1204:           x[gn - nStart] = X[n * coord_dim + d];
1205:         }
1206:         // CGNS nodes use 1-based indexing
1207:         cgsize_t start = nStart + 1, end = nEnd;
1208:         PetscCallCGNS(cgp_coord_write_data(cgv->file_num, base, zone, coord_ids[d], &start, &end, x));
1209:       }
1210:       PetscCall(PetscFree(x));
1211:       PetscCall(VecRestoreArrayRead(coord, &X));
1212:     }

1214:     e_owned = cEnd - cStart;
1215:     if (e_owned > 0) {
1216:       DMPolytopeType cell_type;

1218:       PetscCall(DMPlexGetCellType(dm, cStart, &cell_type));
1219:       for (PetscInt i = cStart, c = 0; i < cEnd; i++) {
1220:         PetscInt closure_dof, *closure_indices, elem_size;

1222:         PetscCall(DMPlexGetClosureIndices(cdm, cdm->localSection, cdm->localSection, i, PETSC_FALSE, &closure_dof, &closure_indices, NULL, NULL));
1223:         elem_size = closure_dof / coord_dim;
1224:         if (!conn) PetscCall(PetscMalloc1(e_owned * elem_size, &conn));
1225:         PetscCall(DMPlexCGNSGetPermutation_Internal(cell_type, closure_dof / coord_dim, &element_type, &perm));
1226:         for (PetscInt j = 0; j < elem_size; j++) conn[c++] = node_l2g[closure_indices[perm[j] * coord_dim] / coord_dim] + 1;
1227:         PetscCall(DMPlexRestoreClosureIndices(cdm, cdm->localSection, cdm->localSection, i, PETSC_FALSE, &closure_dof, &closure_indices, NULL, NULL));
1228:       }
1229:     }

1231:     { // Get global element_type (for ranks that do not have owned elements)
1232:       PetscInt local_element_type, global_element_type;

1234:       local_element_type = e_owned > 0 ? element_type : -1;
1235:       PetscCallMPI(MPIU_Allreduce(&local_element_type, &global_element_type, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)viewer)));
1236:       if (local_element_type != -1) PetscCheck(local_element_type == global_element_type, PETSC_COMM_SELF, PETSC_ERR_SUP, "Ranks with different element types not supported");
1237:       element_type = global_element_type;
1238:     }
1239:     PetscCallMPI(MPIU_Allreduce(&e_owned, &e_global, 1, MPIU_CGSIZE, MPI_SUM, PetscObjectComm((PetscObject)dm)));
1240:     PetscCheck(e_global == num_global_elems, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unexpected number of elements %" PRIdCGSIZE " vs %" PetscInt_FMT, e_global, num_global_elems);
1241:     e_start = 0;
1242:     PetscCallMPI(MPI_Exscan(&e_owned, &e_start, 1, MPIU_CGSIZE, MPI_SUM, PetscObjectComm((PetscObject)dm)));
1243:     PetscCallCGNS(cgp_section_write(cgv->file_num, base, zone, "Elem", element_type, 1, e_global, 0, &section));
1244:     PetscCallCGNS(cgp_elements_write_data(cgv->file_num, base, zone, section, e_start + 1, e_start + e_owned, conn));
1245:     PetscCall(PetscFree(conn));

1247:     cgv->base            = base;
1248:     cgv->zone            = zone;
1249:     cgv->node_l2g        = node_l2g;
1250:     cgv->num_local_nodes = num_local_nodes;
1251:     cgv->nStart          = nStart;
1252:     cgv->nEnd            = nEnd;
1253:     cgv->eStart          = e_start;
1254:     cgv->eEnd            = e_start + e_owned;
1255:     if (1) {
1256:       PetscMPIInt rank;
1257:       int        *efield;
1258:       int         sol, field;
1259:       DMLabel     label;
1260:       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
1261:       PetscCall(PetscMalloc1(e_owned, &efield));
1262:       for (PetscInt i = 0; i < e_owned; i++) efield[i] = rank;
1263:       PetscCallCGNS(cg_sol_write(cgv->file_num, base, zone, "CellInfo", CGNS_ENUMV(CellCenter), &sol));
1264:       PetscCallCGNS(cgp_field_write(cgv->file_num, base, zone, sol, CGNS_ENUMV(Integer), "Rank", &field));
1265:       cgsize_t start = e_start + 1, end = e_start + e_owned;
1266:       PetscCallCGNS(cgp_field_write_data(cgv->file_num, base, zone, sol, field, &start, &end, efield));
1267:       PetscCall(DMGetLabel(dm, "Cell Sets", &label));
1268:       if (label) {
1269:         for (PetscInt c = cStart; c < cEnd; c++) {
1270:           PetscInt value;
1271:           PetscCall(DMLabelGetValue(label, c, &value));
1272:           efield[c - cStart] = value;
1273:         }
1274:         PetscCallCGNS(cgp_field_write(cgv->file_num, base, zone, sol, CGNS_ENUMV(Integer), "CellSet", &field));
1275:         PetscCallCGNS(cgp_field_write_data(cgv->file_num, base, zone, sol, field, &start, &end, efield));
1276:       }
1277:       PetscCall(PetscFree(efield));
1278:     }
1279:   }
1280:   PetscCall(DMDestroy(&colloc_dm));
1281:   PetscFunctionReturn(PETSC_SUCCESS);
1282: }

1284: PetscErrorCode VecView_Plex_Local_CGNS(Vec V, PetscViewer viewer)
1285: {
1286:   PetscViewer_CGNS  *cgv = (PetscViewer_CGNS *)viewer->data;
1287:   DM                 dm;
1288:   PetscSection       section;
1289:   PetscInt           time_step, num_fields, pStart, pEnd, fvGhostStart;
1290:   PetscReal          time, *time_slot;
1291:   size_t            *step_slot;
1292:   const PetscScalar *v;
1293:   char               solution_name[PETSC_MAX_PATH_LEN];
1294:   int                sol;

1296:   PetscFunctionBegin;
1297:   PetscCall(VecGetDM(V, &dm));
1298:   PetscCall(DMGetLocalSection(dm, &section));
1299:   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
1300:   PetscCall(DMPlexGetCellTypeStratum(dm, DM_POLYTOPE_FV_GHOST, &fvGhostStart, NULL));
1301:   if (fvGhostStart >= 0) pEnd = fvGhostStart;

1303:   if (!cgv->node_l2g) PetscCall(DMView(dm, viewer));
1304:   if (!cgv->grid_loc) { // Determine if writing to cell-centers or to nodes
1305:     PetscInt cStart, cEnd;
1306:     PetscInt local_grid_loc, global_grid_loc;

1308:     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1309:     if (fvGhostStart >= 0) cEnd = fvGhostStart;
1310:     if (cgv->num_local_nodes == 0) local_grid_loc = -1;
1311:     else if (cStart == pStart && cEnd == pEnd) local_grid_loc = CGNS_ENUMV(CellCenter);
1312:     else local_grid_loc = CGNS_ENUMV(Vertex);

1314:     PetscCallMPI(MPIU_Allreduce(&local_grid_loc, &global_grid_loc, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)viewer)));
1315:     if (local_grid_loc != -1)
1316:       PetscCheck(local_grid_loc == global_grid_loc, PETSC_COMM_SELF, PETSC_ERR_SUP, "Ranks with different grid locations not supported. Local has %" PetscInt_FMT ", allreduce returned %" PetscInt_FMT, local_grid_loc, global_grid_loc);
1317:     PetscCheck((global_grid_loc == CGNS_ENUMV(CellCenter)) || (global_grid_loc == CGNS_ENUMV(Vertex)), PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Grid location should only be CellCenter (%d) or Vertex(%d), but have %" PetscInt_FMT, CGNS_ENUMV(CellCenter), CGNS_ENUMV(Vertex), global_grid_loc);
1318:     cgv->grid_loc = global_grid_loc;
1319:   }
1320:   if (!cgv->nodal_field) {
1321:     switch (cgv->grid_loc) {
1322:     case CGNS_ENUMV(Vertex): {
1323:       PetscCall(PetscMalloc1(cgv->nEnd - cgv->nStart, &cgv->nodal_field));
1324:     } break;
1325:     case CGNS_ENUMV(CellCenter): {
1326:       PetscCall(PetscMalloc1(cgv->eEnd - cgv->eStart, &cgv->nodal_field));
1327:     } break;
1328:     default:
1329:       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Can only write for Vertex and CellCenter grid locations");
1330:     }
1331:   }
1332:   if (!cgv->output_times) PetscCall(PetscSegBufferCreate(sizeof(PetscReal), 20, &cgv->output_times));
1333:   if (!cgv->output_steps) PetscCall(PetscSegBufferCreate(sizeof(size_t), 20, &cgv->output_steps));

1335:   PetscCall(DMGetOutputSequenceNumber(dm, &time_step, &time));
1336:   if (time_step < 0) {
1337:     time_step = 0;
1338:     time      = 0.;
1339:   }
1340:   PetscCall(PetscSegBufferGet(cgv->output_times, 1, &time_slot));
1341:   *time_slot = time;
1342:   PetscCall(PetscSegBufferGet(cgv->output_steps, 1, &step_slot));
1343:   *step_slot = time_step;
1344:   PetscCall(PetscSNPrintf(solution_name, sizeof solution_name, "FlowSolution%" PetscInt_FMT, time_step));
1345:   PetscCallCGNS(cg_sol_write(cgv->file_num, cgv->base, cgv->zone, solution_name, cgv->grid_loc, &sol));
1346:   PetscCall(VecGetArrayRead(V, &v));
1347:   PetscCall(PetscSectionGetNumFields(section, &num_fields));
1348:   for (PetscInt field = 0; field < num_fields; field++) {
1349:     PetscInt    ncomp;
1350:     const char *field_name;
1351:     PetscCall(PetscSectionGetFieldName(section, field, &field_name));
1352:     PetscCall(PetscSectionGetFieldComponents(section, field, &ncomp));
1353:     for (PetscInt comp = 0; comp < ncomp; comp++) {
1354:       int         cgfield;
1355:       const char *comp_name;
1356:       char        cgns_field_name[32]; // CGNS max field name is 32
1357:       CGNS_ENUMT(DataType_t) datatype;
1358:       PetscCall(PetscSectionGetComponentName(section, field, comp, &comp_name));
1359:       if (ncomp == 1 && comp_name[0] == '0' && comp_name[1] == '\0' && field_name[0] != '\0') PetscCall(PetscStrncpy(cgns_field_name, field_name, sizeof cgns_field_name));
1360:       else if (field_name[0] == '\0') PetscCall(PetscStrncpy(cgns_field_name, comp_name, sizeof cgns_field_name));
1361:       else PetscCall(PetscSNPrintf(cgns_field_name, sizeof cgns_field_name, "%s.%s", field_name, comp_name));
1362:       PetscCall(PetscCGNSDataType(PETSC_SCALAR, &datatype));
1363:       PetscCallCGNS(cgp_field_write(cgv->file_num, cgv->base, cgv->zone, sol, datatype, cgns_field_name, &cgfield));
1364:       for (PetscInt p = pStart, n = 0; p < pEnd; p++) {
1365:         PetscInt off, dof;
1366:         PetscCall(PetscSectionGetFieldDof(section, p, field, &dof));
1367:         if (dof == 0) continue;
1368:         PetscCall(PetscSectionGetFieldOffset(section, p, field, &off));
1369:         for (PetscInt c = comp; c < dof; c += ncomp, n++) {
1370:           switch (cgv->grid_loc) {
1371:           case CGNS_ENUMV(Vertex): {
1372:             PetscInt gn = cgv->node_l2g[n];
1373:             if (gn < cgv->nStart || cgv->nEnd <= gn) continue;
1374:             cgv->nodal_field[gn - cgv->nStart] = v[off + c];
1375:           } break;
1376:           case CGNS_ENUMV(CellCenter): {
1377:             cgv->nodal_field[n] = v[off + c];
1378:           } break;
1379:           default:
1380:             SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Can only pack for Vertex and CellCenter grid locations");
1381:           }
1382:         }
1383:       }
1384:       // CGNS nodes use 1-based indexing
1385:       cgsize_t start, end;
1386:       switch (cgv->grid_loc) {
1387:       case CGNS_ENUMV(Vertex): {
1388:         start = cgv->nStart + 1;
1389:         end   = cgv->nEnd;
1390:       } break;
1391:       case CGNS_ENUMV(CellCenter): {
1392:         start = cgv->eStart + 1;
1393:         end   = cgv->eEnd;
1394:       } break;
1395:       default:
1396:         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Can only write for Vertex and CellCenter grid locations");
1397:       }
1398:       PetscCallCGNS(cgp_field_write_data(cgv->file_num, cgv->base, cgv->zone, sol, cgfield, &start, &end, cgv->nodal_field));
1399:     }
1400:   }
1401:   PetscCall(VecRestoreArrayRead(V, &v));
1402:   PetscCall(PetscViewerCGNSCheckBatch_Internal(viewer));
1403:   PetscFunctionReturn(PETSC_SUCCESS);
1404: }

1406: PetscErrorCode VecLoad_Plex_CGNS_Internal(Vec V, PetscViewer viewer)
1407: {
1408:   MPI_Comm          comm;
1409:   char              buffer[CGIO_MAX_NAME_LENGTH + 1];
1410:   PetscViewer_CGNS *cgv                 = (PetscViewer_CGNS *)viewer->data;
1411:   int               cgid                = cgv->file_num;
1412:   PetscBool         use_parallel_viewer = PETSC_FALSE;
1413:   int               z                   = 1; // Only support one zone
1414:   int               B                   = 1; // Only support one base
1415:   int               numComp;
1416:   PetscInt          V_numComps, mystartv, myendv, myownedv;

1418:   PetscFunctionBeginUser;
1419:   PetscCall(PetscObjectGetComm((PetscObject)V, &comm));

1421:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_plex_cgns_parallel", &use_parallel_viewer, NULL));
1422:   PetscCheck(use_parallel_viewer, comm, PETSC_ERR_USER_INPUT, "Cannot use VecLoad with CGNS file in serial reader; use -dm_plex_cgns_parallel to enable parallel reader");

1424:   { // Get CGNS node ownership information
1425:     int         nbases, nzones;
1426:     PetscInt    NVertices;
1427:     PetscLayout vtx_map;
1428:     cgsize_t    sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */

1430:     PetscCallCGNS(cg_nbases(cgid, &nbases));
1431:     PetscCheck(nbases <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single base, not %d", nbases);
1432:     PetscCallCGNS(cg_nzones(cgid, B, &nzones));
1433:     PetscCheck(nzones == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "limited to one zone %d", (int)nzones);

1435:     PetscCallCGNS(cg_zone_read(cgid, B, z, buffer, sizes));
1436:     NVertices = sizes[0];

1438:     PetscCall(PetscLayoutCreateFromSizes(comm, PETSC_DECIDE, NVertices, 1, &vtx_map));
1439:     PetscCall(PetscLayoutGetRange(vtx_map, &mystartv, &myendv));
1440:     PetscCall(PetscLayoutGetLocalSize(vtx_map, &myownedv));
1441:     PetscCall(PetscLayoutDestroy(&vtx_map));
1442:   }

1444:   { // -- Read data from file into Vec
1445:     PetscScalar *fields = NULL;
1446:     PetscSF      sfNatural;

1448:     { // Check compatibility between sfNatural and the data source and sink
1449:       DM       dm;
1450:       PetscInt nleaves, nroots, V_local_size;

1452:       PetscCall(VecGetDM(V, &dm));
1453:       PetscCall(DMGetNaturalSF(dm, &sfNatural));
1454:       PetscCheck(sfNatural, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "DM of Vec must have sfNatural");
1455:       PetscCall(PetscSFGetGraph(sfNatural, &nroots, &nleaves, NULL, NULL));
1456:       PetscCall(VecGetLocalSize(V, &V_local_size));
1457:       PetscCheck(nleaves == myownedv, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Number of locally owned vertices (% " PetscInt_FMT ") must match number of leaves in sfNatural (% " PetscInt_FMT ")", myownedv, nleaves);
1458:       PetscCheck(V_local_size % nroots == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Local Vec size (% " PetscInt_FMT ") not evenly divisible by number of roots in sfNatural (% " PetscInt_FMT ")", V_local_size, nroots);
1459:       V_numComps = V_local_size / nroots;
1460:     }

1462:     { // Read data into component-major ordering
1463:       int isol, numSols;
1464:       CGNS_ENUMT(DataType_t) datatype;
1465:       double *fields_CGNS;

1467:       PetscCallCGNS(cg_nsols(cgid, B, z, &numSols));
1468:       PetscCall(PetscViewerCGNSGetSolutionFileIndex_Internal(viewer, &isol));
1469:       PetscCallCGNS(cg_nfields(cgid, B, z, isol, &numComp));
1470:       PetscCheck(V_numComps == numComp, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Vec sized for  % " PetscInt_FMT " components per node, but file has %d components per node", V_numComps, numComp);

1472:       cgsize_t range_min[3] = {mystartv + 1, 1, 1};
1473:       cgsize_t range_max[3] = {myendv, 1, 1};
1474:       PetscCall(PetscMalloc1(myownedv * numComp, &fields_CGNS));
1475:       PetscCall(PetscMalloc1(myownedv * numComp, &fields));
1476:       for (int d = 0; d < numComp; ++d) {
1477:         PetscCallCGNS(cg_field_info(cgid, B, z, isol, (d + 1), &datatype, buffer));
1478:         PetscCheck(datatype == CGNS_ENUMV(RealDouble), PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "Field %s in file is not of type double", buffer);
1479:         PetscCallCGNS(cgp_field_read_data(cgid, B, z, isol, (d + 1), range_min, range_max, &fields_CGNS[d * myownedv]));
1480:       }
1481:       for (int d = 0; d < numComp; ++d) {
1482:         for (PetscInt v = 0; v < myownedv; ++v) fields[v * numComp + d] = fields_CGNS[d * myownedv + v];
1483:       }
1484:       PetscCall(PetscFree(fields_CGNS));
1485:     }

1487:     { // Reduce fields into Vec array
1488:       PetscScalar *V_array;
1489:       MPI_Datatype fieldtype;

1491:       PetscCall(VecGetArrayWrite(V, &V_array));
1492:       PetscCallMPI(MPI_Type_contiguous(numComp, MPIU_SCALAR, &fieldtype));
1493:       PetscCallMPI(MPI_Type_commit(&fieldtype));
1494:       PetscCall(PetscSFReduceBegin(sfNatural, fieldtype, fields, V_array, MPI_REPLACE));
1495:       PetscCall(PetscSFReduceEnd(sfNatural, fieldtype, fields, V_array, MPI_REPLACE));
1496:       PetscCallMPI(MPI_Type_free(&fieldtype));
1497:       PetscCall(VecRestoreArrayWrite(V, &V_array));
1498:     }
1499:     PetscCall(PetscFree(fields));
1500:   }
1501:   PetscFunctionReturn(PETSC_SUCCESS);
1502: }