Actual source code: plexcgns2.c

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

439:   PetscFunctionBegin;
440:   PetscAssertPointer(filename, 2);
441:   PetscCall(PetscViewerCGNSRegisterLogEvents_Internal());
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:     PetscCallCGNSOpen(cgp_open(filename, CG_MODE_READ, &cgid), 0, 0);
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:     PetscCallCGNSClose(cgp_close(cgid), 0, 0);
450:   } else {
451:     PetscCallCGNSOpen(cg_open(filename, CG_MODE_READ, &cgid), 0, 0);
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:     PetscCallCGNSClose(cg_close(cgid), 0, 0);
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:     PetscCallCGNSRead(cg_nbases(cgid, &nbases), *dm, 0);
487:     PetscCheck(nbases <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single base, not %d", nbases);
488:     PetscCallCGNSRead(cg_base_read(cgid, B, basename, &dim, &physDim), *dm, 0);
489:     PetscCallCGNSRead(cg_nzones(cgid, B, &nzones), *dm, 0);
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:       PetscCallCGNSRead(cg_zone_read(cgid, B, z, buffer, sizes), *dm, 0);
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:       PetscCallCGNSRead(cg_zone_type(cgid, B, z, &zonetype), *dm, 0);
531:       PetscCheck(zonetype != CGNS_ENUMV(Structured), PETSC_COMM_SELF, PETSC_ERR_LIB, "Can only handle Unstructured zones for CGNS");
532:       PetscCallCGNSRead(cg_nsections(cgid, B, z, &nsections), *dm, 0);
533:       PetscCheck(nsections <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single section, not %d", nsections);
534:       PetscCallCGNSRead(cg_section_read(cgid, B, z, S, buffer, &cellType, &start, &end, &nbndry, &parentFlag), *dm, 0);
535:       if (cellType == CGNS_ENUMV(MIXED)) {
536:         cgsize_t elementDataSize, *elements;
537:         PetscInt off;

539:         PetscCallCGNSRead(cg_ElementDataSize(cgid, B, z, S, &elementDataSize), *dm, 0);
540:         PetscCall(PetscMalloc1(elementDataSize, &elements));
541:         PetscCallCGNSReadData(cg_poly_elements_read(cgid, B, z, S, elements, NULL, NULL), *dm, 0);
542:         for (c_loc = start, off = 0; c_loc <= end; ++c_loc, ++c) {
543:           PetscCall(CGNSElementTypeGetTopologyInfo((CGNS_ENUMT(ElementType_t))elements[off], &ctype, &numCorners, NULL));
544:           PetscCall(CGNSElementTypeGetDiscretizationInfo((CGNS_ENUMT(ElementType_t))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:       PetscCallCGNSRead(cg_section_read(cgid, B, z, S, buffer, &cellType, &start, &end, &nbndry, &parentFlag), *dm, 0);
579:       numc = end - start;
580:       PetscCallCGNSRead(cg_ElementDataSize(cgid, B, z, S, &elementDataSize), *dm, 0);
581:       PetscCall(PetscMalloc2(elementDataSize, &elements, maxCorners, &cone));
582:       PetscCallCGNSReadData(cg_poly_elements_read(cgid, B, z, S, elements, NULL, NULL), *dm, 0);
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((CGNS_ENUMT(ElementType_t))elements[v], NULL, &numCorners, NULL));
587:           PetscCall(CGNSElementTypeGetDiscretizationInfo((CGNS_ENUMT(ElementType_t))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) PetscCall(DMPlexInterpolateInPlace_Internal(*dm));

615:   /* Read coordinates */
616:   PetscCall(DMSetCoordinateDim(*dm, coordDim));
617:   PetscCall(DMGetCoordinateDM(*dm, &cdm));
618:   PetscCall(DMGetLocalSection(cdm, &coordSection));
619:   PetscCall(PetscSectionSetNumFields(coordSection, 1));
620:   PetscCall(PetscSectionSetFieldComponents(coordSection, 0, coordDim));
621:   PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices));
622:   for (v = numCells; v < numCells + numVertices; ++v) {
623:     PetscCall(PetscSectionSetDof(coordSection, v, dim));
624:     PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, coordDim));
625:   }
626:   PetscCall(PetscSectionSetUp(coordSection));

628:   PetscCall(DMCreateLocalVector(cdm, &coordinates));
629:   PetscCall(VecGetArray(coordinates, &coords));
630:   if (rank == 0) {
631:     PetscInt off = 0;
632:     float   *x[3];
633:     int      z, d;

635:     PetscCall(PetscMalloc3(numVertices, &x[0], numVertices, &x[1], numVertices, &x[2]));
636:     for (z = 1; z <= nzones; ++z) {
637:       CGNS_ENUMT(DataType_t) datatype;
638:       cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */
639:       cgsize_t range_min[3] = {1, 1, 1};
640:       cgsize_t range_max[3] = {1, 1, 1};
641:       int      ngrids, ncoords;

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

668:   PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
669:   PetscCall(VecSetBlockSize(coordinates, coordDim));
670:   PetscCall(DMSetCoordinatesLocal(*dm, coordinates));
671:   PetscCall(VecDestroy(&coordinates));

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

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

703:           // List of points:
704:           PetscCallCGNS(cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end"));
705:           PetscCallCGNSRead(cg_gridlocation_read(&gridloc), *dm, 0);
706:           // Range of points: assuming half-open interval
707:           for (c = points[0]; c < points[1]; ++c) {
708:             if (gridloc == CGNS_ENUMV(Vertex)) PetscCall(DMLabelSetValue(label, c - vertStart[z - 1], 1));
709:             else PetscCall(DMLabelSetValue(label, c - cellStart[z - 1], 1));
710:           }
711:         } else if (pointtype == CGNS_ENUMV(PointList)) {
712:           CGNS_ENUMT(GridLocation_t) gridloc;

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

730:   /* Create BC labels at all processes */
731:   for (labelId = labelIdRange[0]; labelId < labelIdRange[1]; ++labelId) {
732:     char       *labelName = buffer;
733:     size_t      len       = sizeof(buffer);
734:     const char *locName;

736:     if (rank == 0) {
737:       PetscCall(DMGetLabelByNum(*dm, labelId, &label));
738:       PetscCall(PetscObjectGetName((PetscObject)label, &locName));
739:       PetscCall(PetscStrncpy(labelName, locName, len));
740:     }
741:     PetscCallMPI(MPI_Bcast(labelName, (PetscMPIInt)len, MPIU_INT, 0, comm));
742:     PetscCallMPI(DMCreateLabel(*dm, labelName));
743:   }
744:   PetscFunctionReturn(PETSC_SUCCESS);
745: }

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

759:   PetscFunctionBegin;
760:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
761:   PetscCallMPI(MPI_Comm_size(comm, &num_proc));
762:   PetscCall(DMCreate(comm, dm));
763:   PetscCall(DMSetType(*dm, DMPLEX));

765:   PetscCallCGNSRead(cg_nbases(cgid, &nbases), *dm, 0);
766:   PetscCheck(nbases <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single base, not %d", nbases);
767:   //  From the CGNS web page                 cell_dim  phys_dim (embedding space in PETSc) CGNS defines as length of spatial vectors/components)
768:   PetscCallCGNSRead(cg_base_read(cgid, B, basename, &dim, &physDim), *dm, 0);
769:   PetscCallCGNSRead(cg_nzones(cgid, B, &nzones), *dm, 0);
770:   PetscCheck(nzones == 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "Parallel reader limited to one zone, not %d", nzones);
771:   {
772:     cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */

774:     PetscCallCGNSRead(cg_zone_read(cgid, B, z, buffer, sizes), *dm, 0);
775:     NVertices = sizes[0];
776:     NCells    = sizes[1];
777:   }

779:   PetscCall(PetscObjectSetName((PetscObject)*dm, basename));
780:   PetscCall(DMSetDimension(*dm, dim));
781:   coordDim = dim;

783:   // 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
784:   // call to get this right but continuing for now with single section, single topology, one zone.
785:   // establish element ranges for my rank
786:   PetscInt    mystarte, myende, mystartv, myendv, myownede, myownedv;
787:   PetscLayout elem_map, vtx_map;
788:   PetscCall(PetscLayoutCreateFromSizes(comm, PETSC_DECIDE, NCells, 1, &elem_map));
789:   PetscCall(PetscLayoutCreateFromSizes(comm, PETSC_DECIDE, NVertices, 1, &vtx_map));
790:   PetscCall(PetscLayoutGetRange(elem_map, &mystarte, &myende));
791:   PetscCall(PetscLayoutGetLocalSize(elem_map, &myownede));
792:   PetscCall(PetscLayoutGetRange(vtx_map, &mystartv, &myendv));
793:   PetscCall(PetscLayoutGetLocalSize(vtx_map, &myownedv));

795:   // -- Build Plex in parallel
796:   DMPolytopeType dm_cell_type = DM_POLYTOPE_UNKNOWN;
797:   PetscInt       pOrder = 1, numClosure = -1;
798:   cgsize_t      *elements;
799:   {
800:     int        nsections;
801:     PetscInt  *elementsQ1, numCorners = -1;
802:     const int *perm;
803:     cgsize_t   start, end; // Throwaway

805:     cg_nsections(cgid, B, z, &nsections);
806:     // Read element connectivity
807:     for (int index_sect = 1; index_sect <= nsections; index_sect++) {
808:       int      nbndry, parentFlag;
809:       PetscInt cell_dim;
810:       CGNS_ENUMT(ElementType_t) cellType;

812:       PetscCallCGNSRead(cg_section_read(cgid, B, z, index_sect, buffer, &cellType, &start, &end, &nbndry, &parentFlag), *dm, 0);

814:       PetscCall(CGNSElementTypeGetTopologyInfo(cellType, &dm_cell_type, &numCorners, &cell_dim));
815:       // Skip over element that are not max dimension (ie. boundary elements)
816:       if (cell_dim != dim) continue;
817:       PetscCall(CGNSElementTypeGetDiscretizationInfo(cellType, &numClosure, &pOrder));
818:       PetscCall(PetscMalloc1(myownede * numClosure, &elements));
819:       PetscCallCGNSReadData(cgp_elements_read_data(cgid, B, z, index_sect, mystarte + 1, myende, elements), *dm, 0);
820:       for (PetscInt v = 0; v < myownede * numClosure; ++v) elements[v] -= 1; // 0 based
821:       break;
822:     }

824:     // Create corners-only connectivity
825:     PetscCall(PetscMalloc1(myownede * numCorners, &elementsQ1));
826:     PetscCall(DMPlexCGNSGetPermutation_Internal(dm_cell_type, numCorners, NULL, &perm));
827:     for (PetscInt e = 0; e < myownede; ++e) {
828:       for (PetscInt v = 0; v < numCorners; ++v) elementsQ1[e * numCorners + perm[v]] = elements[e * numClosure + v];
829:     }

831:     // Build cell-vertex Plex
832:     PetscCall(DMPlexBuildFromCellListParallel(*dm, myownede, myownedv, NVertices, numCorners, elementsQ1, NULL, NULL));
833:     PetscCall(DMViewFromOptions(*dm, NULL, "-corner_dm_view"));
834:     PetscCall(PetscFree(elementsQ1));
835:   }

837:   if (interpolate) PetscCall(DMPlexInterpolateInPlace_Internal(*dm));

839:   // -- Create SF for naive nodal-data read to elements
840:   PetscSF plex_to_cgns_sf;
841:   {
842:     PetscInt     nleaves, num_comp;
843:     PetscInt    *leaf, num_leaves = 0;
844:     PetscInt     cStart, cEnd;
845:     const int   *perm;
846:     PetscSF      cgns_to_local_sf;
847:     PetscSection local_section;
848:     PetscFE      fe;

850:     // sfNatural requires PetscSection to handle DMDistribute, so we use PetscFE to define the section
851:     // Use number of components = 1 to work with just the nodes themselves
852:     PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, 1, dm_cell_type, pOrder, PETSC_DETERMINE, &fe));
853:     PetscCall(PetscObjectSetName((PetscObject)fe, "FE for sfNatural"));
854:     PetscCall(DMAddField(*dm, NULL, (PetscObject)fe));
855:     PetscCall(DMCreateDS(*dm));
856:     PetscCall(PetscFEDestroy(&fe));

858:     PetscCall(DMGetLocalSection(*dm, &local_section));
859:     PetscCall(PetscSectionViewFromOptions(local_section, NULL, "-fe_natural_section_view"));
860:     PetscCall(PetscSectionGetFieldComponents(local_section, 0, &num_comp));
861:     PetscCall(PetscSectionGetStorageSize(local_section, &nleaves));
862:     nleaves /= num_comp;
863:     PetscCall(PetscMalloc1(nleaves, &leaf));
864:     for (PetscInt i = 0; i < nleaves; i++) leaf[i] = -1;

866:     // Get the permutation from CGNS closure numbering to PLEX closure numbering
867:     PetscCall(DMPlexCGNSGetPermutation_Internal(dm_cell_type, numClosure, NULL, &perm));
868:     PetscCall(DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd));
869:     for (PetscInt cell = cStart; cell < cEnd; ++cell) {
870:       PetscInt num_closure_dof, *closure_idx = NULL;

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

878:         PetscInt cgns_idx = elements[cell * numClosure + i];
879:         if (leaf[li] == -1) {
880:           leaf[li] = cgns_idx;
881:           num_leaves++;
882:         } else PetscAssert(leaf[li] == cgns_idx, PETSC_COMM_SELF, PETSC_ERR_PLIB, "leaf does not match previously set");
883:       }
884:       PetscCall(DMPlexRestoreClosureIndices(*dm, local_section, local_section, cell, PETSC_FALSE, &num_closure_dof, &closure_idx, NULL, NULL));
885:     }
886:     PetscAssert(num_leaves == nleaves, PETSC_COMM_SELF, PETSC_ERR_PLIB, "leaf count in closure does not match nleaves");
887:     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)*dm), &cgns_to_local_sf));
888:     PetscCall(PetscSFSetGraphLayout(cgns_to_local_sf, vtx_map, nleaves, NULL, PETSC_USE_POINTER, leaf));
889:     PetscCall(PetscObjectSetName((PetscObject)cgns_to_local_sf, "CGNS to Plex SF"));
890:     PetscCall(PetscSFViewFromOptions(cgns_to_local_sf, NULL, "-CGNStoPlex_sf_view"));
891:     PetscCall(PetscFree(leaf));
892:     PetscCall(PetscFree(elements));

