1 #define PETSCDM_DLL 2 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "DMPlexCreateFluentFromFile" 6 /*@C 7 DMPlexCreateFluentFromFile - Create a DMPlex mesh from a Fluent mesh file 8 9 + comm - The MPI communicator 10 . filename - Name of the Fluent mesh file 11 - interpolate - Create faces and edges in the mesh 12 13 Output Parameter: 14 . dm - The DM object representing the mesh 15 16 Level: beginner 17 18 .seealso: DMPlexCreateFromFile(), DMPlexCreateFluent(), DMPlexCreate() 19 @*/ 20 PetscErrorCode DMPlexCreateFluentFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 21 { 22 PetscViewer viewer; 23 PetscErrorCode ierr; 24 25 PetscFunctionBegin; 26 /* Create file viewer and build plex */ 27 ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr); 28 ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr); 29 ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr); 30 ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr); 31 ierr = DMPlexCreateFluent(comm, viewer, interpolate, dm);CHKERRQ(ierr); 32 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 33 PetscFunctionReturn(0); 34 } 35 36 #undef __FUNCT__ 37 #define __FUNCT__ "DMPlexCreateFluent_ReadString" 38 PetscErrorCode DMPlexCreateFluent_ReadString(PetscViewer viewer, char *buffer, char delim) 39 { 40 PetscInt ret, i = 0; 41 PetscErrorCode ierr; 42 43 PetscFunctionBegin; 44 do {ierr = PetscViewerRead(viewer, &(buffer[i++]), 1, &ret, PETSC_CHAR);CHKERRQ(ierr);} 45 while (ret > 0 && buffer[i-1] != '\0' && buffer[i-1] != delim); 46 buffer[i] = '\0'; 47 PetscFunctionReturn(0); 48 } 49 50 #undef __FUNCT__ 51 #define __FUNCT__ "DMPlexCreateFluent_ReadValues" 52 PetscErrorCode DMPlexCreateFluent_ReadValues(PetscViewer viewer, void *data, PetscInt count, PetscDataType dtype, PetscBool binary) 53 { 54 int fdes=0; 55 FILE *file; 56 PetscInt i; 57 PetscErrorCode ierr; 58 59 PetscFunctionBegin; 60 if (binary) { 61 /* Extract raw file descriptor to read binary block */ 62 ierr = PetscViewerASCIIGetPointer(viewer, &file);CHKERRQ(ierr); 63 fflush(file); fdes = fileno(file); 64 } 65 66 if (!binary && dtype == PETSC_INT) { 67 char cbuf[256]; 68 int ibuf, snum; 69 /* Parse hexadecimal ascii integers */ 70 for (i = 0; i < count; i++) { 71 ierr = PetscViewerRead(viewer, cbuf, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 72 snum = sscanf(cbuf, "%x", &ibuf); 73 if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 74 ((PetscInt*)data)[i] = (PetscInt)ibuf; 75 } 76 } else if (binary && dtype == PETSC_INT) { 77 /* Always read 32-bit ints and cast to PetscInt */ 78 int *ibuf; 79 ierr = PetscMalloc1(count, &ibuf);CHKERRQ(ierr); 80 ierr = PetscBinaryRead(fdes, ibuf, count, PETSC_ENUM);CHKERRQ(ierr); 81 ierr = PetscByteSwap(ibuf, PETSC_ENUM, count);CHKERRQ(ierr); 82 for (i = 0; i < count; i++) ((PetscInt*)data)[i] = (PetscInt)(ibuf[i]); 83 ierr = PetscFree(ibuf);CHKERRQ(ierr); 84 85 } else if (binary && dtype == PETSC_SCALAR) { 86 float *fbuf; 87 /* Always read 32-bit floats and cast to PetscScalar */ 88 ierr = PetscMalloc1(count, &fbuf);CHKERRQ(ierr); 89 ierr = PetscBinaryRead(fdes, fbuf, count, PETSC_FLOAT);CHKERRQ(ierr); 90 ierr = PetscByteSwap(fbuf, PETSC_FLOAT, count);CHKERRQ(ierr); 91 for (i = 0; i < count; i++) ((PetscScalar*)data)[i] = (PetscScalar)(fbuf[i]); 92 ierr = PetscFree(fbuf);CHKERRQ(ierr); 93 } else { 94 ierr = PetscViewerASCIIRead(viewer, data, count, NULL, dtype);CHKERRQ(ierr); 95 } 96 PetscFunctionReturn(0); 97 } 98 99 #undef __FUNCT__ 100 #define __FUNCT__ "DMPlexCreateFluent_ReadSection" 101 PetscErrorCode DMPlexCreateFluent_ReadSection(PetscViewer viewer, FluentSection *s) 102 { 103 char buffer[PETSC_MAX_PATH_LEN]; 104 int snum; 105 PetscErrorCode ierr; 106 107 PetscFunctionBegin; 108 /* Fast-forward to next section and derive its index */ 109 ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr); 110 ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ' ');CHKERRQ(ierr); 111 snum = sscanf(buffer, "%d", &(s->index)); 112 /* If we can't match an index return -1 to signal end-of-file */ 113 if (snum < 1) {s->index = -1; PetscFunctionReturn(0);} 114 115 if (s->index == 0) { /* Comment */ 116 ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 117 118 } else if (s->index == 2) { /* Dimension */ 119 ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 120 snum = sscanf(buffer, "%d", &(s->nd)); 121 if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 122 123 } else if (s->index == 10 || s->index == 2010) { /* Vertices */ 124 ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 125 snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd)); 126 if (snum != 5) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 127 if (s->zoneID > 0) { 128 PetscInt numCoords = s->last - s->first + 1; 129 ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr); 130 ierr = PetscMalloc1(s->nd*numCoords, (PetscScalar**)&s->data);CHKERRQ(ierr); 131 ierr = DMPlexCreateFluent_ReadValues(viewer, s->data, s->nd*numCoords, PETSC_SCALAR, s->index==2010 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 132 ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 133 } 134 ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 135 136 } else if (s->index == 12 || s->index == 2012) { /* Cells */ 137 ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 138 snum = sscanf(buffer, "(%x", &(s->zoneID)); 139 if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 140 if (s->zoneID == 0) { /* Header section */ 141 snum = sscanf(buffer, "(%x %x %x %d)", &(s->zoneID), &(s->first), &(s->last), &(s->nd)); 142 if (snum != 4) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 143 } else { /* Data section */ 144 snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd)); 145 if (snum != 5) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 146 if (s->nd == 0) { 147 /* Read cell type definitions for mixed cells */ 148 PetscInt numCells = s->last - s->first + 1; 149 ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr); 150 ierr = PetscMalloc1(numCells, (PetscInt**)&s->data);CHKERRQ(ierr); 151 ierr = DMPlexCreateFluent_ReadValues(viewer, s->data, numCells, PETSC_INT, s->index==2012 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 152 ierr = PetscFree(s->data);CHKERRQ(ierr); 153 ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 154 } 155 } 156 ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 157 158 } else if (s->index == 13 || s->index == 2013) { /* Faces */ 159 ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 160 snum = sscanf(buffer, "(%x", &(s->zoneID)); 161 if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 162 if (s->zoneID == 0) { /* Header section */ 163 snum = sscanf(buffer, "(%x %x %x %d)", &(s->zoneID), &(s->first), &(s->last), &(s->nd)); 164 if (snum != 4) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 165 } else { /* Data section */ 166 PetscInt f, numEntries, numFaces, numFaceVert; 167 snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd)); 168 if (snum != 5) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 169 ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr); 170 switch (s->nd) { 171 case 0: numEntries = PETSC_DETERMINE; break; 172 case 2: numEntries = 2 + 2; break; /* linear */ 173 case 3: numEntries = 2 + 3; break; /* triangular */ 174 case 4: numEntries = 2 + 4; break; /* quadrilateral */ 175 default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown face type in Fluent file"); 176 } 177 numFaces = s->last-s->first + 1; 178 if (numEntries != PETSC_DETERMINE) { 179 /* Allocate space only if we already know the size of the block */ 180 ierr = PetscMalloc1(numEntries*numFaces, (PetscInt**)&s->data);CHKERRQ(ierr); 181 } 182 for (f = 0; f < numFaces; f++) { 183 if (s->nd == 0) { 184 /* Determine the size of the block for "mixed" facets */ 185 ierr = DMPlexCreateFluent_ReadValues(viewer, &numFaceVert, 1, PETSC_INT, s->index==2013 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 186 if (numEntries == PETSC_DETERMINE) { 187 numEntries = numFaceVert + 2; 188 ierr = PetscMalloc1(numEntries*numFaces, (PetscInt**)&s->data);CHKERRQ(ierr); 189 } else { 190 if (numEntries != numFaceVert + 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No support for mixed faces in Fluent files"); 191 } 192 } 193 ierr = DMPlexCreateFluent_ReadValues(viewer, &(((PetscInt*)s->data)[f*numEntries]), numEntries, PETSC_INT, s->index==2013 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 194 } 195 s->nd = numEntries - 2; 196 ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 197 } 198 ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 199 200 } else { /* Unknown section type */ 201 PetscInt depth = 1; 202 do { 203 /* Match parentheses when parsing unknown sections */ 204 do {ierr = PetscViewerRead(viewer, &(buffer[0]), 1, NULL, PETSC_CHAR);CHKERRQ(ierr);} 205 while (buffer[0] != '(' && buffer[0] != ')'); 206 if (buffer[0] == '(') depth++; 207 if (buffer[0] == ')') depth--; 208 } while (depth > 0); 209 ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '\n');CHKERRQ(ierr); 210 } 211 PetscFunctionReturn(0); 212 } 213 214 #undef __FUNCT__ 215 #define __FUNCT__ "DMPlexCreateFluent" 216 /*@C 217 DMPlexCreateFluent - Create a DMPlex mesh from a Fluent mesh file. 218 219 Collective on comm 220 221 Input Parameters: 222 + comm - The MPI communicator 223 . viewer - The Viewer associated with a Fluent mesh file 224 - interpolate - Create faces and edges in the mesh 225 226 Output Parameter: 227 . dm - The DM object representing the mesh 228 229 Note: http://aerojet.engr.ucdavis.edu/fluenthelp/html/ug/node1490.htm 230 231 Level: beginner 232 233 .keywords: mesh, fluent, case 234 .seealso: DMPLEX, DMCreate() 235 @*/ 236 PetscErrorCode DMPlexCreateFluent(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 237 { 238 PetscMPIInt rank; 239 PetscInt c, f, v, dim = PETSC_DETERMINE, numCells = 0, numVertices = 0, numCellVertices = PETSC_DETERMINE; 240 PetscInt numFaces = PETSC_DETERMINE, numFaceEntries = PETSC_DETERMINE, numFaceVertices = PETSC_DETERMINE; 241 PetscInt *faces = NULL, *cellVertices, *faceZoneIDs = NULL; 242 PetscInt d, coordSize; 243 PetscScalar *coords, *coordsIn = NULL; 244 PetscSection coordSection; 245 Vec coordinates; 246 PetscErrorCode ierr; 247 248 PetscFunctionBegin; 249 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 250 251 if (!rank) { 252 FluentSection s; 253 do { 254 ierr = DMPlexCreateFluent_ReadSection(viewer, &s);CHKERRQ(ierr); 255 if (s.index == 2) { /* Dimension */ 256 dim = s.nd; 257 258 } else if (s.index == 10 || s.index == 2010) { /* Vertices */ 259 if (s.zoneID == 0) numVertices = s.last; 260 else { 261 if (coordsIn) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Currently no support for multiple coordinate sets in Fluent files"); 262 coordsIn = (PetscScalar *) s.data; 263 } 264 265 } else if (s.index == 12 || s.index == 2012) { /* Cells */ 266 if (s.zoneID == 0) numCells = s.last; 267 else { 268 switch (s.nd) { 269 case 0: numCellVertices = PETSC_DETERMINE; break; 270 case 1: numCellVertices = 3; break; /* triangular */ 271 case 2: numCellVertices = 4; break; /* tetrahedral */ 272 case 3: numCellVertices = 4; break; /* quadrilateral */ 273 case 4: numCellVertices = 8; break; /* hexahedral */ 274 case 5: numCellVertices = 5; break; /* pyramid */ 275 case 6: numCellVertices = 6; break; /* wedge */ 276 default: numCellVertices = PETSC_DETERMINE; 277 } 278 } 279 280 } else if (s.index == 13 || s.index == 2013) { /* Facets */ 281 if (s.zoneID == 0) { /* Header section */ 282 numFaces = s.last - s.first + 1; 283 if (s.nd == 0 || s.nd == 5) numFaceVertices = PETSC_DETERMINE; 284 else numFaceVertices = s.nd; 285 } else { /* Data section */ 286 if (numFaceVertices != PETSC_DETERMINE && s.nd != numFaceVertices) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mixed facets in Fluent files are not supported"); 287 if (numFaces < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No header section for facets in Fluent file"); 288 if (numFaceVertices == PETSC_DETERMINE) numFaceVertices = s.nd; 289 numFaceEntries = numFaceVertices + 2; 290 if (!faces) {ierr = PetscMalloc1(numFaces*numFaceEntries, &faces);CHKERRQ(ierr);} 291 if (!faceZoneIDs) {ierr = PetscMalloc1(numFaces, &faceZoneIDs);CHKERRQ(ierr);} 292 ierr = PetscMemcpy(&(faces[(s.first-1)*numFaceEntries]), s.data, (s.last-s.first+1)*numFaceEntries*sizeof(PetscInt));CHKERRQ(ierr); 293 /* Record the zoneID for each face set */ 294 for (f = s.first -1; f < s.last; f++) faceZoneIDs[f] = s.zoneID; 295 ierr = PetscFree(s.data);CHKERRQ(ierr); 296 } 297 } 298 } while (s.index >= 0); 299 } 300 ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 301 if (dim < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Fluent file does not include dimension"); 302 303 /* Allocate cell-vertex mesh */ 304 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 305 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 306 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 307 ierr = DMPlexSetChart(*dm, 0, numCells + numVertices);CHKERRQ(ierr); 308 if (!rank) { 309 if (numCells < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown number of cells in Fluent file"); 310 /* If no cell type was given we assume simplices */ 311 if (numCellVertices == PETSC_DETERMINE) numCellVertices = numFaceVertices + 1; 312 for (c = 0; c < numCells; ++c) {ierr = DMPlexSetConeSize(*dm, c, numCellVertices);CHKERRQ(ierr);} 313 } 314 ierr = DMSetUp(*dm);CHKERRQ(ierr); 315 316 if (!rank) { 317 /* Derive cell-vertex list from face-vertex and face-cell maps */ 318 ierr = PetscMalloc1(numCells*numCellVertices, &cellVertices);CHKERRQ(ierr); 319 for (c = 0; c < numCells*numCellVertices; c++) cellVertices[c] = -1; 320 for (f = 0; f < numFaces; f++) { 321 PetscInt *cell; 322 const PetscInt cl = faces[f*numFaceEntries + numFaceVertices]; 323 const PetscInt cr = faces[f*numFaceEntries + numFaceVertices + 1]; 324 const PetscInt *face = &(faces[f*numFaceEntries]); 325 326 if (cl > 0) { 327 cell = &(cellVertices[(cl-1) * numCellVertices]); 328 for (v = 0; v < numFaceVertices; v++) { 329 PetscBool found = PETSC_FALSE; 330 for (c = 0; c < numCellVertices; c++) { 331 if (cell[c] < 0) break; 332 if (cell[c] == face[v]-1 + numCells) {found = PETSC_TRUE; break;} 333 } 334 if (!found) cell[c] = face[v]-1 + numCells; 335 } 336 } 337 if (cr > 0) { 338 cell = &(cellVertices[(cr-1) * numCellVertices]); 339 for (v = 0; v < numFaceVertices; v++) { 340 PetscBool found = PETSC_FALSE; 341 for (c = 0; c < numCellVertices; c++) { 342 if (cell[c] < 0) break; 343 if (cell[c] == face[v]-1 + numCells) {found = PETSC_TRUE; break;} 344 } 345 if (!found) cell[c] = face[v]-1 + numCells; 346 } 347 } 348 } 349 for (c = 0; c < numCells; c++) { 350 ierr = DMPlexSetCone(*dm, c, &(cellVertices[c*numCellVertices]));CHKERRQ(ierr); 351 } 352 } 353 ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 354 ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 355 if (interpolate) { 356 DM idm = NULL; 357 358 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 359 ierr = DMDestroy(dm);CHKERRQ(ierr); 360 *dm = idm; 361 } 362 363 if (!rank) { 364 PetscInt fi, joinSize, meetSize, *fverts, cells[2]; 365 const PetscInt *join, *meet; 366 ierr = PetscMalloc1(numFaceVertices, &fverts);CHKERRQ(ierr); 367 /* Mark facets by finding the full join of all adjacent vertices */ 368 for (f = 0; f < numFaces; f++) { 369 const PetscInt cl = faces[f*numFaceEntries + numFaceVertices] - 1; 370 const PetscInt cr = faces[f*numFaceEntries + numFaceVertices + 1] - 1; 371 if (cl > 0 && cr > 0) { 372 /* If we know both adjoining cells we can use a single-level meet */ 373 cells[0] = cl; cells[1] = cr; 374 ierr = DMPlexGetMeet(*dm, 2, cells, &meetSize, &meet);CHKERRQ(ierr); 375 if (meetSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for Fluent face %d", f); 376 ierr = DMSetLabelValue(*dm, "Face Sets", meet[0], faceZoneIDs[f]);CHKERRQ(ierr); 377 ierr = DMPlexRestoreMeet(*dm, numFaceVertices, fverts, &meetSize, &meet);CHKERRQ(ierr); 378 } else { 379 for (fi = 0; fi < numFaceVertices; fi++) fverts[fi] = faces[f*numFaceEntries + fi] + numCells - 1; 380 ierr = DMPlexGetFullJoin(*dm, numFaceVertices, fverts, &joinSize, &join);CHKERRQ(ierr); 381 if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for Fluent face %d", f); 382 ierr = DMSetLabelValue(*dm, "Face Sets", join[0], faceZoneIDs[f]);CHKERRQ(ierr); 383 ierr = DMPlexRestoreJoin(*dm, numFaceVertices, fverts, &joinSize, &join);CHKERRQ(ierr); 384 } 385 } 386 ierr = PetscFree(fverts);CHKERRQ(ierr); 387 } 388 389 /* Read coordinates */ 390 ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 391 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 392 ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 393 ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 394 for (v = numCells; v < numCells+numVertices; ++v) { 395 ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 396 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 397 } 398 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 399 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 400 ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr); 401 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 402 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 403 ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 404 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 405 if (!rank) { 406 for (v = 0; v < numVertices; ++v) { 407 for (d = 0; d < dim; ++d) { 408 coords[v*dim+d] = coordsIn[v*dim+d]; 409 } 410 } 411 } 412 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 413 ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 414 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 415 if (!rank) { 416 ierr = PetscFree(cellVertices);CHKERRQ(ierr); 417 ierr = PetscFree(faces);CHKERRQ(ierr); 418 ierr = PetscFree(faceZoneIDs);CHKERRQ(ierr); 419 ierr = PetscFree(coordsIn);CHKERRQ(ierr); 420 } 421 PetscFunctionReturn(0); 422 } 423