Actual source code: plexmed.c

petsc-3.14.6 2021-03-30
Report Typos and Errors
  1: #define PETSCDM_DLL
  2: #include <petsc/private/dmpleximpl.h>

  4: #if defined(PETSC_HAVE_MED)
  5: #include <med.h>
  6: #endif

  8: /*@C
  9:   DMPlexCreateMedFromFile - Create a DMPlex mesh from a (Salome-)Med file.

 11: + comm        - The MPI communicator
 12: . filename    - Name of the .med file
 13: - interpolate - Create faces and edges in the mesh

 15:   Output Parameter:
 16: . dm  - The DM object representing the mesh

 18:   Note: https://www.salome-platform.org/user-section/about/med, http://docs.salome-platform.org/latest/MED_index.html

 20:   Level: beginner

 22: .seealso: DMPlexCreateFromFile(), DMPlexCreateGmsh(), DMPlexCreate()
 23: @*/
 24: PetscErrorCode DMPlexCreateMedFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
 25: {
 26:   PetscMPIInt     rank, size;
 27: #if defined(PETSC_HAVE_MED)
 28:   med_idt         fileID;
 29:   PetscInt        i, ngeo, spaceDim, meshDim;
 30:   PetscInt        numVertices = 0, numCells = 0, c, numCorners, numCellsLocal, numVerticesLocal;
 31:   med_int        *medCellList;
 32:   PetscInt       *cellList;
 33:   med_float      *coordinates = NULL;
 34:   PetscReal      *vertcoords = NULL;
 35:   PetscLayout     vLayout, cLayout;
 36:   const PetscInt *vrange, *crange;
 37:   PetscSF         sfVertices;
 38:   char           *axisname, *unitname, meshname[MED_NAME_SIZE+1], geotypename[MED_NAME_SIZE+1];
 39:   char            meshdescription[MED_COMMENT_SIZE+1], dtunit[MED_SNAME_SIZE+1];
 40:   med_sorting_type sortingtype;
 41:   med_mesh_type   meshtype;
 42:   med_axis_type   axistype;
 43:   med_bool        coordinatechangement, geotransformation, hdfok, medok;
 44:   med_geometry_type geotype[2], cellID, facetID;
 45:   med_filter      vfilter = MED_FILTER_INIT;
 46:   med_filter      cfilter = MED_FILTER_INIT;
 47:   med_err         mederr;
 48:   med_int         major, minor, release;
 49: #endif
 50:   PetscErrorCode  ierr;

 53:   MPI_Comm_rank(comm, &rank);
 54:   MPI_Comm_size(comm, &size);
 55: #if defined(PETSC_HAVE_MED)
 56:   mederr = MEDfileCompatibility(filename, &hdfok, &medok);
 57:   if (mederr) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Cannot determine MED file compatibility: %s", filename);
 58:   if (!hdfok) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Not a compatible HDF format: %s", filename);
 59:   if (!medok) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Not a compatible MED format: %s", filename);

 61:   fileID = MEDfileOpen(filename, MED_ACC_RDONLY);
 62:   if (fileID < 0) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Unable to open .med mesh file: %s", filename);
 63:   mederr = MEDfileNumVersionRd(fileID, &major, &minor, &release);
 64:   if (MEDnMesh (fileID) < 1) SETERRQ4(comm, PETSC_ERR_ARG_WRONG, "No meshes found in .med v%d.%d.%d mesh file: %s", major, minor, release, filename);
 65:   spaceDim = MEDmeshnAxis(fileID, 1);
 66:   if (spaceDim < 1) SETERRQ4(comm, PETSC_ERR_ARG_WRONG, "Mesh of unknown space dimension found in .med v%d.%d.%d mesh file: %s", major, minor, release, filename);
 67:   /* Read general mesh information */
 68:   PetscMalloc1(MED_SNAME_SIZE*spaceDim+1, &axisname);
 69:   PetscMalloc1(MED_SNAME_SIZE*spaceDim+1, &unitname);
 70:   {
 71:     med_int medMeshDim, medNstep;
 72:     med_int medSpaceDim = spaceDim;

 74:     MEDmeshInfo(fileID, 1, meshname, &medSpaceDim, &medMeshDim, &meshtype, meshdescription,
 75:                        dtunit, &sortingtype, &medNstep, &axistype, axisname, unitname);
 76:     spaceDim = medSpaceDim;
 77:     meshDim  = medMeshDim;
 78:   }
 79:   PetscFree(axisname);
 80:   PetscFree(unitname);
 81:   /* Partition mesh coordinates */
 82:   numVertices = MEDmeshnEntity(fileID, meshname, MED_NO_DT, MED_NO_IT, MED_NODE, MED_NO_GEOTYPE,
 83:                                MED_COORDINATE, MED_NO_CMODE,&coordinatechangement, &geotransformation);
 84:   PetscLayoutCreate(comm, &vLayout);
 85:   PetscLayoutSetSize(vLayout, numVertices);
 86:   PetscLayoutSetBlockSize(vLayout, 1);
 87:   PetscLayoutSetUp(vLayout);
 88:   PetscLayoutGetRanges(vLayout, &vrange);
 89:   numVerticesLocal = vrange[rank+1]-vrange[rank];
 90:   MEDfilterBlockOfEntityCr(fileID, numVertices, 1, spaceDim, MED_ALL_CONSTITUENT, MED_FULL_INTERLACE, MED_COMPACT_STMODE,
 91:                                   MED_NO_PROFILE, vrange[rank]+1, 1, numVerticesLocal, 1, 1, &vfilter);
 92:   /* Read mesh coordinates */
 93:   if (numVertices < 0) SETERRQ4(comm, PETSC_ERR_ARG_WRONG, "No nodes found in .med v%d.%d.%d mesh file: %s", major, minor, release, filename);
 94:   PetscMalloc1(numVerticesLocal*spaceDim, &coordinates);
 95:   MEDmeshNodeCoordinateAdvancedRd(fileID, meshname, MED_NO_DT, MED_NO_IT, &vfilter, coordinates);
 96:   /* Read the types of entity sets in the mesh */
 97:   ngeo = MEDmeshnEntity(fileID, meshname, MED_NO_DT, MED_NO_IT, MED_CELL,MED_GEO_ALL, MED_CONNECTIVITY,
 98:                         MED_NODAL, &coordinatechangement, &geotransformation);
 99:   if (ngeo < 1) SETERRQ4(comm, PETSC_ERR_ARG_WRONG, "No cells found in .med v%d.%d.%d mesh file: %s", major, minor, release, filename);
100:   if (ngeo > 2) SETERRQ4(comm, PETSC_ERR_ARG_WRONG, "Currently no support for hybrid meshes in .med v%d.%d.%d mesh file: %s", major, minor, release, filename);
101:   MEDmeshEntityInfo(fileID, meshname, MED_NO_DT, MED_NO_IT, MED_CELL, 1, geotypename, &(geotype[0]));
102:   if (ngeo > 1) {MEDmeshEntityInfo(fileID, meshname, MED_NO_DT, MED_NO_IT, MED_CELL, 2, geotypename, &(geotype[1]));}
103:   else geotype[1] = 0;
104:   /* Determine topological dim and set ID for cells */
105:   cellID = geotype[0]/100 > geotype[1]/100 ? 0 : 1;
106:   facetID = geotype[0]/100 > geotype[1]/100 ? 1 : 0;
107:   meshDim = geotype[cellID] / 100;
108:   numCorners = geotype[cellID] % 100;
109:   /* Partition cells */
110:   numCells = MEDmeshnEntity(fileID, meshname, MED_NO_DT, MED_NO_IT, MED_CELL, geotype[cellID],
111:                             MED_CONNECTIVITY, MED_NODAL,&coordinatechangement, &geotransformation);
112:   PetscLayoutCreate(comm, &cLayout);
113:   PetscLayoutSetSize(cLayout, numCells);
114:   PetscLayoutSetBlockSize(cLayout, 1);
115:   PetscLayoutSetUp(cLayout);
116:   PetscLayoutGetRanges(cLayout, &crange);
117:   numCellsLocal = crange[rank+1]-crange[rank];
118:   MEDfilterBlockOfEntityCr(fileID, numCells, 1, numCorners, MED_ALL_CONSTITUENT, MED_FULL_INTERLACE, MED_COMPACT_STMODE,
119:                                   MED_NO_PROFILE, crange[rank]+1, 1, numCellsLocal, 1, 1, &cfilter);
120:   /* Read cell connectivity */
121:   if (numCells < 0) SETERRQ4(comm, PETSC_ERR_ARG_WRONG, "No cells found in .med v%d.%d.%d mesh file: %s", major, minor, release, filename);
122:   PetscMalloc1(numCellsLocal*numCorners, &medCellList);
123:   MEDmeshElementConnectivityAdvancedRd(fileID, meshname, MED_NO_DT, MED_NO_IT, MED_CELL, geotype[cellID],
124:                                               MED_NODAL, &cfilter, medCellList);
125:   if (sizeof(med_int) > sizeof(PetscInt)) SETERRQ2(comm, PETSC_ERR_ARG_SIZ, "Size of PetscInt %zd less than  size of med_int %zd. Reconfigure PETSc --with-64-bit-indices=1", sizeof(PetscInt), sizeof(med_int));
126:   PetscMalloc1(numCellsLocal*numCorners, &cellList);
127:   for (i = 0; i < numCellsLocal*numCorners; i++) {
128:     cellList[i] = ((PetscInt) medCellList[i]) - 1; /* Correct entity counting */
129:   }
130:   PetscFree(medCellList);
131:   /* Generate the DM */
132:   if (sizeof(med_float) == sizeof(PetscReal)) {
133:     vertcoords = (PetscReal *) coordinates;
134:   } else {
135:     PetscMalloc1(numVerticesLocal*spaceDim, &vertcoords);
136:     for (i = 0; i < numVerticesLocal*spaceDim; i++) vertcoords[i] = (PetscReal) coordinates[i];
137:   }
138:   /* Account for cell inversion */
139:   for (c = 0; c < numCellsLocal; ++c) {
140:     PetscInt *pcone = &cellList[c*numCorners];

142:     if (meshDim == 3) {
143:       /* Hexahedra are inverted */
144:       if (numCorners == 8) {
145:         PetscInt tmp = pcone[4+1];
146:         pcone[4+1] = pcone[4+3];
147:         pcone[4+3] = tmp;
148:       }
149:     }
150:   }
151:   DMPlexCreateFromCellListParallelPetsc(comm, meshDim, numCellsLocal, numVerticesLocal, numVertices, numCorners, interpolate, cellList, spaceDim, vertcoords, &sfVertices, dm);
152:   if (sizeof(med_float) == sizeof(PetscReal)) {
153:     vertcoords = NULL;
154:   } else {
155:     PetscFree(vertcoords);
156:   }
157:   if (ngeo > 1) {
158:     PetscInt        numFacets = 0, numFacetsLocal, numFacetCorners, numFacetsRendezvous;
159:     PetscInt        c, f, v, vStart, joinSize, vertices[8];
160:     PetscInt       *facetList, *facetListRendezvous, *facetIDs, *facetIDsRendezvous, *facetListRemote, *facetIDsRemote;
161:     const PetscInt *frange, *join;
162:     PetscLayout     fLayout;
163:     med_filter      ffilter = MED_FILTER_INIT, fidfilter = MED_FILTER_INIT;
164:     PetscSection    facetSectionRemote, facetSectionIDsRemote;
165:     /* Partition facets */
166:     numFacets = MEDmeshnEntity(fileID, meshname, MED_NO_DT, MED_NO_IT, MED_CELL, geotype[facetID],
167:                                MED_CONNECTIVITY, MED_NODAL,&coordinatechangement, &geotransformation);
168:     numFacetCorners = geotype[facetID] % 100;
169:     PetscLayoutCreate(comm, &fLayout);
170:     PetscLayoutSetSize(fLayout, numFacets);
171:     PetscLayoutSetBlockSize(fLayout, 1);
172:     PetscLayoutSetUp(fLayout);
173:     PetscLayoutGetRanges(fLayout, &frange);
174:     numFacetsLocal = frange[rank+1]-frange[rank];
175:     MEDfilterBlockOfEntityCr(fileID, numFacets, 1, numFacetCorners, MED_ALL_CONSTITUENT, MED_FULL_INTERLACE, MED_COMPACT_STMODE,
176:                                     MED_NO_PROFILE, frange[rank]+1, 1, numFacetsLocal, 1, 1, &ffilter);
177:     MEDfilterBlockOfEntityCr(fileID, numFacets, 1, 1, MED_ALL_CONSTITUENT, MED_FULL_INTERLACE, MED_COMPACT_STMODE,
178:                                     MED_NO_PROFILE, frange[rank]+1, 1, numFacetsLocal, 1, 1, &fidfilter);
179:     DMPlexGetDepthStratum(*dm, 0, &vStart, NULL);
180:     PetscMalloc1(numFacetsLocal, &facetIDs);
181:     PetscMalloc1(numFacetsLocal*numFacetCorners, &facetList);

183:     /* Read facet connectivity */
184:     {
185:       med_int *medFacetList;

187:       PetscMalloc1(numFacetsLocal*numFacetCorners, &medFacetList);
188:       MEDmeshElementConnectivityAdvancedRd(fileID, meshname, MED_NO_DT, MED_NO_IT, MED_CELL, geotype[facetID], MED_NODAL, &ffilter, medFacetList);
189:       for (i = 0; i < numFacetsLocal*numFacetCorners; i++) {
190:         facetList[i] = ((PetscInt) medFacetList[i]) - 1 ; /* Correct entity counting */
191:       }
192:       PetscFree(medFacetList);
193:     }

195:     /* Read facet IDs */
196:     {
197:       med_int *medFacetIDs;

199:       PetscMalloc1(numFacetsLocal, &medFacetIDs);
200:       MEDmeshEntityAttributeAdvancedRd(fileID, meshname, MED_FAMILY_NUMBER, MED_NO_DT, MED_NO_IT, MED_CELL, geotype[facetID], &fidfilter, medFacetIDs);
201:       for (i = 0; i < numFacetsLocal; i++) {
202:         facetIDs[i] = (PetscInt) medFacetIDs[i];
203:       }
204:       PetscFree(medFacetIDs);
205:     }

207:     /* Send facets and IDs to a rendezvous partition that is based on the initial vertex partitioning. */
208:     {
209:       PetscInt           r;
210:       DMLabel            lblFacetRendezvous, lblFacetMigration;
211:       PetscSection       facetSection, facetSectionRendezvous;
212:       PetscSF            sfProcess, sfFacetMigration;
213:       const PetscSFNode *remoteVertices;
214:       DMLabelCreate(PETSC_COMM_SELF, "Facet Rendezvous", &lblFacetRendezvous);
215:       DMLabelCreate(PETSC_COMM_SELF, "Facet Migration", &lblFacetMigration);
216:       PetscSFGetGraph(sfVertices, NULL, NULL, NULL, &remoteVertices);
217:       for (f = 0; f < numFacetsLocal; f++) {
218:         for (v = 0; v < numFacetCorners; v++) {
219:           /* Find vertex owner on rendezvous partition and mark in label */
220:           const PetscInt vertex = facetList[f*numFacetCorners+v];
221:           r = rank; while (vrange[r] > vertex) r--; while (vrange[r + 1] < vertex) r++;
222:           DMLabelSetValue(lblFacetRendezvous, f, r);
223:         }
224:       }
225:       /* Build a global process SF */
226:       PetscSFCreate(comm,&sfProcess);
227:       PetscSFSetGraphWithPattern(sfProcess,NULL,PETSCSF_PATTERN_ALLTOALL);
228:       PetscObjectSetName((PetscObject) sfProcess, "Process SF");
229:       /* Convert facet rendezvous label into SF for migration */
230:       DMPlexPartitionLabelInvert(*dm, lblFacetRendezvous, sfProcess, lblFacetMigration);
231:       DMPlexPartitionLabelCreateSF(*dm, lblFacetMigration, &sfFacetMigration);
232:       /* Migrate facet connectivity data */
233:       PetscSectionCreate(comm, &facetSection);
234:       PetscSectionSetChart(facetSection, 0, numFacetsLocal);
235:       for (f = 0; f < numFacetsLocal; f++) {PetscSectionSetDof(facetSection, f, numFacetCorners);}
236:       PetscSectionSetUp(facetSection);
237:       PetscSectionCreate(comm, &facetSectionRendezvous);
238:       DMPlexDistributeData(*dm, sfFacetMigration, facetSection, MPIU_INT, facetList, facetSectionRendezvous, (void**) &facetListRendezvous);
239:       /* Migrate facet IDs */
240:       PetscSFGetGraph(sfFacetMigration, NULL, &numFacetsRendezvous, NULL, NULL);
241:       PetscMalloc1(numFacetsRendezvous, &facetIDsRendezvous);
242:       PetscSFBcastBegin(sfFacetMigration, MPIU_INT, facetIDs, facetIDsRendezvous);
243:       PetscSFBcastEnd(sfFacetMigration, MPIU_INT, facetIDs, facetIDsRendezvous);
244:       /* Clean up */
245:       DMLabelDestroy(&lblFacetRendezvous);
246:       DMLabelDestroy(&lblFacetMigration);
247:       PetscSFDestroy(&sfProcess);
248:       PetscSFDestroy(&sfFacetMigration);
249:       PetscSectionDestroy(&facetSection);
250:       PetscSectionDestroy(&facetSectionRendezvous);
251:     }

253:     /* On the rendevouz partition we build a vertex-wise section/array of facets and IDs. */
254:     {
255:       PetscInt               sizeVertexFacets, offset, sizeFacetIDsRemote;
256:       PetscInt              *vertexFacets, *vertexIdx, *vertexFacetIDs;
257:       PetscSection           facetSectionVertices, facetSectionIDs;
258:       ISLocalToGlobalMapping ltogVertexNumbering;
259:       PetscSectionCreate(comm, &facetSectionVertices);
260:       PetscSectionSetChart(facetSectionVertices, 0, numVerticesLocal);
261:       PetscSectionCreate(comm, &facetSectionIDs);
262:       PetscSectionSetChart(facetSectionIDs, 0, numVerticesLocal);
263:       for (f = 0; f < numFacetsRendezvous*numFacetCorners; f++) {
264:         const PetscInt vertex = facetListRendezvous[f];
265:         if (vrange[rank] <= vertex && vertex < vrange[rank+1]) {
266:           PetscSectionAddDof(facetSectionIDs, vertex-vrange[rank], 1);
267:           PetscSectionAddDof(facetSectionVertices, vertex-vrange[rank], numFacetCorners);
268:         }
269:       }
270:       PetscSectionSetUp(facetSectionVertices);
271:       PetscSectionSetUp(facetSectionIDs);
272:       PetscSectionGetStorageSize(facetSectionVertices, &sizeVertexFacets);
273:       PetscSectionGetStorageSize(facetSectionVertices, &sizeFacetIDsRemote);
274:       PetscMalloc1(sizeVertexFacets, &vertexFacets);
275:       PetscMalloc1(sizeFacetIDsRemote, &vertexFacetIDs);
276:       PetscCalloc1(numVerticesLocal, &vertexIdx);
277:       for (f = 0; f < numFacetsRendezvous; f++) {
278:         for (c = 0; c < numFacetCorners; c++) {
279:           const PetscInt vertex = facetListRendezvous[f*numFacetCorners+c];
280:           if (vrange[rank] <= vertex && vertex < vrange[rank+1]) {
281:             /* Flip facet connectivities and IDs to a vertex-wise layout */
282:             PetscSectionGetOffset(facetSectionVertices, vertex-vrange[rank], &offset);
283:             offset += vertexIdx[vertex-vrange[rank]] * numFacetCorners;
284:             PetscArraycpy(&(vertexFacets[offset]), &(facetListRendezvous[f*numFacetCorners]), numFacetCorners);
285:             PetscSectionGetOffset(facetSectionIDs, vertex-vrange[rank], &offset);
286:             offset += vertexIdx[vertex-vrange[rank]];
287:             vertexFacetIDs[offset] = facetIDsRendezvous[f];
288:             vertexIdx[vertex-vrange[rank]]++;
289:           }
290:         }
291:       }
292:       /* Distribute the vertex-wise facet connectivities over the vertexSF */
293:       PetscSectionCreate(comm, &facetSectionRemote);
294:       DMPlexDistributeData(*dm, sfVertices, facetSectionVertices, MPIU_INT, vertexFacets, facetSectionRemote, (void**) &facetListRemote);
295:       PetscSectionCreate(comm, &facetSectionIDsRemote);
296:       DMPlexDistributeData(*dm, sfVertices, facetSectionIDs, MPIU_INT, vertexFacetIDs, facetSectionIDsRemote, (void**) &facetIDsRemote);
297:       /* Convert facet connectivities to local vertex numbering */
298:       PetscSectionGetStorageSize(facetSectionRemote, &sizeVertexFacets);
299:       ISLocalToGlobalMappingCreateSF(sfVertices, vrange[rank], &ltogVertexNumbering);
300:       ISGlobalToLocalMappingApplyBlock(ltogVertexNumbering, IS_GTOLM_MASK, sizeVertexFacets, facetListRemote, NULL, facetListRemote);
301:       /* Clean up */
302:       PetscFree(vertexFacets);
303:       PetscFree(vertexIdx);
304:       PetscFree(vertexFacetIDs);
305:       PetscSectionDestroy(&facetSectionVertices);
306:       PetscSectionDestroy(&facetSectionIDs);
307:       ISLocalToGlobalMappingDestroy(&ltogVertexNumbering);
308:     }
309:     {
310:       PetscInt offset, dof;
311:       DMLabel lblFaceSets;
312:       PetscBool verticesLocal;
313:       /* Identify and mark facets locally with facet joins */
314:       DMCreateLabel(*dm, "Face Sets");
315:       DMGetLabel(*dm, "Face Sets", &lblFaceSets);
316:       /* We need to set a new default value here, since -1 is a legitimate facet ID */
317:       DMLabelSetDefaultValue(lblFaceSets, -666666666);
318:       for (v = 0; v < numVerticesLocal; v++) {
319:         PetscSectionGetOffset(facetSectionRemote, v, &offset);
320:         PetscSectionGetDof(facetSectionRemote, v, &dof);
321:         for (f = 0; f < dof; f += numFacetCorners) {
322:           for (verticesLocal = PETSC_TRUE, c = 0; c < numFacetCorners; ++c) {
323:             if (facetListRemote[offset+f+c] < 0) {verticesLocal = PETSC_FALSE; break;}
324:             vertices[c] = vStart + facetListRemote[offset+f+c];
325:           }
326:           if (verticesLocal) {
327:             DMPlexGetFullJoin(*dm, numFacetCorners, (const PetscInt*)vertices, &joinSize, &join);
328:             if (joinSize == 1) {
329:               DMLabelSetValue(lblFaceSets, join[0], facetIDsRemote[(offset+f) / numFacetCorners]);
330:             }
331:             DMPlexRestoreJoin(*dm, numFacetCorners, (const PetscInt*)vertices, &joinSize, &join);
332:           }
333:         }
334:       }
335:     }
336:     PetscFree(facetList);
337:     PetscFree(facetListRendezvous);
338:     PetscFree(facetListRemote);
339:     PetscFree(facetIDs);
340:     PetscFree(facetIDsRendezvous);
341:     PetscFree(facetIDsRemote);
342:     PetscLayoutDestroy(&fLayout);
343:     PetscSectionDestroy(&facetSectionRemote);
344:     PetscSectionDestroy(&facetSectionIDsRemote);
345:   }
346:   MEDfileClose(fileID);
347:   PetscFree(coordinates);
348:   PetscFree(cellList);
349:   PetscLayoutDestroy(&vLayout);
350:   PetscLayoutDestroy(&cLayout);
351:   PetscSFDestroy(&sfVertices);
352:   return(0);
353: #else
354:   SETERRQ(comm, PETSC_ERR_SUP, "This method requires Med mesh reader support. Reconfigure using --download-med");
355: #endif
356: }