894:     { // Convert cgns_to_local to global_to_cgns
895:       PetscSF sectionsf, cgns_to_global_sf;

897:       PetscCall(DMGetSectionSF(*dm, &sectionsf));
898:       PetscCall(PetscSFComposeInverse(cgns_to_local_sf, sectionsf, &cgns_to_global_sf));
899:       PetscCall(PetscSFDestroy(&cgns_to_local_sf));
900:       PetscCall(PetscSFCreateInverseSF(cgns_to_global_sf, &plex_to_cgns_sf));
901:       PetscCall(PetscObjectSetName((PetscObject)plex_to_cgns_sf, "Global Plex to CGNS"));
902:       PetscCall(PetscSFDestroy(&cgns_to_global_sf));
903:     }
904:   }

906:   { // -- Set coordinates for DM
907:     PetscScalar *coords;
908:     float       *x[3];
909:     double      *xd[3];
910:     PetscBool    read_with_double;
911:     PetscFE      cfe;

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

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

921:       PetscCallCGNSRead(cg_coord_info(cgid, B, z, 1, &datatype, buffer), *dm, 0);
922:       read_with_double = datatype == CGNS_ENUMV(RealDouble) ? PETSC_TRUE : PETSC_FALSE;
923:     }

925:     // Read coords from file and set into component-major ordering
926:     if (read_with_double) PetscCall(PetscMalloc3(myownedv, &xd[0], myownedv, &xd[1], myownedv, &xd[2]));
927:     else PetscCall(PetscMalloc3(myownedv, &x[0], myownedv, &x[1], myownedv, &x[2]));
928:     PetscCall(PetscMalloc1(myownedv * coordDim, &coords));
929:     {
930:       cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */
931:       cgsize_t range_min[3] = {mystartv + 1, 1, 1};
932:       cgsize_t range_max[3] = {myendv, 1, 1};
933:       int      ngrids, ncoords;

935:       PetscCallCGNSRead(cg_zone_read(cgid, B, z, buffer, sizes), *dm, 0);
936:       PetscCallCGNSRead(cg_ngrids(cgid, B, z, &ngrids), *dm, 0);
937:       PetscCheck(ngrids <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single grid, not %d", ngrids);
938:       PetscCallCGNSRead(cg_ncoords(cgid, B, z, &ncoords), *dm, 0);
939:       PetscCheck(ncoords == coordDim, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a coordinate array for each dimension, not %d", ncoords);
940:       if (read_with_double) {
941:         for (int d = 0; d < coordDim; ++d) PetscCallCGNSReadData(cgp_coord_read_data(cgid, B, z, (d + 1), range_min, range_max, xd[d]), *dm, 0);
942:         if (coordDim >= 1) {
943:           for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 0] = xd[0][v];
944:         }
945:         if (coordDim >= 2) {
946:           for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 1] = xd[1][v];
947:         }
948:         if (coordDim >= 3) {
949:           for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 2] = xd[2][v];
950:         }
951:       } else {
952:         for (int d = 0; d < coordDim; ++d) PetscCallCGNSReadData(cgp_coord_read_data(cgid, 1, z, (d + 1), range_min, range_max, x[d]), *dm, 0);
953:         if (coordDim >= 1) {
954:           for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 0] = x[0][v];
955:         }
956:         if (coordDim >= 2) {
957:           for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 1] = x[1][v];
958:         }
959:         if (coordDim >= 3) {
960:           for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 2] = x[2][v];
961:         }
962:       }
963:     }
964:     if (read_with_double) PetscCall(PetscFree3(xd[0], xd[1], xd[2]));
965:     else PetscCall(PetscFree3(x[0], x[1], x[2]));

