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