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
15: PetscErrorCode DMPlexCreateCGNSFromFile_Internal(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
16: {
17: PetscMPIInt rank;
18: int cgid = -1;
20: PetscFunctionBegin;
21: PetscAssertPointer(filename, 2);
22: PetscCallMPI(MPI_Comm_rank(comm, &rank));
23: if (rank == 0) {
24: PetscCallCGNS(cg_open(filename, CG_MODE_READ, &cgid));
25: PetscCheck(cgid > 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "cg_open(\"%s\",...) did not return a valid file ID", filename);
26: }
27: PetscCall(DMPlexCreateCGNS(comm, cgid, interpolate, dm));
28: if (rank == 0) PetscCallCGNS(cg_close(cgid));
29: PetscFunctionReturn(PETSC_SUCCESS);
30: }
32: PetscErrorCode DMPlexCreateCGNS_Internal(MPI_Comm comm, PetscInt cgid, PetscBool interpolate, DM *dm)
33: {
34: PetscMPIInt num_proc, rank;
35: DM cdm;
36: DMLabel label;
37: PetscSection coordSection;
38: Vec coordinates;
39: PetscScalar *coords;
40: PetscInt *cellStart, *vertStart, v;
41: PetscInt labelIdRange[2], labelId;
42: /* Read from file */
43: char basename[CGIO_MAX_NAME_LENGTH + 1];
44: char buffer[CGIO_MAX_NAME_LENGTH + 1];
45: int dim = 0, physDim = 0, coordDim = 0, numVertices = 0, numCells = 0;
46: int nzones = 0;
48: PetscFunctionBegin;
49: PetscCallMPI(MPI_Comm_rank(comm, &rank));
50: PetscCallMPI(MPI_Comm_size(comm, &num_proc));
51: PetscCall(DMCreate(comm, dm));
52: PetscCall(DMSetType(*dm, DMPLEX));
54: /* Open CGNS II file and read basic information on rank 0, then broadcast to all processors */
55: if (rank == 0) {
56: int nbases, z;
58: PetscCallCGNS(cg_nbases(cgid, &nbases));
59: PetscCheck(nbases <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single base, not %d", nbases);
60: PetscCallCGNS(cg_base_read(cgid, 1, basename, &dim, &physDim));
61: PetscCallCGNS(cg_nzones(cgid, 1, &nzones));
62: PetscCall(PetscCalloc2(nzones + 1, &cellStart, nzones + 1, &vertStart));
63: for (z = 1; z <= nzones; ++z) {
64: cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */
66: PetscCallCGNS(cg_zone_read(cgid, 1, z, buffer, sizes));
67: numVertices += sizes[0];
68: numCells += sizes[1];
69: cellStart[z] += sizes[1] + cellStart[z - 1];
70: vertStart[z] += sizes[0] + vertStart[z - 1];
71: }
72: for (z = 1; z <= nzones; ++z) vertStart[z] += numCells;
73: coordDim = dim;
74: }
75: PetscCallMPI(MPI_Bcast(basename, CGIO_MAX_NAME_LENGTH + 1, MPI_CHAR, 0, comm));
76: PetscCallMPI(MPI_Bcast(&dim, 1, MPI_INT, 0, comm));
77: PetscCallMPI(MPI_Bcast(&coordDim, 1, MPI_INT, 0, comm));
78: PetscCallMPI(MPI_Bcast(&nzones, 1, MPI_INT, 0, comm));
80: PetscCall(PetscObjectSetName((PetscObject)*dm, basename));
81: PetscCall(DMSetDimension(*dm, dim));
82: PetscCall(DMCreateLabel(*dm, "celltype"));
83: PetscCall(DMPlexSetChart(*dm, 0, numCells + numVertices));
85: /* Read zone information */
86: if (rank == 0) {
87: int z, c, c_loc;
89: /* Read the cell set connectivity table and build mesh topology
90: CGNS standard requires that cells in a zone be numbered sequentially and be pairwise disjoint. */
91: /* First set sizes */
92: for (z = 1, c = 0; z <= nzones; ++z) {
93: CGNS_ENUMT(ZoneType_t) zonetype;
94: int nsections;
95: CGNS_ENUMT(ElementType_t) cellType;
96: cgsize_t start, end;
97: int nbndry, parentFlag;
98: PetscInt numCorners;
99: DMPolytopeType ctype;
101: PetscCallCGNS(cg_zone_type(cgid, 1, z, &zonetype));
102: PetscCheck(zonetype != CGNS_ENUMV(Structured), PETSC_COMM_SELF, PETSC_ERR_LIB, "Can only handle Unstructured zones for CGNS");
103: PetscCallCGNS(cg_nsections(cgid, 1, z, &nsections));
104: PetscCheck(nsections <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single section, not %d", nsections);
105: PetscCallCGNS(cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag));
106: /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */
107: if (cellType == CGNS_ENUMV(MIXED)) {
108: cgsize_t elementDataSize, *elements;
109: PetscInt off;
111: PetscCallCGNS(cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize));
112: PetscCall(PetscMalloc1(elementDataSize, &elements));
113: PetscCallCGNS(cg_poly_elements_read(cgid, 1, z, 1, elements, NULL, NULL));
114: for (c_loc = start, off = 0; c_loc <= end; ++c_loc, ++c) {
115: switch (elements[off]) {
116: case CGNS_ENUMV(BAR_2):
117: numCorners = 2;
118: ctype = DM_POLYTOPE_SEGMENT;
119: break;
120: case CGNS_ENUMV(TRI_3):
121: numCorners = 3;
122: ctype = DM_POLYTOPE_TRIANGLE;
123: break;
124: case CGNS_ENUMV(QUAD_4):
125: numCorners = 4;
126: ctype = DM_POLYTOPE_QUADRILATERAL;
127: break;
128: case CGNS_ENUMV(TETRA_4):
129: numCorners = 4;
130: ctype = DM_POLYTOPE_TETRAHEDRON;
131: break;
132: case CGNS_ENUMV(PENTA_6):
133: numCorners = 6;
134: ctype = DM_POLYTOPE_TRI_PRISM;
135: break;
136: case CGNS_ENUMV(HEXA_8):
137: numCorners = 8;
138: ctype = DM_POLYTOPE_HEXAHEDRON;
139: break;
140: default:
141: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int)elements[off]);
142: }
143: PetscCall(DMPlexSetConeSize(*dm, c, numCorners));
144: PetscCall(DMPlexSetCellType(*dm, c, ctype));
145: off += numCorners + 1;
146: }
147: PetscCall(PetscFree(elements));
148: } else {
149: switch (cellType) {
150: case CGNS_ENUMV(BAR_2):
151: numCorners = 2;
152: ctype = DM_POLYTOPE_SEGMENT;
153: break;
154: case CGNS_ENUMV(TRI_3):
155: numCorners = 3;
156: ctype = DM_POLYTOPE_TRIANGLE;
157: break;
158: case CGNS_ENUMV(QUAD_4):
159: numCorners = 4;
160: ctype = DM_POLYTOPE_QUADRILATERAL;
161: break;
162: case CGNS_ENUMV(TETRA_4):
163: numCorners = 4;
164: ctype = DM_POLYTOPE_TETRAHEDRON;
165: break;
166: case CGNS_ENUMV(PENTA_6):
167: numCorners = 6;
168: ctype = DM_POLYTOPE_TRI_PRISM;
169: break;
170: case CGNS_ENUMV(HEXA_8):
171: numCorners = 8;
172: ctype = DM_POLYTOPE_HEXAHEDRON;
173: break;
174: default:
175: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int)cellType);
176: }
177: for (c_loc = start; c_loc <= end; ++c_loc, ++c) {
178: PetscCall(DMPlexSetConeSize(*dm, c, numCorners));
179: PetscCall(DMPlexSetCellType(*dm, c, ctype));
180: }
181: }
182: }
183: for (v = numCells; v < numCells + numVertices; ++v) PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT));
184: }
186: PetscCall(DMSetUp(*dm));
188: PetscCall(DMCreateLabel(*dm, "zone"));
189: if (rank == 0) {
190: int z, c, c_loc, v_loc;
192: PetscCall(DMGetLabel(*dm, "zone", &label));
193: for (z = 1, c = 0; z <= nzones; ++z) {
194: CGNS_ENUMT(ElementType_t) cellType;
195: cgsize_t elementDataSize, *elements, start, end;
196: int nbndry, parentFlag;
197: PetscInt *cone, numc, numCorners, maxCorners = 27;
199: PetscCallCGNS(cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag));
200: numc = end - start;
201: /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */
202: PetscCallCGNS(cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize));
203: PetscCall(PetscMalloc2(elementDataSize, &elements, maxCorners, &cone));
204: PetscCallCGNS(cg_poly_elements_read(cgid, 1, z, 1, elements, NULL, NULL));
205: if (cellType == CGNS_ENUMV(MIXED)) {
206: /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */
207: for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) {
208: switch (elements[v]) {
209: case CGNS_ENUMV(BAR_2):
210: numCorners = 2;
211: break;
212: case CGNS_ENUMV(TRI_3):
213: numCorners = 3;
214: break;
215: case CGNS_ENUMV(QUAD_4):
216: numCorners = 4;
217: break;
218: case CGNS_ENUMV(TETRA_4):
219: numCorners = 4;
220: break;
221: case CGNS_ENUMV(PENTA_6):
222: numCorners = 6;
223: break;
224: case CGNS_ENUMV(HEXA_8):
225: numCorners = 8;
226: break;
227: default:
228: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int)elements[v]);
229: }
230: ++v;
231: for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) cone[v_loc] = elements[v] + numCells - 1;
232: PetscCall(DMPlexReorderCell(*dm, c, cone));
233: PetscCall(DMPlexSetCone(*dm, c, cone));
234: PetscCall(DMLabelSetValue(label, c, z));
235: }
236: } else {
237: switch (cellType) {
238: case CGNS_ENUMV(BAR_2):
239: numCorners = 2;
240: break;
241: case CGNS_ENUMV(TRI_3):
242: numCorners = 3;
243: break;
244: case CGNS_ENUMV(QUAD_4):
245: numCorners = 4;
246: break;
247: case CGNS_ENUMV(TETRA_4):
248: numCorners = 4;
249: break;
250: case CGNS_ENUMV(PENTA_6):
251: numCorners = 6;
252: break;
253: case CGNS_ENUMV(HEXA_8):
254: numCorners = 8;
255: break;
256: default:
257: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int)cellType);
258: }
259: /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */
260: for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) {
261: for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) cone[v_loc] = elements[v] + numCells - 1;
262: PetscCall(DMPlexReorderCell(*dm, c, cone));
263: PetscCall(DMPlexSetCone(*dm, c, cone));
264: PetscCall(DMLabelSetValue(label, c, z));
265: }
266: }
267: PetscCall(PetscFree2(elements, cone));
268: }
269: }
271: PetscCall(DMPlexSymmetrize(*dm));
272: PetscCall(DMPlexStratify(*dm));
273: if (interpolate) {
274: DM idm;
276: PetscCall(DMPlexInterpolate(*dm, &idm));
277: PetscCall(DMDestroy(dm));
278: *dm = idm;
279: }
281: /* Read coordinates */
282: PetscCall(DMSetCoordinateDim(*dm, coordDim));
283: PetscCall(DMGetCoordinateDM(*dm, &cdm));
284: PetscCall(DMGetLocalSection(cdm, &coordSection));
285: PetscCall(PetscSectionSetNumFields(coordSection, 1));
286: PetscCall(PetscSectionSetFieldComponents(coordSection, 0, coordDim));
287: PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices));
288: for (v = numCells; v < numCells + numVertices; ++v) {
289: PetscCall(PetscSectionSetDof(coordSection, v, dim));
290: PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, coordDim));
291: }
292: PetscCall(PetscSectionSetUp(coordSection));
294: PetscCall(DMCreateLocalVector(cdm, &coordinates));
295: PetscCall(VecGetArray(coordinates, &coords));
296: if (rank == 0) {
297: PetscInt off = 0;
298: float *x[3];
299: int z, d;
301: PetscCall(PetscMalloc3(numVertices, &x[0], numVertices, &x[1], numVertices, &x[2]));
302: for (z = 1; z <= nzones; ++z) {
303: CGNS_ENUMT(DataType_t) datatype;
304: cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */
305: cgsize_t range_min[3] = {1, 1, 1};
306: cgsize_t range_max[3] = {1, 1, 1};
307: int ngrids, ncoords;
309: PetscCallCGNS(cg_zone_read(cgid, 1, z, buffer, sizes));
310: range_max[0] = sizes[0];
311: PetscCallCGNS(cg_ngrids(cgid, 1, z, &ngrids));
312: PetscCheck(ngrids <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single grid, not %d", ngrids);
313: PetscCallCGNS(cg_ncoords(cgid, 1, z, &ncoords));
314: PetscCheck(ncoords == coordDim, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a coordinate array for each dimension, not %d", ncoords);
315: for (d = 0; d < coordDim; ++d) {
316: PetscCallCGNS(cg_coord_info(cgid, 1, z, 1 + d, &datatype, buffer));
317: PetscCallCGNS(cg_coord_read(cgid, 1, z, buffer, CGNS_ENUMV(RealSingle), range_min, range_max, x[d]));
318: }
319: if (coordDim >= 1) {
320: for (v = 0; v < sizes[0]; ++v) coords[(v + off) * coordDim + 0] = x[0][v];
321: }
322: if (coordDim >= 2) {
323: for (v = 0; v < sizes[0]; ++v) coords[(v + off) * coordDim + 1] = x[1][v];
324: }
325: if (coordDim >= 3) {
326: for (v = 0; v < sizes[0]; ++v) coords[(v + off) * coordDim + 2] = x[2][v];
327: }
328: off += sizes[0];
329: }
330: PetscCall(PetscFree3(x[0], x[1], x[2]));
331: }
332: PetscCall(VecRestoreArray(coordinates, &coords));
334: PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
335: PetscCall(VecSetBlockSize(coordinates, coordDim));
336: PetscCall(DMSetCoordinatesLocal(*dm, coordinates));
337: PetscCall(VecDestroy(&coordinates));
339: /* Read boundary conditions */
340: PetscCall(DMGetNumLabels(*dm, &labelIdRange[0]));
341: if (rank == 0) {
342: CGNS_ENUMT(BCType_t) bctype;
343: CGNS_ENUMT(DataType_t) datatype;
344: CGNS_ENUMT(PointSetType_t) pointtype;
345: cgsize_t *points;
346: PetscReal *normals;
347: int normal[3];
348: char *bcname = buffer;
349: cgsize_t npoints, nnormals;
350: int z, nbc, bc, c, ndatasets;
352: for (z = 1; z <= nzones; ++z) {
353: PetscCallCGNS(cg_nbocos(cgid, 1, z, &nbc));
354: for (bc = 1; bc <= nbc; ++bc) {
355: PetscCallCGNS(cg_boco_info(cgid, 1, z, bc, bcname, &bctype, &pointtype, &npoints, normal, &nnormals, &datatype, &ndatasets));
356: PetscCall(DMCreateLabel(*dm, bcname));
357: PetscCall(DMGetLabel(*dm, bcname, &label));
358: PetscCall(PetscMalloc2(npoints, &points, nnormals, &normals));
359: PetscCallCGNS(cg_boco_read(cgid, 1, z, bc, points, (void *)normals));
360: if (pointtype == CGNS_ENUMV(ElementRange)) {
361: /* Range of cells: assuming half-open interval since the documentation sucks */
362: for (c = points[0]; c < points[1]; ++c) PetscCall(DMLabelSetValue(label, c - cellStart[z - 1], 1));
363: } else if (pointtype == CGNS_ENUMV(ElementList)) {
364: /* List of cells */
365: for (c = 0; c < npoints; ++c) PetscCall(DMLabelSetValue(label, points[c] - cellStart[z - 1], 1));
366: } else if (pointtype == CGNS_ENUMV(PointRange)) {
367: CGNS_ENUMT(GridLocation_t) gridloc;
369: /* List of points: Oh please, someone get the CGNS developers away from a computer. This is unconscionable. */
370: PetscCallCGNS(cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end"));
371: PetscCallCGNS(cg_gridlocation_read(&gridloc));
372: /* Range of points: assuming half-open interval since the documentation sucks */
373: for (c = points[0]; c < points[1]; ++c) {
374: if (gridloc == CGNS_ENUMV(Vertex)) PetscCall(DMLabelSetValue(label, c - vertStart[z - 1], 1));
375: else PetscCall(DMLabelSetValue(label, c - cellStart[z - 1], 1));
376: }
377: } else if (pointtype == CGNS_ENUMV(PointList)) {
378: CGNS_ENUMT(GridLocation_t) gridloc;
380: /* List of points: Oh please, someone get the CGNS developers away from a computer. This is unconscionable. */
381: PetscCallCGNS(cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end"));
382: PetscCallCGNS(cg_gridlocation_read(&gridloc));
383: for (c = 0; c < npoints; ++c) {
384: if (gridloc == CGNS_ENUMV(Vertex)) PetscCall(DMLabelSetValue(label, points[c] - vertStart[z - 1], 1));
385: else PetscCall(DMLabelSetValue(label, points[c] - cellStart[z - 1], 1));
386: }
387: } else SETERRQ(comm, PETSC_ERR_SUP, "Unsupported point set type %d", (int)pointtype);
388: PetscCall(PetscFree2(points, normals));
389: }
390: }
391: PetscCall(PetscFree2(cellStart, vertStart));
392: }
393: PetscCall(DMGetNumLabels(*dm, &labelIdRange[1]));
394: PetscCallMPI(MPI_Bcast(labelIdRange, 2, MPIU_INT, 0, comm));
396: /* Create BC labels at all processes */
397: for (labelId = labelIdRange[0]; labelId < labelIdRange[1]; ++labelId) {
398: char *labelName = buffer;
399: size_t len = sizeof(buffer);
400: const char *locName;
402: if (rank == 0) {
403: PetscCall(DMGetLabelByNum(*dm, labelId, &label));
404: PetscCall(PetscObjectGetName((PetscObject)label, &locName));
405: PetscCall(PetscStrncpy(labelName, locName, len));
406: }
407: PetscCallMPI(MPI_Bcast(labelName, (PetscMPIInt)len, MPIU_INT, 0, comm));
408: PetscCallMPI(DMCreateLabel(*dm, labelName));
409: }
410: PetscFunctionReturn(PETSC_SUCCESS);
411: }
413: // Permute plex closure ordering to CGNS
414: static PetscErrorCode DMPlexCGNSGetPermutation_Internal(DMPolytopeType cell_type, PetscInt closure_size, CGNS_ENUMT(ElementType_t) * element_type, const int **perm)
415: {
416: // https://cgns.github.io/CGNS_docs_current/sids/conv.html#unst_example
417: static const int bar_2[2] = {0, 1};
418: static const int bar_3[3] = {1, 2, 0};
419: static const int bar_4[4] = {2, 3, 0, 1};
420: static const int bar_5[5] = {3, 4, 0, 1, 2};
421: static const int tri_3[3] = {0, 1, 2};
422: static const int tri_6[6] = {3, 4, 5, 0, 1, 2};
423: static const int tri_10[10] = {7, 8, 9, 1, 2, 3, 4, 5, 6, 0};
424: static const int quad_4[4] = {0, 1, 2, 3};
425: static const int quad_9[9] = {
426: 5, 6, 7, 8, // vertices
427: 1, 2, 3, 4, // edges
428: 0, // center
429: };
430: static const int quad_16[] = {
431: 12, 13, 14, 15, // vertices
432: 4, 5, 6, 7, 8, 9, 10, 11, // edges
433: 0, 1, 3, 2, // centers
434: };
435: static const int quad_25[] = {
436: 21, 22, 23, 24, // vertices
437: 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19, 20, // edges
438: 0, 1, 2, 5, 8, 7, 6, 3, 4, // centers
439: };
440: static const int tetra_4[4] = {0, 2, 1, 3};
441: static const int tetra_10[10] = {6, 8, 7, 9, 2, 1, 0, 3, 5, 4};
442: static const int tetra_20[20] = {
443: 16, 18, 17, 19, // vertices
444: 9, 8, 7, 6, 5, 4, // bottom edges
445: 10, 11, 14, 15, 13, 12, // side edges
446: 0, 2, 3, 1, // faces
447: };
448: static const int hexa_8[8] = {0, 3, 2, 1, 4, 5, 6, 7};
449: static const int hexa_27[27] = {
450: 19, 22, 21, 20, 23, 24, 25, 26, // vertices
451: 10, 9, 8, 7, // bottom edges
452: 16, 15, 18, 17, // mid edges
453: 11, 12, 13, 14, // top edges
454: 1, 3, 5, 4, 6, 2, // faces
455: 0, // center
456: };
457: static const int hexa_64[64] = {
458: // 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
459: 56, 59, 58, 57, 60, 61, 62, 63, // vertices
460: 39, 38, 37, 36, 35, 34, 33, 32, // bottom edges
461: 51, 50, 48, 49, 52, 53, 55, 54, // mid edges; Paraview needs edge 21-22 swapped with 23-24
462: 40, 41, 42, 43, 44, 45, 46, 47, // top edges
463: 8, 10, 11, 9, // z-minus face
464: 16, 17, 19, 18, // y-minus face
465: 24, 25, 27, 26, // x-plus face
466: 20, 21, 23, 22, // y-plus face
467: 30, 28, 29, 31, // x-minus face
468: 12, 13, 15, 14, // z-plus face
469: 0, 1, 3, 2, 4, 5, 7, 6, // center
470: };
472: PetscFunctionBegin;
473: *element_type = CGNS_ENUMV(ElementTypeNull);
474: *perm = NULL;
475: switch (cell_type) {
476: case DM_POLYTOPE_SEGMENT:
477: switch (closure_size) {
478: case 2:
479: *element_type = CGNS_ENUMV(BAR_2);
480: *perm = bar_2;
481: case 3:
482: *element_type = CGNS_ENUMV(BAR_3);
483: *perm = bar_3;
484: case 4:
485: *element_type = CGNS_ENUMV(BAR_4);
486: *perm = bar_4;
487: break;
488: case 5:
489: *element_type = CGNS_ENUMV(BAR_5);
490: *perm = bar_5;
491: break;
492: default:
493: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
494: }
495: break;
496: case DM_POLYTOPE_TRIANGLE:
497: switch (closure_size) {
498: case 3:
499: *element_type = CGNS_ENUMV(TRI_3);
500: *perm = tri_3;
501: break;
502: case 6:
503: *element_type = CGNS_ENUMV(TRI_6);
504: *perm = tri_6;
505: break;
506: case 10:
507: *element_type = CGNS_ENUMV(TRI_10);
508: *perm = tri_10;
509: break;
510: default:
511: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
512: }
513: break;
514: case DM_POLYTOPE_QUADRILATERAL:
515: switch (closure_size) {
516: case 4:
517: *element_type = CGNS_ENUMV(QUAD_4);
518: *perm = quad_4;
519: break;
520: case 9:
521: *element_type = CGNS_ENUMV(QUAD_9);
522: *perm = quad_9;
523: break;
524: case 16:
525: *element_type = CGNS_ENUMV(QUAD_16);
526: *perm = quad_16;
527: break;
528: case 25:
529: *element_type = CGNS_ENUMV(QUAD_25);
530: *perm = quad_25;
531: break;
532: default:
533: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
534: }
535: break;
536: case DM_POLYTOPE_TETRAHEDRON:
537: switch (closure_size) {
538: case 4:
539: *element_type = CGNS_ENUMV(TETRA_4);
540: *perm = tetra_4;
541: break;
542: case 10:
543: *element_type = CGNS_ENUMV(TETRA_10);
544: *perm = tetra_10;
545: break;
546: case 20:
547: *element_type = CGNS_ENUMV(TETRA_20);
548: *perm = tetra_20;
549: break;
550: default:
551: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
552: }
553: break;
554: case DM_POLYTOPE_HEXAHEDRON:
555: switch (closure_size) {
556: case 8:
557: *element_type = CGNS_ENUMV(HEXA_8);
558: *perm = hexa_8;
559: break;
560: case 27:
561: *element_type = CGNS_ENUMV(HEXA_27);
562: *perm = hexa_27;
563: break;
564: case 64:
565: *element_type = CGNS_ENUMV(HEXA_64);
566: *perm = hexa_64;
567: break;
568: default:
569: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
570: }
571: break;
572: default:
573: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
574: }
575: PetscFunctionReturn(PETSC_SUCCESS);
576: }
578: // node_l2g must be freed
579: static PetscErrorCode DMPlexCreateNodeNumbering(DM dm, PetscInt *num_local_nodes, PetscInt *num_global_nodes, PetscInt *nStart, PetscInt *nEnd, const PetscInt **node_l2g)
580: {
581: PetscSection local_section;
582: PetscSF point_sf;
583: PetscInt pStart, pEnd, spStart, spEnd, *points, nleaves, ncomp, *nodes;
584: PetscMPIInt comm_size;
585: const PetscInt *ilocal, field = 0;
587: PetscFunctionBegin;
588: *num_local_nodes = -1;
589: *num_global_nodes = -1;
590: *nStart = -1;
591: *nEnd = -1;
592: *node_l2g = NULL;
594: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &comm_size));
595: PetscCall(DMGetLocalSection(dm, &local_section));
596: PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
597: PetscCall(PetscSectionGetChart(local_section, &spStart, &spEnd));
598: PetscCall(DMGetPointSF(dm, &point_sf));
599: if (comm_size == 1) nleaves = 0;
600: else PetscCall(PetscSFGetGraph(point_sf, NULL, &nleaves, &ilocal, NULL));
601: PetscCall(PetscSectionGetFieldComponents(local_section, field, &ncomp));
603: PetscInt local_node = 0, owned_node = 0, owned_start = 0;
604: for (PetscInt p = spStart, leaf = 0; p < spEnd; p++) {
605: PetscInt dof;
606: PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof));
607: PetscAssert(dof % ncomp == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Field dof %" PetscInt_FMT " must be divisible by components %" PetscInt_FMT, dof, ncomp);
608: local_node += dof / ncomp;
609: if (leaf < nleaves && p == ilocal[leaf]) { // skip points owned by a different process
610: leaf++;
611: } else {
612: owned_node += dof / ncomp;
613: }
614: }
615: PetscCallMPI(MPI_Exscan(&owned_node, &owned_start, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
616: PetscCall(PetscMalloc1(pEnd - pStart, &points));
617: owned_node = 0;
618: for (PetscInt p = spStart, leaf = 0; p < spEnd; p++) {
619: if (leaf < nleaves && p == ilocal[leaf]) { // skip points owned by a different process
620: points[p - pStart] = -1;
621: leaf++;
622: continue;
623: }
624: PetscInt dof, offset;
625: PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof));
626: PetscCall(PetscSectionGetFieldOffset(local_section, p, field, &offset));
627: PetscAssert(offset % ncomp == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Field offset %" PetscInt_FMT " must be divisible by components %" PetscInt_FMT, offset, ncomp);
628: points[p - pStart] = owned_start + owned_node;
629: owned_node += dof / ncomp;
630: }
631: if (comm_size > 1) {
632: PetscCall(PetscSFBcastBegin(point_sf, MPIU_INT, points, points, MPI_REPLACE));
633: PetscCall(PetscSFBcastEnd(point_sf, MPIU_INT, points, points, MPI_REPLACE));
634: }
636: // Set up global indices for each local node
637: PetscCall(PetscMalloc1(local_node, &nodes));
638: for (PetscInt p = spStart; p < spEnd; p++) {
639: PetscInt dof, offset;
640: PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof));
641: PetscCall(PetscSectionGetFieldOffset(local_section, p, field, &offset));
642: for (PetscInt n = 0; n < dof / ncomp; n++) nodes[offset / ncomp + n] = points[p - pStart] + n;
643: }
644: PetscCall(PetscFree(points));
645: *num_local_nodes = local_node;
646: *nStart = owned_start;
647: *nEnd = owned_start + owned_node;
648: PetscCall(MPIU_Allreduce(&owned_node, num_global_nodes, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
649: *node_l2g = nodes;
650: PetscFunctionReturn(PETSC_SUCCESS);
651: }
653: PetscErrorCode DMView_PlexCGNS(DM dm, PetscViewer viewer)
654: {
655: PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;
656: PetscInt topo_dim, coord_dim, num_global_elems;
657: PetscInt cStart, cEnd, num_local_nodes, num_global_nodes, nStart, nEnd;
658: const PetscInt *node_l2g;
659: Vec coord;
660: DM colloc_dm, cdm;
661: PetscMPIInt size;
662: const char *dm_name;
663: int base, zone;
664: cgsize_t isize[3];
666: PetscFunctionBegin;
667: if (!cgv->file_num) {
668: PetscInt time_step;
669: PetscCall(DMGetOutputSequenceNumber(dm, &time_step, NULL));
670: PetscCall(PetscViewerCGNSFileOpen_Internal(viewer, time_step));
671: }
672: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
673: PetscCall(DMGetDimension(dm, &topo_dim));
674: PetscCall(DMGetCoordinateDim(dm, &coord_dim));
675: PetscCall(PetscObjectGetName((PetscObject)dm, &dm_name));
676: PetscCallCGNS(cg_base_write(cgv->file_num, dm_name, topo_dim, coord_dim, &base));
677: PetscCallCGNS(cg_goto(cgv->file_num, base, NULL));
678: PetscCallCGNS(cg_dataclass_write(CGNS_ENUMV(NormalizedByDimensional)));
680: {
681: PetscFE fe, fe_coord;
682: PetscClassId ds_id;
683: PetscDualSpace dual_space, dual_space_coord;
684: PetscInt num_fields, field_order = -1, field_order_coord;
685: PetscBool is_simplex;
686: PetscCall(DMGetNumFields(dm, &num_fields));
687: if (num_fields > 0) {
688: PetscCall(DMGetField(dm, 0, NULL, (PetscObject *)&fe));
689: PetscCall(PetscObjectGetClassId((PetscObject)fe, &ds_id));
690: if (ds_id != PETSCFE_CLASSID) {
691: fe = NULL;
692: if (ds_id == PETSCFV_CLASSID) field_order = -1; // use whatever is present for coords; field will be CellCenter
693: else field_order = 1; // assume vertex-based linear elements
694: }
695: } else fe = NULL;
696: if (fe) {
697: PetscCall(PetscFEGetDualSpace(fe, &dual_space));
698: PetscCall(PetscDualSpaceGetOrder(dual_space, &field_order));
699: }
700: PetscCall(DMGetCoordinateDM(dm, &cdm));
701: PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *)&fe_coord));
702: {
703: PetscClassId id;
704: PetscCall(PetscObjectGetClassId((PetscObject)fe_coord, &id));
705: if (id != PETSCFE_CLASSID) fe_coord = NULL;
706: }
707: if (fe_coord) {
708: PetscCall(PetscFEGetDualSpace(fe_coord, &dual_space_coord));
709: PetscCall(PetscDualSpaceGetOrder(dual_space_coord, &field_order_coord));
710: } else field_order_coord = 1;
711: if (field_order > 0 && field_order != field_order_coord) {
712: PetscInt quadrature_order = field_order;
713: PetscCall(DMClone(dm, &colloc_dm));
714: { // Inform the new colloc_dm that it is a coordinate DM so isoperiodic affine corrections can be applied
715: const PetscSF *face_sfs;
716: PetscInt num_face_sfs;
717: PetscCall(DMPlexGetIsoperiodicFaceSF(dm, &num_face_sfs, &face_sfs));
718: PetscCall(DMPlexSetIsoperiodicFaceSF(colloc_dm, num_face_sfs, (PetscSF *)face_sfs));
719: if (face_sfs) colloc_dm->periodic.setup = DMPeriodicCoordinateSetUp_Internal;
720: }
721: PetscCall(DMPlexIsSimplex(dm, &is_simplex));
722: PetscCall(PetscFECreateLagrange(PetscObjectComm((PetscObject)dm), topo_dim, coord_dim, is_simplex, field_order, quadrature_order, &fe));
723: PetscCall(DMSetCoordinateDisc(colloc_dm, fe, PETSC_TRUE));
724: PetscCall(PetscFEDestroy(&fe));
725: } else {
726: PetscCall(PetscObjectReference((PetscObject)dm));
727: colloc_dm = dm;
728: }
729: }
730: PetscCall(DMGetCoordinateDM(colloc_dm, &cdm));
731: PetscCall(DMPlexCreateNodeNumbering(cdm, &num_local_nodes, &num_global_nodes, &nStart, &nEnd, &node_l2g));
732: PetscCall(DMGetCoordinatesLocal(colloc_dm, &coord));
733: PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
734: num_global_elems = cEnd - cStart;
735: PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &num_global_elems, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
736: isize[0] = num_global_nodes;
737: isize[1] = num_global_elems;
738: isize[2] = 0;
739: PetscCallCGNS(cg_zone_write(cgv->file_num, base, "Zone", isize, CGNS_ENUMV(Unstructured), &zone));
741: {
742: const PetscScalar *X;
743: PetscScalar *x;
744: int coord_ids[3];
746: PetscCall(VecGetArrayRead(coord, &X));
747: for (PetscInt d = 0; d < coord_dim; d++) {
748: const double exponents[] = {0, 1, 0, 0, 0};
749: char coord_name[64];
750: PetscCall(PetscSNPrintf(coord_name, sizeof coord_name, "Coordinate%c", 'X' + (int)d));
751: PetscCallCGNS(cgp_coord_write(cgv->file_num, base, zone, CGNS_ENUMV(RealDouble), coord_name, &coord_ids[d]));
752: PetscCallCGNS(cg_goto(cgv->file_num, base, "Zone_t", zone, "GridCoordinates", 0, coord_name, 0, NULL));
753: PetscCallCGNS(cg_exponents_write(CGNS_ENUMV(RealDouble), exponents));
754: }
756: DMPolytopeType cell_type;
757: int section;
758: cgsize_t e_owned, e_global, e_start, *conn = NULL;
759: const int *perm;
760: CGNS_ENUMT(ElementType_t) element_type = CGNS_ENUMV(ElementTypeNull);
761: {
762: PetscCall(PetscMalloc1(nEnd - nStart, &x));
763: for (PetscInt d = 0; d < coord_dim; d++) {
764: for (PetscInt n = 0; n < num_local_nodes; n++) {
765: PetscInt gn = node_l2g[n];
766: if (gn < nStart || nEnd <= gn) continue;
767: x[gn - nStart] = X[n * coord_dim + d];
768: }
769: // CGNS nodes use 1-based indexing
770: cgsize_t start = nStart + 1, end = nEnd;
771: PetscCallCGNS(cgp_coord_write_data(cgv->file_num, base, zone, coord_ids[d], &start, &end, x));
772: }
773: PetscCall(PetscFree(x));
774: PetscCall(VecRestoreArrayRead(coord, &X));
775: }
777: PetscCall(DMPlexGetCellType(dm, cStart, &cell_type));
778: for (PetscInt i = cStart, c = 0; i < cEnd; i++) {
779: PetscInt closure_dof, *closure_indices, elem_size;
780: PetscCall(DMPlexGetClosureIndices(cdm, cdm->localSection, cdm->localSection, i, PETSC_FALSE, &closure_dof, &closure_indices, NULL, NULL));
781: elem_size = closure_dof / coord_dim;
782: if (!conn) PetscCall(PetscMalloc1((cEnd - cStart) * elem_size, &conn));
783: PetscCall(DMPlexCGNSGetPermutation_Internal(cell_type, closure_dof / coord_dim, &element_type, &perm));
784: for (PetscInt j = 0; j < elem_size; j++) conn[c++] = node_l2g[closure_indices[perm[j] * coord_dim] / coord_dim] + 1;
785: PetscCall(DMPlexRestoreClosureIndices(cdm, cdm->localSection, cdm->localSection, i, PETSC_FALSE, &closure_dof, &closure_indices, NULL, NULL));
786: }
787: e_owned = cEnd - cStart;
788: PetscCall(MPIU_Allreduce(&e_owned, &e_global, 1, MPIU_CGSIZE, MPI_SUM, PetscObjectComm((PetscObject)dm)));
789: 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);
790: e_start = 0;
791: PetscCallMPI(MPI_Exscan(&e_owned, &e_start, 1, MPIU_CGSIZE, MPI_SUM, PetscObjectComm((PetscObject)dm)));
792: PetscCallCGNS(cgp_section_write(cgv->file_num, base, zone, "Elem", element_type, 1, e_global, 0, §ion));
793: PetscCallCGNS(cgp_elements_write_data(cgv->file_num, base, zone, section, e_start + 1, e_start + e_owned, conn));
794: PetscCall(PetscFree(conn));
796: cgv->base = base;
797: cgv->zone = zone;
798: cgv->node_l2g = node_l2g;
799: cgv->num_local_nodes = num_local_nodes;
800: cgv->nStart = nStart;
801: cgv->nEnd = nEnd;
802: cgv->eStart = e_start;
803: cgv->eEnd = e_start + e_owned;
804: if (1) {
805: PetscMPIInt rank;
806: int *efield;
807: int sol, field;
808: DMLabel label;
809: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
810: PetscCall(PetscMalloc1(e_owned, &efield));
811: for (PetscInt i = 0; i < e_owned; i++) efield[i] = rank;
812: PetscCallCGNS(cg_sol_write(cgv->file_num, base, zone, "CellInfo", CGNS_ENUMV(CellCenter), &sol));
813: PetscCallCGNS(cgp_field_write(cgv->file_num, base, zone, sol, CGNS_ENUMV(Integer), "Rank", &field));
814: cgsize_t start = e_start + 1, end = e_start + e_owned;
815: PetscCallCGNS(cgp_field_write_data(cgv->file_num, base, zone, sol, field, &start, &end, efield));
816: PetscCall(DMGetLabel(dm, "Cell Sets", &label));
817: if (label) {
818: for (PetscInt c = cStart; c < cEnd; c++) {
819: PetscInt value;
820: PetscCall(DMLabelGetValue(label, c, &value));
821: efield[c - cStart] = value;
822: }
823: PetscCallCGNS(cgp_field_write(cgv->file_num, base, zone, sol, CGNS_ENUMV(Integer), "CellSet", &field));
824: PetscCallCGNS(cgp_field_write_data(cgv->file_num, base, zone, sol, field, &start, &end, efield));
825: }
826: PetscCall(PetscFree(efield));
827: }
828: }
829: PetscCall(DMDestroy(&colloc_dm));
830: PetscFunctionReturn(PETSC_SUCCESS);
831: }
833: static PetscErrorCode PetscCGNSDataType(PetscDataType pd, CGNS_ENUMT(DataType_t) * cd)
834: {
835: PetscFunctionBegin;
836: switch (pd) {
837: case PETSC_FLOAT:
838: *cd = CGNS_ENUMV(RealSingle);
839: break;
840: case PETSC_DOUBLE:
841: *cd = CGNS_ENUMV(RealDouble);
842: break;
843: case PETSC_COMPLEX:
844: *cd = CGNS_ENUMV(ComplexDouble);
845: break;
846: default:
847: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Data type %s", PetscDataTypes[pd]);
848: }
849: PetscFunctionReturn(PETSC_SUCCESS);
850: }
852: PetscErrorCode VecView_Plex_Local_CGNS(Vec V, PetscViewer viewer)
853: {
854: PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;
855: DM dm;
856: PetscSection section;
857: PetscInt time_step, num_fields, pStart, pEnd, cStart, cEnd;
858: PetscReal time, *time_slot;
859: size_t *step_slot;
860: const PetscScalar *v;
861: char solution_name[PETSC_MAX_PATH_LEN];
862: int sol;
864: PetscFunctionBegin;
865: PetscCall(VecGetDM(V, &dm));
866: if (!cgv->node_l2g) PetscCall(DMView(dm, viewer));
867: if (!cgv->nodal_field) PetscCall(PetscMalloc1(PetscMax(cgv->nEnd - cgv->nStart, cgv->eEnd - cgv->eStart), &cgv->nodal_field));
868: if (!cgv->output_times) PetscCall(PetscSegBufferCreate(sizeof(PetscReal), 20, &cgv->output_times));
869: if (!cgv->output_steps) PetscCall(PetscSegBufferCreate(sizeof(size_t), 20, &cgv->output_steps));
871: PetscCall(DMGetOutputSequenceNumber(dm, &time_step, &time));
872: if (time_step < 0) {
873: time_step = 0;
874: time = 0.;
875: }
876: PetscCall(PetscSegBufferGet(cgv->output_times, 1, &time_slot));
877: *time_slot = time;
878: PetscCall(PetscSegBufferGet(cgv->output_steps, 1, &step_slot));
879: *step_slot = time_step;
880: PetscCall(PetscSNPrintf(solution_name, sizeof solution_name, "FlowSolution%" PetscInt_FMT, time_step));
881: PetscCall(DMGetLocalSection(dm, §ion));
882: PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
883: PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
884: CGNS_ENUMT(GridLocation_t) grid_loc = CGNS_ENUMV(Vertex);
885: if (cStart == pStart && cEnd == pEnd) grid_loc = CGNS_ENUMV(CellCenter);
886: PetscCallCGNS(cg_sol_write(cgv->file_num, cgv->base, cgv->zone, solution_name, grid_loc, &sol));
887: PetscCall(VecGetArrayRead(V, &v));
888: PetscCall(PetscSectionGetNumFields(section, &num_fields));
889: for (PetscInt field = 0; field < num_fields; field++) {
890: PetscInt ncomp;
891: const char *field_name;
892: PetscCall(PetscSectionGetFieldName(section, field, &field_name));
893: PetscCall(PetscSectionGetFieldComponents(section, field, &ncomp));
894: for (PetscInt comp = 0; comp < ncomp; comp++) {
895: int cgfield;
896: const char *comp_name;
897: char cgns_field_name[32]; // CGNS max field name is 32
898: CGNS_ENUMT(DataType_t) datatype;
899: PetscCall(PetscSectionGetComponentName(section, field, comp, &comp_name));
900: if (ncomp == 1 && comp_name[0] == '0' && comp_name[1] == '\0') PetscCall(PetscStrncpy(cgns_field_name, field_name, sizeof cgns_field_name));
901: else if (field_name[0] == '\0') PetscCall(PetscStrncpy(cgns_field_name, comp_name, sizeof cgns_field_name));
902: else PetscCall(PetscSNPrintf(cgns_field_name, sizeof cgns_field_name, "%s.%s", field_name, comp_name));
903: PetscCall(PetscCGNSDataType(PETSC_SCALAR, &datatype));
904: PetscCallCGNS(cgp_field_write(cgv->file_num, cgv->base, cgv->zone, sol, datatype, cgns_field_name, &cgfield));
905: for (PetscInt p = pStart, n = 0; p < pEnd; p++) {
906: PetscInt off, dof;
907: PetscCall(PetscSectionGetFieldDof(section, p, field, &dof));
908: if (dof == 0) continue;
909: PetscCall(PetscSectionGetFieldOffset(section, p, field, &off));
910: for (PetscInt c = comp; c < dof; c += ncomp, n++) {
911: switch (grid_loc) {
912: case CGNS_ENUMV(Vertex): {
913: PetscInt gn = cgv->node_l2g[n];
914: if (gn < cgv->nStart || cgv->nEnd <= gn) continue;
915: cgv->nodal_field[gn - cgv->nStart] = v[off + c];
916: } break;
917: case CGNS_ENUMV(CellCenter): {
918: cgv->nodal_field[n] = v[off + c];
919: } break;
920: default:
921: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Can only pack for Vertex and CellCenter grid locations");
922: }
923: }
924: }
925: // CGNS nodes use 1-based indexing
926: cgsize_t start = cgv->nStart + 1, end = cgv->nEnd;
927: if (grid_loc == CGNS_ENUMV(CellCenter)) {
928: start = cgv->eStart + 1;
929: end = cgv->eEnd;
930: }
931: PetscCallCGNS(cgp_field_write_data(cgv->file_num, cgv->base, cgv->zone, sol, cgfield, &start, &end, cgv->nodal_field));
932: }
933: }
934: PetscCall(VecRestoreArrayRead(V, &v));
935: PetscCall(PetscViewerCGNSCheckBatch_Internal(viewer));
936: PetscFunctionReturn(PETSC_SUCCESS);
937: }