967:     { // Reduce CGNS-ordered coordinate nodes to Plex ordering, and set DM's coordinates
968:       Vec          coord_global;
969:       MPI_Datatype unit;
970:       PetscScalar *coord_global_array;
971:       DM           cdm;

973:       PetscCall(DMGetCoordinateDM(*dm, &cdm));
974:       PetscCall(DMCreateGlobalVector(cdm, &coord_global));
975:       PetscCall(VecGetArrayWrite(coord_global, &coord_global_array));
976:       PetscCallMPI(MPI_Type_contiguous(coordDim, MPIU_SCALAR, &unit));
977:       PetscCallMPI(MPI_Type_commit(&unit));
978:       PetscCall(PetscSFReduceBegin(plex_to_cgns_sf, unit, coords, coord_global_array, MPI_REPLACE));
979:       PetscCall(PetscSFReduceEnd(plex_to_cgns_sf, unit, coords, coord_global_array, MPI_REPLACE));
980:       PetscCall(VecRestoreArrayWrite(coord_global, &coord_global_array));
981:       PetscCallMPI(MPI_Type_free(&unit));
982:       PetscCall(DMSetCoordinates(*dm, coord_global));
983:       PetscCall(VecDestroy(&coord_global));
984:     }
985:     PetscCall(PetscFree(coords));
986:   }

