Actual source code: plexmed.c
petsc-3.12.5 2020-03-29
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: int *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(int)) {
126: cellList = (int *) medCellList;
127: medCellList = NULL;
128: } else {
129: PetscMalloc1(numCellsLocal*numCorners, &cellList);
130: for (i = 0; i < numCellsLocal*numCorners; i++) cellList[i] = medCellList[i];
131: PetscFree(medCellList);
132: }
133: for (i = 0; i < numCellsLocal*numCorners; i++) cellList[i]--; /* Correct entity counting */
134: /* Generate the DM */
135: if (sizeof(med_float) == sizeof(PetscReal)) {
136: vertcoords = (PetscReal *) coordinates;
137: } else {
138: PetscMalloc1(numVerticesLocal*spaceDim, &vertcoords);
139: for (i = 0; i < numVerticesLocal*spaceDim; i++) vertcoords[i] = coordinates[i];
140: }
141: /* Account for cell inversion */
142: for (c = 0; c < numCellsLocal; ++c) {
143: PetscInt *pcone = &cellList[c*numCorners];
145: if (meshDim == 3) {
146: /* Hexahedra are inverted */
147: if (numCorners == 8) {
148: PetscInt tmp = pcone[4+1];
149: pcone[4+1] = pcone[4+3];
150: pcone[4+3] = tmp;
151: }
152: }
153: }
154: DMPlexCreateFromCellListParallel(comm, meshDim, numCellsLocal, numVerticesLocal, numCorners, interpolate, cellList, spaceDim, coordinates, &sfVertices, dm);
155: if (sizeof(med_float) == sizeof(PetscReal)) {
156: vertcoords = NULL;
157: } else {
158: PetscFree(vertcoords);
159: }
161: if (ngeo > 1) {
162: PetscInt numFacets = 0, numFacetsLocal, numFacetCorners, numFacetsRendezvous;
163: PetscInt c, f, v, vStart, joinSize, vertices[8];
164: med_int *medFacetList;
165: PetscInt *facetList, *facetListRendezvous, *facetIDs, *facetIDsRendezvous, *facetListRemote, *facetIDsRemote;
166: const PetscInt *frange, *join;
167: PetscLayout fLayout;
168: med_filter ffilter = MED_FILTER_INIT, fidfilter = MED_FILTER_INIT;
169: PetscSection facetSectionRemote, facetSectionIDsRemote;
170: /* Partition facets */
171: numFacets = MEDmeshnEntity(fileID, meshname, MED_NO_DT, MED_NO_IT, MED_CELL, geotype[facetID],
172: MED_CONNECTIVITY, MED_NODAL,&coordinatechangement, &geotransformation);
173: numFacetCorners = geotype[facetID] % 100;
174: PetscLayoutCreate(comm, &fLayout);
175: PetscLayoutSetSize(fLayout, numFacets);
176: PetscLayoutSetBlockSize(fLayout, 1);
177: PetscLayoutSetUp(fLayout);
178: PetscLayoutGetRanges(fLayout, &frange);
179: numFacetsLocal = frange[rank+1]-frange[rank];
180: MEDfilterBlockOfEntityCr(fileID, numFacets, 1, numFacetCorners, MED_ALL_CONSTITUENT, MED_FULL_INTERLACE, MED_COMPACT_STMODE,
181: MED_NO_PROFILE, frange[rank]+1, 1, numFacetsLocal, 1, 1, &ffilter);
182: MEDfilterBlockOfEntityCr(fileID, numFacets, 1, 1, MED_ALL_CONSTITUENT, MED_FULL_INTERLACE, MED_COMPACT_STMODE,
183: MED_NO_PROFILE, frange[rank]+1, 1, numFacetsLocal, 1, 1, &fidfilter);
184: /* Read facet connectivity */
185: DMPlexGetDepthStratum(*dm, 0, &vStart, NULL);
186: PetscMalloc1(numFacetsLocal*numFacetCorners, &medFacetList);
187: PetscCalloc1(numFacetsLocal, &facetIDs);
188: MEDmeshElementConnectivityAdvancedRd(fileID, meshname, MED_NO_DT, MED_NO_IT, MED_CELL, geotype[facetID],
189: MED_NODAL, &ffilter, medFacetList);
190: if (sizeof(med_int) == sizeof(PetscInt)) {
191: facetList = (PetscInt *) medFacetList;
192: medFacetList = NULL;
193: } else {
194: PetscMalloc1(numFacetsLocal*numFacetCorners, &facetList);
195: for (i = 0; i < numFacetsLocal*numFacetCorners; i++) facetList[i] = medFacetList[i];
196: PetscFree(medFacetList);
197: }
198: for (i = 0; i < numFacetsLocal*numFacetCorners; i++) facetList[i]--; /* Correct entity counting */
199: /* Read facet IDs */
200: MEDmeshEntityAttributeAdvancedRd(fileID, meshname, MED_FAMILY_NUMBER, MED_NO_DT, MED_NO_IT, MED_CELL,
201: geotype[facetID], &fidfilter, facetIDs);
202: {
203: /* Send facets and IDs to a rendezvous partition that is based on the initial vertex partitioning. */
204: PetscInt r;
205: DMLabel lblFacetRendezvous, lblFacetMigration;
206: PetscSection facetSection, facetSectionRendezvous;
207: PetscSF sfProcess, sfFacetMigration;
208: const PetscSFNode *remoteVertices;
209: DMLabelCreate(PETSC_COMM_SELF, "Facet Rendezvous", &lblFacetRendezvous);
210: DMLabelCreate(PETSC_COMM_SELF, "Facet Migration", &lblFacetMigration);
211: PetscSFGetGraph(sfVertices, NULL, NULL, NULL, &remoteVertices);
212: for (f = 0; f < numFacetsLocal; f++) {
213: for (v = 0; v < numFacetCorners; v++) {
214: /* Find vertex owner on rendezvous partition and mark in label */
215: const PetscInt vertex = facetList[f*numFacetCorners+v];
216: r = rank; while (vrange[r] > vertex) r--; while (vrange[r + 1] < vertex) r++;
217: DMLabelSetValue(lblFacetRendezvous, f, r);
218: }
219: }
220: /* Build a global process SF */
221: PetscSFCreate(comm,&sfProcess);
222: PetscSFSetGraphWithPattern(sfProcess,NULL,PETSCSF_PATTERN_ALLTOALL);
223: PetscObjectSetName((PetscObject) sfProcess, "Process SF");
224: /* Convert facet rendezvous label into SF for migration */
225: DMPlexPartitionLabelInvert(*dm, lblFacetRendezvous, sfProcess, lblFacetMigration);
226: DMPlexPartitionLabelCreateSF(*dm, lblFacetMigration, &sfFacetMigration);
227: /* Migrate facet connectivity data */
228: PetscSectionCreate(comm, &facetSection);
229: PetscSectionSetChart(facetSection, 0, numFacetsLocal);
230: for (f = 0; f < numFacetsLocal; f++) {PetscSectionSetDof(facetSection, f, numFacetCorners);}
231: PetscSectionSetUp(facetSection);
232: PetscSectionCreate(comm, &facetSectionRendezvous);
233: DMPlexDistributeData(*dm, sfFacetMigration, facetSection, MPIU_INT, facetList, facetSectionRendezvous, (void**) &facetListRendezvous);
234: /* Migrate facet IDs */
235: PetscSFGetGraph(sfFacetMigration, NULL, &numFacetsRendezvous, NULL, NULL);
236: PetscMalloc1(numFacetsRendezvous, &facetIDsRendezvous);
237: PetscSFBcastBegin(sfFacetMigration, MPIU_INT, facetIDs, facetIDsRendezvous);
238: PetscSFBcastEnd(sfFacetMigration, MPIU_INT, facetIDs, facetIDsRendezvous);
239: /* Clean up */
240: DMLabelDestroy(&lblFacetRendezvous);
241: DMLabelDestroy(&lblFacetMigration);
242: PetscSFDestroy(&sfProcess);
243: PetscSFDestroy(&sfFacetMigration);
244: PetscSectionDestroy(&facetSection);
245: PetscSectionDestroy(&facetSectionRendezvous);
246: }
247: {
248: /* On the rendevouz partition we build a vertex-wise section/array of facets and IDs. */
249: PetscInt sizeVertexFacets, offset, sizeFacetIDsRemote;
250: PetscInt *vertexFacets, *vertexIdx, *vertexFacetIDs;
251: PetscSection facetSectionVertices, facetSectionIDs;
252: ISLocalToGlobalMapping ltogVertexNumbering;
253: PetscSectionCreate(comm, &facetSectionVertices);
254: PetscSectionSetChart(facetSectionVertices, 0, numVerticesLocal);
255: PetscSectionCreate(comm, &facetSectionIDs);
256: PetscSectionSetChart(facetSectionIDs, 0, numVerticesLocal);
257: for (f = 0; f < numFacetsRendezvous*numFacetCorners; f++) {
258: const PetscInt vertex = facetListRendezvous[f];
259: if (vrange[rank] <= vertex && vertex < vrange[rank+1]) {
260: PetscSectionAddDof(facetSectionIDs, vertex-vrange[rank], 1);
261: PetscSectionAddDof(facetSectionVertices, vertex-vrange[rank], numFacetCorners);
262: }
263: }
264: PetscSectionSetUp(facetSectionVertices);
265: PetscSectionSetUp(facetSectionIDs);
266: PetscSectionGetStorageSize(facetSectionVertices, &sizeVertexFacets);
267: PetscSectionGetStorageSize(facetSectionVertices, &sizeFacetIDsRemote);
268: PetscMalloc1(sizeVertexFacets, &vertexFacets);
269: PetscMalloc1(sizeFacetIDsRemote, &vertexFacetIDs);
270: PetscCalloc1(numVerticesLocal, &vertexIdx);
271: for (f = 0; f < numFacetsRendezvous; f++) {
272: for (c = 0; c < numFacetCorners; c++) {
273: const PetscInt vertex = facetListRendezvous[f*numFacetCorners+c];
274: if (vrange[rank] <= vertex && vertex < vrange[rank+1]) {
275: /* Flip facet connectivities and IDs to a vertex-wise layout */
276: PetscSectionGetOffset(facetSectionVertices, vertex-vrange[rank], &offset);
277: offset += vertexIdx[vertex-vrange[rank]] * numFacetCorners;
278: PetscArraycpy(&(vertexFacets[offset]), &(facetListRendezvous[f*numFacetCorners]), numFacetCorners);
279: PetscSectionGetOffset(facetSectionIDs, vertex-vrange[rank], &offset);
280: offset += vertexIdx[vertex-vrange[rank]];
281: vertexFacetIDs[offset] = facetIDsRendezvous[f];
282: vertexIdx[vertex-vrange[rank]]++;
283: }
284: }
285: }
286: /* Distribute the vertex-wise facet connectivities over the vertexSF */
287: PetscSectionCreate(comm, &facetSectionRemote);
288: DMPlexDistributeData(*dm, sfVertices, facetSectionVertices, MPIU_INT, vertexFacets, facetSectionRemote, (void**) &facetListRemote);
289: PetscSectionCreate(comm, &facetSectionIDsRemote);
290: DMPlexDistributeData(*dm, sfVertices, facetSectionIDs, MPIU_INT, vertexFacetIDs, facetSectionIDsRemote, (void**) &facetIDsRemote);
291: /* Convert facet connectivities to local vertex numbering */
292: PetscSectionGetStorageSize(facetSectionRemote, &sizeVertexFacets);
293: ISLocalToGlobalMappingCreateSF(sfVertices, vrange[rank], <ogVertexNumbering);
294: ISGlobalToLocalMappingApplyBlock(ltogVertexNumbering, IS_GTOLM_MASK, sizeVertexFacets, facetListRemote, NULL, facetListRemote);
295: /* Clean up */
296: PetscFree(vertexFacets);
297: PetscFree(vertexIdx);
298: PetscFree(vertexFacetIDs);
299: PetscSectionDestroy(&facetSectionVertices);
300: PetscSectionDestroy(&facetSectionIDs);
301: ISLocalToGlobalMappingDestroy(<ogVertexNumbering);
302: }
303: {
304: PetscInt offset, dof;
305: DMLabel lblFaceSets;
306: PetscBool verticesLocal;
307: /* Identify and mark facets locally with facet joins */
308: DMCreateLabel(*dm, "Face Sets");
309: DMGetLabel(*dm, "Face Sets", &lblFaceSets);
310: /* We need to set a new default value here, since -1 is a legitimate facet ID */
311: DMLabelSetDefaultValue(lblFaceSets, -666666666);
312: for (v = 0; v < numVerticesLocal; v++) {
313: PetscSectionGetOffset(facetSectionRemote, v, &offset);
314: PetscSectionGetDof(facetSectionRemote, v, &dof);
315: for (f = 0; f < dof; f += numFacetCorners) {
316: for (verticesLocal = PETSC_TRUE, c = 0; c < numFacetCorners; ++c) {
317: if (facetListRemote[offset+f+c] < 0) {verticesLocal = PETSC_FALSE; break;}
318: vertices[c] = vStart + facetListRemote[offset+f+c];
319: }
320: if (verticesLocal) {
321: DMPlexGetFullJoin(*dm, numFacetCorners, (const PetscInt*)vertices, &joinSize, &join);
322: if (joinSize == 1) {
323: DMLabelSetValue(lblFaceSets, join[0], facetIDsRemote[(offset+f) / numFacetCorners]);
324: }
325: DMPlexRestoreJoin(*dm, numFacetCorners, (const PetscInt*)vertices, &joinSize, &join);
326: }
327: }
328: }
329: }
330: PetscFree(facetList);
331: PetscFree(facetListRendezvous);
332: PetscFree(facetListRemote);
333: PetscFree(facetIDs);
334: PetscFree(facetIDsRendezvous);
335: PetscFree(facetIDsRemote);
336: PetscLayoutDestroy(&fLayout);
337: PetscSectionDestroy(&facetSectionRemote);
338: PetscSectionDestroy(&facetSectionIDsRemote);
339: }
340: MEDfileClose(fileID);
341: PetscFree(coordinates);
342: PetscFree(cellList);
343: PetscLayoutDestroy(&vLayout);
344: PetscLayoutDestroy(&cLayout);
345: PetscSFDestroy(&sfVertices);
346: #else
347: SETERRQ(comm, PETSC_ERR_SUP, "This method requires Med mesh reader support. Reconfigure using --download-med");
348: #endif
349: return(0);
350: }