988:   // -- Set sfNatural for solution vectors in CGNS file
989:   // NOTE: We set sfNatural to be the map between the original CGNS ordering of nodes and the Plex ordering of nodes.
990:   PetscCall(PetscSFViewFromOptions(plex_to_cgns_sf, NULL, "-sfNatural_init_view"));
991:   PetscCall(DMSetNaturalSF(*dm, plex_to_cgns_sf));
992:   PetscCall(DMSetUseNatural(*dm, PETSC_TRUE));
993:   PetscCall(PetscSFDestroy(&plex_to_cgns_sf));

995:   PetscCall(PetscLayoutDestroy(&elem_map));
996:   PetscCall(PetscLayoutDestroy(&vtx_map));
997:   PetscFunctionReturn(PETSC_SUCCESS);
998: }

1000: // node_l2g must be freed
1001: static PetscErrorCode DMPlexCreateNodeNumbering(DM dm, PetscInt *num_local_nodes, PetscInt *num_global_nodes, PetscInt *nStart, PetscInt *nEnd, const PetscInt **node_l2g)
1002: {
1003:   PetscSection    local_section;
1004:   PetscSF         point_sf;
1005:   PetscInt        pStart, pEnd, spStart, spEnd, *points, nleaves, ncomp, *nodes;
1006:   PetscMPIInt     comm_size;
1007:   const PetscInt *ilocal, field = 0;

1009:   PetscFunctionBegin;
1010:   *num_local_nodes  = -1;
1011:   *num_global_nodes = -1;
1012:   *nStart           = -1;
1013:   *nEnd             = -1;
1014:   *node_l2g         = NULL;

1016:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &comm_size));
1017:   PetscCall(DMGetLocalSection(dm, &local_section));
1018:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
1019:   PetscCall(PetscSectionGetChart(local_section, &spStart, &spEnd));
1020:   PetscCall(DMGetPointSF(dm, &point_sf));
1021:   if (comm_size == 1) nleaves = 0;
1022:   else PetscCall(PetscSFGetGraph(point_sf, NULL, &nleaves, &ilocal, NULL));
1023:   PetscCall(PetscSectionGetFieldComponents(local_section, field, &ncomp));

1025:   PetscInt local_node = 0, owned_node = 0, owned_start = 0;
1026:   for (PetscInt p = spStart, leaf = 0; p < spEnd; p++) {
1027:     PetscInt dof;
1028:     PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof));
1029:     PetscAssert(dof % ncomp == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Field dof %" PetscInt_FMT " must be divisible by components %" PetscInt_FMT, dof, ncomp);
1030:     local_node += dof / ncomp;
1031:     if (leaf < nleaves && p == ilocal[leaf]) { // skip points owned by a different process
1032:       leaf++;
1033:     } else {
1034:       owned_node += dof / ncomp;
1035:     }
1036:   }
1037:   PetscCallMPI(MPI_Exscan(&owned_node, &owned_start, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
1038:   PetscCall(PetscMalloc1(pEnd - pStart, &points));
1039:   owned_node = 0;
1040:   for (PetscInt p = spStart, leaf = 0; p < spEnd; p++) {
1041:     if (leaf < nleaves && p == ilocal[leaf]) { // skip points owned by a different process
1042:       points[p - pStart] = -1;
1043:       leaf++;
1044:       continue;
1045:     }
1046:     PetscInt dof, offset;
1047:     PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof));
1048:     PetscCall(PetscSectionGetFieldOffset(local_section, p, field, &offset));
1049:     PetscAssert(offset % ncomp == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Field offset %" PetscInt_FMT " must be divisible by components %" PetscInt_FMT, offset, ncomp);
1050:     points[p - pStart] = owned_start + owned_node;
1051:     owned_node += dof / ncomp;
1052:   }
1053:   if (comm_size > 1) {
1054:     PetscCall(PetscSFBcastBegin(point_sf, MPIU_INT, points, points, MPI_REPLACE));
1055:     PetscCall(PetscSFBcastEnd(point_sf, MPIU_INT, points, points, MPI_REPLACE));
1056:   }

1058:   // Set up global indices for each local node
1059:   PetscCall(PetscMalloc1(local_node, &nodes));
1060:   for (PetscInt p = spStart; p < spEnd; p++) {
1061:     PetscInt dof, offset;
1062:     PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof));
1063:     PetscCall(PetscSectionGetFieldOffset(local_section, p, field, &offset));
1064:     for (PetscInt n = 0; n < dof / ncomp; n++) nodes[offset / ncomp + n] = points[p - pStart] + n;
1065:   }
1066:   PetscCall(PetscFree(points));
1067:   *num_local_nodes = local_node;
1068:   *nStart          = owned_start;
1069:   *nEnd            = owned_start + owned_node;
1070:   PetscCallMPI(MPIU_Allreduce(&owned_node, num_global_nodes, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
1071:   *node_l2g = nodes;
1072:   PetscFunctionReturn(PETSC_SUCCESS);
1073: }

1075: PetscErrorCode DMView_PlexCGNS(DM dm, PetscViewer viewer)
1076: {
1077:   PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;
1078:   PetscInt          fvGhostStart;
1079:   PetscInt          topo_dim, coord_dim, num_global_elems;
1080:   PetscInt          cStart, cEnd, num_local_nodes, num_global_nodes, nStart, nEnd;
1081:   const PetscInt   *node_l2g;
1082:   Vec               coord;
1083:   DM                colloc_dm, cdm;
1084:   PetscMPIInt       size;
1085:   const char       *dm_name;
1086:   int               base, zone;
1087:   cgsize_t          isize[3];

1089:   PetscFunctionBegin;
1090:   if (!cgv->file_num) {
1091:     PetscInt time_step;
1092:     PetscCall(DMGetOutputSequenceNumber(dm, &time_step, NULL));
1093:     PetscCall(PetscViewerCGNSFileOpen_Internal(viewer, time_step));
1094:   }
1095:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
1096:   PetscCall(DMGetDimension(dm, &topo_dim));
1097:   PetscCall(DMGetCoordinateDim(dm, &coord_dim));
1098:   PetscCall(PetscObjectGetName((PetscObject)dm, &dm_name));
1099:   PetscCallCGNSWrite(cg_base_write(cgv->file_num, dm_name, topo_dim, coord_dim, &base), dm, viewer);
1100:   PetscCallCGNS(cg_goto(cgv->file_num, base, NULL));
1101:   PetscCallCGNSWrite(cg_dataclass_write(CGNS_ENUMV(NormalizedByDimensional)), dm, viewer);

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

1166:   {
1167:     const PetscScalar *X;
1168:     PetscScalar       *x;
1169:     int                coord_ids[3];

1171:     PetscCall(VecGetArrayRead(coord, &X));
1172:     for (PetscInt d = 0; d < coord_dim; d++) {
1173:       const double exponents[] = {0, 1, 0, 0, 0};
1174:       char         coord_name[64];
1175:       PetscCall(PetscSNPrintf(coord_name, sizeof coord_name, "Coordinate%c", 'X' + (int)d));
1176:       PetscCallCGNSWrite(cgp_coord_write(cgv->file_num, base, zone, CGNS_ENUMV(RealDouble), coord_name, &coord_ids[d]), dm, viewer);
1177:       PetscCallCGNS(cg_goto(cgv->file_num, base, "Zone_t", zone, "GridCoordinates", 0, coord_name, 0, NULL));
1178:       PetscCallCGNSWrite(cg_exponents_write(CGNS_ENUMV(RealDouble), exponents), dm, viewer);
1179:     }

1181:     int        section;
1182:     cgsize_t   e_owned, e_global, e_start, *conn = NULL;
1183:     const int *perm;
1184:     CGNS_ENUMT(ElementType_t) element_type = CGNS_ENUMV(ElementTypeNull);
1185:     {
1186:       PetscCall(PetscMalloc1(nEnd - nStart, &x));
1187:       for (PetscInt d = 0; d < coord_dim; d++) {
1188:         for (PetscInt n = 0; n < num_local_nodes; n++) {
1189:           PetscInt gn = node_l2g[n];
1190:           if (gn < nStart || nEnd <= gn) continue;
1191:           x[gn - nStart] = X[n * coord_dim + d];
1192:         }
1193:         // CGNS nodes use 1-based indexing
1194:         cgsize_t start = nStart + 1, end = nEnd;
1195:         PetscCallCGNSWriteData(cgp_coord_write_data(cgv->file_num, base, zone, coord_ids[d], &start, &end, x), dm, viewer);
1196:       }
1197:       PetscCall(PetscFree(x));
1198:       PetscCall(VecRestoreArrayRead(coord, &X));
1199:     }

1201:     e_owned = cEnd - cStart;
1202:     if (e_owned > 0) {
1203:       DMPolytopeType cell_type;

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

1209:         PetscCall(DMPlexGetClosureIndices(cdm, cdm->localSection, cdm->localSection, i, PETSC_FALSE, &closure_dof, &closure_indices, NULL, NULL));
1210:         elem_size = closure_dof / coord_dim;
1211:         if (!conn) PetscCall(PetscMalloc1(e_owned * elem_size, &conn));
1212:         PetscCall(DMPlexCGNSGetPermutation_Internal(cell_type, closure_dof / coord_dim, &element_type, &perm));
1213:         for (PetscInt j = 0; j < elem_size; j++) conn[c++] = node_l2g[closure_indices[perm[j] * coord_dim] / coord_dim] + 1;
1214:         PetscCall(DMPlexRestoreClosureIndices(cdm, cdm->localSection, cdm->localSection, i, PETSC_FALSE, &closure_dof, &closure_indices, NULL, NULL));
1215:       }
1216:     }

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

1221:       local_element_type = e_owned > 0 ? element_type : -1;
1222:       PetscCallMPI(MPIU_Allreduce(&local_element_type, &global_element_type, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)viewer)));
1223:       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");
1224:       element_type = (CGNS_ENUMT(ElementType_t))global_element_type;
1225:     }
1226:     PetscCallMPI(MPIU_Allreduce(&e_owned, &e_global, 1, MPIU_CGSIZE, MPI_SUM, PetscObjectComm((PetscObject)dm)));
1227:     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);
1228:     e_start = 0;
1229:     PetscCallMPI(MPI_Exscan(&e_owned, &e_start, 1, MPIU_CGSIZE, MPI_SUM, PetscObjectComm((PetscObject)dm)));
1230:     PetscCallCGNSWrite(cgp_section_write(cgv->file_num, base, zone, "Elem", element_type, 1, e_global, 0, &section), dm, viewer);
1231:     PetscCallCGNSWriteData(cgp_elements_write_data(cgv->file_num, base, zone, section, e_start + 1, e_start + e_owned, conn), dm, viewer);
1232:     PetscCall(PetscFree(conn));

1234:     cgv->base            = base;
1235:     cgv->zone            = zone;
1236:     cgv->node_l2g        = node_l2g;
1237:     cgv->num_local_nodes = num_local_nodes;
1238:     cgv->nStart          = nStart;
1239:     cgv->nEnd            = nEnd;
1240:     cgv->eStart          = e_start;
1241:     cgv->eEnd            = e_start + e_owned;
1242:     if (1) {
1243:       PetscMPIInt rank;
1244:       int        *efield;
1245:       int         sol, field;
1246:       DMLabel     label;
1247:       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
1248:       PetscCall(PetscMalloc1(e_owned, &efield));
1249:       for (PetscInt i = 0; i < e_owned; i++) efield[i] = rank;
1250:       PetscCallCGNSWrite(cg_sol_write(cgv->file_num, base, zone, "CellInfo", CGNS_ENUMV(CellCenter), &sol), dm, viewer);
1251:       PetscCallCGNSWrite(cgp_field_write(cgv->file_num, base, zone, sol, CGNS_ENUMV(Integer), "Rank", &field), dm, viewer);
1252:       cgsize_t start = e_start + 1, end = e_start + e_owned;
1253:       PetscCallCGNSWriteData(cgp_field_write_data(cgv->file_num, base, zone, sol, field, &start, &end, efield), dm, viewer);
1254:       PetscCall(DMGetLabel(dm, "Cell Sets", &label));
1255:       if (label) {
1256:         for (PetscInt c = cStart; c < cEnd; c++) {
1257:           PetscInt value;
1258:           PetscCall(DMLabelGetValue(label, c, &value));
1259:           efield[c - cStart] = value;
1260:         }
1261:         PetscCallCGNSWrite(cgp_field_write(cgv->file_num, base, zone, sol, CGNS_ENUMV(Integer), "CellSet", &field), dm, viewer);
1262:         PetscCallCGNSWriteData(cgp_field_write_data(cgv->file_num, base, zone, sol, field, &start, &end, efield), dm, viewer);
1263:       }
1264:       PetscCall(PetscFree(efield));
1265:     }
1266:   }
1267:   PetscCall(DMDestroy(&colloc_dm));
1268:   PetscFunctionReturn(PETSC_SUCCESS);
1269: }

1271: PetscErrorCode VecView_Plex_Local_CGNS(Vec V, PetscViewer viewer)
1272: {
1273:   PetscViewer_CGNS  *cgv = (PetscViewer_CGNS *)viewer->data;
1274:   DM                 dm;
1275:   PetscSection       section;
1276:   PetscInt           time_step, num_fields, pStart, pEnd, fvGhostStart;
1277:   PetscReal          time, *time_slot;
1278:   size_t            *step_slot;
1279:   const PetscScalar *v;
1280:   char               solution_name[PETSC_MAX_PATH_LEN];
1281:   int                sol;

1283:   PetscFunctionBegin;
1284:   PetscCall(VecGetDM(V, &dm));
1285:   PetscCall(DMGetLocalSection(dm, &section));
1286:   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
1287:   PetscCall(DMPlexGetCellTypeStratum(dm, DM_POLYTOPE_FV_GHOST, &fvGhostStart, NULL));
1288:   if (fvGhostStart >= 0) pEnd = fvGhostStart;

1290:   if (!cgv->node_l2g) PetscCall(DMView(dm, viewer));
1291:   if (!cgv->grid_loc) { // Determine if writing to cell-centers or to nodes
1292:     PetscInt cStart, cEnd;
1293:     PetscInt local_grid_loc, global_grid_loc;

1295:     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1296:     if (fvGhostStart >= 0) cEnd = fvGhostStart;
1297:     if (cgv->num_local_nodes == 0) local_grid_loc = -1;
1298:     else if (cStart == pStart && cEnd == pEnd) local_grid_loc = CGNS_ENUMV(CellCenter);
1299:     else local_grid_loc = CGNS_ENUMV(Vertex);

1301:     PetscCallMPI(MPIU_Allreduce(&local_grid_loc, &global_grid_loc, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)viewer)));
1302:     if (local_grid_loc != -1)
1303:       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);
1304:     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);
1305:     cgv->grid_loc = (CGNS_ENUMT(GridLocation_t))global_grid_loc;
1306:   }
1307:   if (!cgv->nodal_field) {
1308:     switch (cgv->grid_loc) {
1309:     case CGNS_ENUMV(Vertex): {
1310:       PetscCall(PetscMalloc1(cgv->nEnd - cgv->nStart, &cgv->nodal_field));
1311:     } break;
1312:     case CGNS_ENUMV(CellCenter): {
1313:       PetscCall(PetscMalloc1(cgv->eEnd - cgv->eStart, &cgv->nodal_field));
1314:     } break;
1315:     default:
1316:       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Can only write for Vertex and CellCenter grid locations");
1317:     }
1318:   }
1319:   if (!cgv->output_times) PetscCall(PetscSegBufferCreate(sizeof(PetscReal), 20, &cgv->output_times));
1320:   if (!cgv->output_steps) PetscCall(PetscSegBufferCreate(sizeof(size_t), 20, &cgv->output_steps));

1322:   PetscCall(DMGetOutputSequenceNumber(dm, &time_step, &time));
1323:   if (time_step < 0) {
1324:     time_step = 0;
1325:     time      = 0.;
1326:   }
1327:   // Avoid "Duplicate child name found" error by not writing an already-written solution.
1328:   // This usually occurs when a solution is written and then diverges on the very next timestep.
1329:   if (time_step == cgv->previous_output_step) PetscFunctionReturn(PETSC_SUCCESS);

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

1397: PetscErrorCode VecLoad_Plex_CGNS_Internal(Vec V, PetscViewer viewer)
1398: {
1399:   MPI_Comm          comm;
1400:   char              buffer[CGIO_MAX_NAME_LENGTH + 1];
1401:   PetscViewer_CGNS *cgv                 = (PetscViewer_CGNS *)viewer->data;
1402:   int               cgid                = cgv->file_num;
1403:   PetscBool         use_parallel_viewer = PETSC_FALSE;
1404:   int               z                   = 1; // Only support one zone
1405:   int               B                   = 1; // Only support one base
1406:   int               numComp;
1407:   PetscInt          V_numComps, mystartv, myendv, myownedv;

1409:   PetscFunctionBeginUser;
1410:   PetscCall(PetscObjectGetComm((PetscObject)V, &comm));

1412:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_plex_cgns_parallel", &use_parallel_viewer, NULL));
1413:   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");

1415:   { // Get CGNS node ownership information
1416:     int         nbases, nzones;
1417:     PetscInt    NVertices;
1418:     PetscLayout vtx_map;
1419:     cgsize_t    sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */

1421:     PetscCallCGNSRead(cg_nbases(cgid, &nbases), V, viewer);
1422:     PetscCheck(nbases <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single base, not %d", nbases);
1423:     PetscCallCGNSRead(cg_nzones(cgid, B, &nzones), V, viewer);
1424:     PetscCheck(nzones == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "limited to one zone %d", (int)nzones);

1426:     PetscCallCGNSRead(cg_zone_read(cgid, B, z, buffer, sizes), V, viewer);
1427:     NVertices = sizes[0];

1429:     PetscCall(PetscLayoutCreateFromSizes(comm, PETSC_DECIDE, NVertices, 1, &vtx_map));
1430:     PetscCall(PetscLayoutGetRange(vtx_map, &mystartv, &myendv));
1431:     PetscCall(PetscLayoutGetLocalSize(vtx_map, &myownedv));
1432:     PetscCall(PetscLayoutDestroy(&vtx_map));
1433:   }

1435:   { // -- Read data from file into Vec
1436:     PetscScalar *fields = NULL;
1437:     PetscSF      sfNatural;

1439:     { // Check compatibility between sfNatural and the data source and sink
1440:       DM       dm;
1441:       PetscInt nleaves, nroots, V_local_size;

1443:       PetscCall(VecGetDM(V, &dm));
1444:       PetscCall(DMGetNaturalSF(dm, &sfNatural));
1445:       PetscCheck(sfNatural, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "DM of Vec must have sfNatural");
1446:       PetscCall(PetscSFGetGraph(sfNatural, &nroots, &nleaves, NULL, NULL));
1447:       PetscCall(VecGetLocalSize(V, &V_local_size));
1448:       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);
1449:       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);
1450:       V_numComps = V_local_size / nroots;
1451:     }

1453:     { // Read data into component-major ordering
1454:       int isol, numSols;
1455:       CGNS_ENUMT(DataType_t) datatype;
1456:       double *fields_CGNS;

1458:       PetscCallCGNSRead(cg_nsols(cgid, B, z, &numSols), V, viewer);
1459:       PetscCall(PetscViewerCGNSGetSolutionFileIndex_Internal(viewer, &isol));
1460:       PetscCallCGNSRead(cg_nfields(cgid, B, z, isol, &numComp), V, viewer);
1461:       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);

1463:       cgsize_t range_min[3] = {mystartv + 1, 1, 1};
1464:       cgsize_t range_max[3] = {myendv, 1, 1};
1465:       PetscCall(PetscMalloc1(myownedv * numComp, &fields_CGNS));
1466:       PetscCall(PetscMalloc1(myownedv * numComp, &fields));
1467:       for (int d = 0; d < numComp; ++d) {
1468:         PetscCallCGNSRead(cg_field_info(cgid, B, z, isol, (d + 1), &datatype, buffer), V, viewer);
1469:         PetscCheck(datatype == CGNS_ENUMV(RealDouble), PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "Field %s in file is not of type double", buffer);
1470:         PetscCallCGNSReadData(cgp_field_read_data(cgid, B, z, isol, (d + 1), range_min, range_max, &fields_CGNS[d * myownedv]), V, viewer);
1471:       }
1472:       for (int d = 0; d < numComp; ++d) {
1473:         for (PetscInt v = 0; v < myownedv; ++v) fields[v * numComp + d] = fields_CGNS[d * myownedv + v];
1474:       }
1475:       PetscCall(PetscFree(fields_CGNS));
1476:     }

1478:     { // Reduce fields into Vec array
1479:       PetscScalar *V_array;
1480:       MPI_Datatype fieldtype;

1482:       PetscCall(VecGetArrayWrite(V, &V_array));
1483:       PetscCallMPI(MPI_Type_contiguous(numComp, MPIU_SCALAR, &fieldtype));
1484:       PetscCallMPI(MPI_Type_commit(&fieldtype));
1485:       PetscCall(PetscSFReduceBegin(sfNatural, fieldtype, fields, V_array, MPI_REPLACE));
1486:       PetscCall(PetscSFReduceEnd(sfNatural, fieldtype, fields, V_array, MPI_REPLACE));
1487:       PetscCallMPI(MPI_Type_free(&fieldtype));
1488:       PetscCall(VecRestoreArrayWrite(V, &V_array));
1489:     }
1490:     PetscCall(PetscFree(fields));
1491:   }
1492:   PetscFunctionReturn(PETSC_SUCCESS);
1493: }