void MESH_Renumber(Mesh_ptr mesh, int renum_type, MType mtype) { MVertex_ptr mv, v0=NULL; MEdge_ptr me, e0=NULL; MFace_ptr mf, f0=NULL; MRegion_ptr mr, r0=NULL; int idx, idx2, idx3; int i, j; int done; MAttrib_ptr vidatt; List_ptr vlist; double xyz[3]; double rval; int bandwidth, maxbandwidth1, maxbandwidth2; double avebandwidth1, avebandwidth2; void *pval; if (renum_type == 0) { if (mtype == MVERTEX || mtype == MALLTYPE) { int nv = 0; idx = 0; while ((mv = MESH_Next_Vertex(mesh,&idx))) MV_Set_ID(mv,++nv); } if (mtype == MEDGE || mtype == MALLTYPE) { int ne = 0; idx = 0; while ((me = MESH_Next_Edge(mesh,&idx))) ME_Set_ID(me,++ne); } if (mtype == MFACE || mtype == MALLTYPE) { int nf = 0; idx = 0; while ((mf = MESH_Next_Face(mesh,&idx))) MF_Set_ID(mf,++nf); } if (mtype == MREGION || mtype == MALLTYPE) { int nr = 0; idx = 0; while ((mr = MESH_Next_Region(mesh,&idx))) MR_Set_ID(mr,++nr); } } else if (renum_type == 1) { double minx, miny, minz; int minid, maxid; int *nadj, *newmap, *adj, *offset, nconn; int nalloc, depth, maxwidth; #ifdef MSTK_USE_MARKERS int mkid = MSTK_GetMarker(); #else MAttrib_ptr mkatt = MAttrib_New(mesh, "mkatt", INT, MALLTYPE); #endif if (mtype == MVERTEX || mtype == MALLTYPE) { int nv = MESH_Num_Vertices(mesh); /* Compute a graph of vertex connections across elements (faces for surface meshes, regions for solid meshes */ /* Start with the vertex in the lower leftmost corner */ minx = miny = minz = 1.0e+12; v0 = NULL; idx = 0; while ((mv = MESH_Next_Vertex(mesh,&idx))) { MV_Coords(mv,xyz); if (xyz[0] <= minx && xyz[1] <= miny && xyz[2] <= minz) { minx = xyz[0]; miny = xyz[1]; minz = xyz[2]; v0 = mv; } } nadj = (int *) malloc(nv*sizeof(int)); nalloc = nv*5; adj = (int *) malloc(nalloc*sizeof(int)); if (!MESH_Num_Regions(mesh)) { int nentries = 0; i = 0; idx = 0; while ((mv = MESH_Next_Vertex(mesh,&idx))) { List_ptr vfaces, adjvlist; MFace_ptr vf; MVertex_ptr adjv; adjvlist = List_New(0); vfaces = MV_Faces(mv); idx2 = 0; while ((vf = List_Next_Entry(vfaces,&idx2))) { List_ptr fverts = MF_Vertices(vf,1,0); idx3 = 0; while ((adjv = List_Next_Entry(fverts,&idx3))) { if (adjv != mv) { int vmarked; #ifdef MSTK_USE_MARKERS vmarked = MEnt_IsMarked(adjv,mkid); #else MEnt_Get_AttVal(adjv, mkatt, &vmarked, &rval, &pval); #endif if (!vmarked) { List_Add(adjvlist,adjv); #ifdef MSTK_USE_MARKERS MEnt_Mark(adjv,mkid); #else MEnt_Set_AttVal(adjv, mkatt, 1, 0.0, NULL); #endif } } } List_Delete(fverts); } List_Delete(vfaces); #ifdef MSTK_USE_MARKERS List_Unmark(adjvlist,mkid); #endif nadj[i] = List_Num_Entries(adjvlist); if (nentries+nadj[i] > nalloc) { nalloc *= 2; adj = (int *) realloc(adj,nalloc*sizeof(int)); } idx2 = 0; while ((adjv = List_Next_Entry(adjvlist,&idx2))) adj[nentries++] = MV_ID(adjv)-1; List_Delete(adjvlist); i++; } } else { int nentries = 0; i = 0; idx = 0; while ((mv = MESH_Next_Vertex(mesh,&idx))) { List_ptr vregions, adjvlist; MRegion_ptr vr; MVertex_ptr adjv; adjvlist = List_New(0); vregions = MV_Regions(mv); idx2 = 0; while ((vr = List_Next_Entry(vregions,&idx2))) { List_ptr rverts = MR_Vertices(vr); idx3 = 0; while ((adjv = List_Next_Entry(rverts,&idx3))) { if (adjv != mv) { int vmarked; #ifdef MSTK_USE_MARKERS vmarked = MEnt_IsMarked(adjv,mkid); #else MEnt_Get_AttVal(adjv, mkatt, &vmarked, &rval, &pval); #endif if (!vmarked) { List_Add(adjvlist,adjv); #ifdef MSTK_USE_MARKERS MEnt_Mark(adjv,mkid); #else MEnt_Set_AttVal(adjv, mkatt, 1, 0.0, NULL); #endif } } } List_Delete(rverts); } List_Delete(vregions); #ifdef MSTK_USE_MARKERS List_Unmark(adjvlist,mkid); #endif nadj[i] = List_Num_Entries(adjvlist); if (nentries+nadj[i] > nalloc) { nalloc *= 2; adj = (int *) realloc(adj,nalloc*sizeof(int)); } idx2 = 0; while ((adjv = List_Next_Entry(adjvlist,&idx2))) adj[nentries++] = MV_ID(adjv)-1; List_Delete(adjvlist); i++; } } /* Compute offsets into adj array */ offset = (int *) malloc(nv*sizeof(int)); offset[0] = 0; for (i = 1; i < nv; i++) offset[i] = offset[i-1] + nadj[i-1]; /* Compute maximum bandwidth before renumbering */ maxbandwidth1 = 0; avebandwidth1 = 0; for (i = 0; i < nv; i++) { int off = offset[i]; int curid = i; for (j = 0; j < nadj[i]; j++) { int adjid = adj[off+j]; int diff = abs(adjid-curid); maxbandwidth1 = (diff > maxbandwidth1) ? diff : maxbandwidth1; avebandwidth1 += diff; nconn++; } } nconn = offset[nv-1]+nadj[nv-1]; avebandwidth1 /= nconn; fprintf(stderr, "Ave vertex ID difference on elements before renumbering: %-lf\n", avebandwidth1); fprintf(stderr, "Max vertex ID difference on elements before renumbering: %-d\n", maxbandwidth1); fprintf(stderr,"\n"); newmap = (int *) malloc(nv*sizeof(int)); Graph_Renumber_GPS(nv, MV_ID(v0)-1, nadj, adj, newmap, &depth, &maxwidth); /* Compute bandwidth after renumbering */ maxbandwidth2 = 0; avebandwidth2 = 0; for (i = 0; i < nv; i++) { int off = offset[i]; int curid = newmap[i]; for (j = 0; j < nadj[i]; j++) { int adjid = newmap[adj[off+j]]; int diff = abs(adjid-curid); maxbandwidth2 = (diff > maxbandwidth2) ? diff : maxbandwidth2; avebandwidth2 += diff; nconn++; } } nconn = offset[nv-1]+nadj[nv-1]; avebandwidth2 /= nconn; if (maxbandwidth2 < maxbandwidth1 && avebandwidth2 < avebandwidth1) { /* Renumber */ idx = 0; i = 0; while ((mv = MESH_Next_Vertex(mesh,&idx))) { MV_Set_ID(mv,newmap[i]+1); i++; } fprintf(stderr, "Ave vertex ID difference on elements after renumbering: %-lf\n", avebandwidth2); fprintf(stderr, "Max vertex ID difference on elements after renumbering: %-d\n", maxbandwidth2); } else { nv = 0; idx = 0; while ((mv = MESH_Next_Vertex(mesh,&idx))) MV_Set_ID(mv,++nv); fprintf(stderr,"Bandwidth did not improve. Keeping old numbering with gaps eliminated\n"); } fprintf(stderr,"\n\n\n"); free(nadj); free(adj); free(offset); free(newmap); } /* Reorder edges according to a breadth first algorithm applied to edges (differs from RCM in that it does not add adjacent nodes in ascending order of their valence) */ if (mtype == MEDGE || mtype == MALLTYPE) { int ne = MESH_Num_Edges(mesh); MEdge_ptr ve; List_ptr elist; /************************* renumbering code ****************************/ ne = MESH_Num_Edges(mesh); if (mtype == MALLTYPE) { /* RCM algorithm already applied on the vertices. Use an edge connected to the starting vertex as the first edge */ List_ptr vedges = MV_Edges(v0); e0 = List_Entry(vedges,0); List_Delete(vedges); } else { /* Find the edge whose mid point is a minimum point */ minx = miny = minz = 1.0e+12; e0 = NULL; idx = 0; while ((me = MESH_Next_Edge(mesh,&idx))) { double exyz[2][3]; MV_Coords(ME_Vertex(me,0),exyz[0]); MV_Coords(ME_Vertex(me,1),exyz[1]); xyz[0] = (exyz[0][0]+exyz[1][0])/2.0; xyz[1] = (exyz[0][1]+exyz[1][1])/2.0; xyz[2] = (exyz[0][2]+exyz[1][2])/2.0; if (xyz[0] < minx && xyz[1] < miny && xyz[2] < minz) { minx = xyz[0]; miny = xyz[1]; minz = xyz[2]; e0 = me; } } } nadj = (int *) malloc(ne*sizeof(int)); nalloc = ne*5; adj = (int *) malloc(nalloc*sizeof(int)); if (!MESH_Num_Regions(mesh)) { int nentries = 0; i = 0; idx = 0; while ((me = MESH_Next_Edge(mesh,&idx))) { List_ptr efaces, adjelist; MFace_ptr ef; MEdge_ptr adje; adjelist = List_New(0); efaces = ME_Faces(me); idx2 = 0; while ((ef = List_Next_Entry(efaces,&idx2))) { List_ptr fedges = MF_Edges(ef,1,0); idx3 = 0; while ((adje = List_Next_Entry(fedges,&idx3))) { if (adje != me) { int emarked; #ifdef MSTK_USE_MARKERS emarked = MEnt_IsMarked(adje,mkid); #else MEnt_Get_AttVal(adje, mkatt, &emarked, &rval, &pval); #endif if (!emarked) { List_Add(adjelist,adje); #ifdef MSTK_USE_MARKERS MEnt_Mark(adje,mkid); #else MEnt_Set_AttVal(adje, mkatt, 1, 0.0, NULL); #endif } } } List_Delete(fedges); } List_Delete(efaces); #ifdef MSTK_USE_MARKERS List_Unmark(adjelist,mkid); #endif nadj[i] = List_Num_Entries(adjelist); if (nentries+nadj[i] > nalloc) { nalloc *= 2; adj = (int *) realloc(adj,nalloc*sizeof(int)); } idx2 = 0; while ((adje = List_Next_Entry(adjelist,&idx2))) adj[nentries++] = ME_ID(adje)-1; List_Delete(adjelist); i++; } } else { int nentries = 0; i = 0; idx = 0; while ((me = MESH_Next_Edge(mesh,&idx))) { List_ptr eregions, adjelist; MRegion_ptr er; MEdge_ptr adje; adjelist = List_New(0); eregions = ME_Regions(me); idx2 = 0; while ((er = List_Next_Entry(eregions,&idx2))) { List_ptr redges = MR_Edges(er); idx3 = 0; while ((adje = List_Next_Entry(redges,&idx3))) { if (adje != me) { int emarked; #ifdef MSTK_USE_MARKERS emarked = MEnt_IsMarked(adje,mkid); #else MEnt_Get_AttVal(adje, mkatt, &emarked, &rval, &pval); #endif if (!emarked) { List_Add(adjelist,adje); #ifdef MSTK_USE_MARKERS MEnt_Mark(adje,mkid); #else MEnt_Set_AttVal(adje, mkatt, 1, 0.0, NULL); #endif } } } List_Delete(redges); } List_Delete(eregions); #ifdef MSTK_USE_MARKERS List_Unmark(adjelist,mkid); #endif nadj[i] = List_Num_Entries(adjelist); if (nentries+nadj[i] > nalloc) { nalloc *= 2; adj = (int *) realloc(adj,nalloc*sizeof(int)); } idx2 = 0; while ((adje = List_Next_Entry(adjelist,&idx2))) adj[nentries++] = ME_ID(adje)-1; List_Delete(adjelist); i++; } } /* Compute offsets into adj array */ offset = (int *) malloc(ne*sizeof(int)); offset[0] = 0; for (i = 1; i < ne; i++) offset[i] = offset[i-1] + nadj[i-1]; /* Compute maximum bandwidth before renumbering */ maxbandwidth1 = 0; avebandwidth1 = 0; for (i = 0; i < ne; i++) { int off = offset[i]; int curid = i; for (j = 0; j < nadj[i]; j++) { int adjid = adj[off+j]; int diff = abs(adjid-curid); maxbandwidth1 = (diff > maxbandwidth1) ? diff : maxbandwidth1; avebandwidth1 += diff; nconn++; } } nconn = offset[ne-1]+nadj[ne-1]; avebandwidth1 /= nconn; fprintf(stderr, "Ave edge ID difference on elements before renumbering: %-lf\n", avebandwidth1); fprintf(stderr, "Max edge ID difference on elements before renumbering: %-d\n", maxbandwidth1); fprintf(stderr,"\n"); /* Call Graph Renumbering algorithm */ newmap = (int *) malloc(ne*sizeof(int)); Graph_Renumber_GPS(ne, ME_ID(e0)-1, nadj, adj, newmap, &depth, &maxwidth); /* Compute bandwidth after renumbering */ maxbandwidth2 = 0; avebandwidth2 = 0; for (i = 0; i < ne; i++) { int off = offset[i]; int curid = newmap[i]; for (j = 0; j < nadj[i]; j++) { int adjid = newmap[adj[off+j]]; int diff = abs(adjid-curid); maxbandwidth2 = (diff > maxbandwidth2) ? diff : maxbandwidth2; avebandwidth2 += diff; nconn++; } } nconn = offset[ne-1]+nadj[ne-1]; avebandwidth2 /= nconn; if (maxbandwidth2 < maxbandwidth1 && avebandwidth2 < avebandwidth1) { /* Renumber */ idx = 0; i = 0; while ((me = MESH_Next_Edge(mesh,&idx))) { ME_Set_ID(me,newmap[i]+1); i++; } fprintf(stderr, "Ave edge ID difference on elements after renumbering: %-lf\n", avebandwidth2); fprintf(stderr, "Max edge ID difference on elements after renumbering: %-d\n", maxbandwidth2); } else { ne = 0; idx = 0; while ((me = MESH_Next_Edge(mesh,&idx))) ME_Set_ID(me,++ne); fprintf(stderr,"Bandwidth did not improve. Keeping old numbering with gaps eliminated\n"); } fprintf(stderr,"\n\n\n"); free(nadj); free(adj); free(offset); free(newmap); } /* Reorder faces according to a breadth first algorithm applied to edges (differs from RCM in that it does not add adjacent graph nodes in ascending order of their valence) */ if (mtype == MFACE || mtype == MALLTYPE) { int nf = MESH_Num_Faces(mesh); if (mtype == MALLTYPE) { /* RCM algorithm already applied on the vertices. Use an edge connected to the starting vertex as the first edge */ List_ptr vfaces = MV_Faces(v0); f0 = List_Entry(vfaces,0); List_Delete(vfaces); } else { /* Find the face whose mid point is a minimum point */ minx = miny = minz = 1.0e+12; f0 = NULL; idx = 0; while ((mf = MESH_Next_Face(mesh,&idx))) { double fxyz[MAXPV2][3]; int nfv; MF_Coords(mf,&nfv,fxyz); xyz[0] = fxyz[0][0]; xyz[1] = fxyz[0][1]; xyz[2] = fxyz[0][2]; for (i = 1; i < nfv; i++) { xyz[0] += fxyz[i][0]; xyz[1] += fxyz[i][1]; xyz[2] += fxyz[i][2]; } xyz[0] /= nfv; xyz[1] /= nfv; xyz[2] /= nfv; if (xyz[0] < minx && xyz[1] < miny && xyz[2] < minz) { minx = xyz[0]; miny = xyz[1]; minz = xyz[2]; f0 = mf; } } } nadj = (int *) malloc(nf*sizeof(int)); nalloc = nf*5; adj = (int *) malloc(nalloc*sizeof(int)); if (!MESH_Num_Regions(mesh)) { int nentries = 0; i = 0; idx = 0; while ((mf = MESH_Next_Face(mesh,&idx))) { List_ptr vfaces, fverts, adjflist; MFace_ptr vf, adjf; MVertex_ptr fv; adjflist = List_New(0); fverts = MF_Vertices(mf,1,0); idx2 = 0; while ((fv = List_Next_Entry(fverts,&idx2))) { List_ptr vfaces = MV_Faces(fv); idx3 = 0; while ((adjf = List_Next_Entry(vfaces,&idx3))) { if (adjf != mf) { int fmarked; #ifdef MSTK_USE_MARKERS fmarked = MEnt_IsMarked(adjf,mkid); #else MEnt_Get_AttVal(adjf, mkatt, &fmarked, &rval, &pval); #endif if (fmarked) { List_Add(adjflist,adjf); #ifdef MSTK_USE_MARKERS MEnt_Mark(adjf,mkid); #else MEnt_Set_AttVal(adjf, mkatt, 1, 0.0, NULL); #endif } } } List_Delete(vfaces); } List_Delete(fverts); #ifdef MSTK_USE_MARKERS List_Unmark(adjflist,mkid); #endif nadj[i] = List_Num_Entries(adjflist); if (nentries+nadj[i] > nalloc) { nalloc *= 2; adj = (int *) realloc(adj,nalloc*sizeof(int)); } idx2 = 0; while ((adjf = List_Next_Entry(adjflist,&idx2))) adj[nentries++] = MF_ID(adjf)-1; List_Delete(adjflist); i++; } } else { int nentries = 0; i = 0; idx = 0; while ((mf = MESH_Next_Face(mesh,&idx))) { List_ptr fregions, adjflist; MRegion_ptr fr; MFace_ptr adjf; adjflist = List_New(0); fregions = MF_Regions(mf); idx2 = 0; while ((fr = List_Next_Entry(fregions,&idx2))) { List_ptr rfaces = MR_Faces(fr); idx3 = 0; while ((adjf = List_Next_Entry(rfaces,&idx3))) { if (adjf != mf) { int fmarked; #ifdef MSTK_USE_MARKERS fmarked = MEnt_IsMarked(adjf,mkid); #else MEnt_Get_AttVal(adjf, mkatt, &fmarked, &rval, &pval); #endif if (fmarked) { List_Add(adjflist,adjf); #ifdef MSTK_USE_MARKERS MEnt_Mark(adjf,mkid); #else MEnt_Set_AttVal(adjf, mkatt, 1, 0.0, NULL); #endif } } } List_Delete(rfaces); } List_Delete(fregions); #ifdef MSTK_USE_MARKERS List_Unmark(adjflist,mkid); #endif nadj[i] = List_Num_Entries(adjflist); if (nentries+nadj[i] > nalloc) { nalloc *= 2; adj = (int *) realloc(adj,nalloc*sizeof(int)); } idx2 = 0; while ((adjf = List_Next_Entry(adjflist,&idx2))) adj[nentries++] = MF_ID(adjf)-1; List_Delete(adjflist); i++; } } /* Compute offsets into adj array */ offset = (int *) malloc(nf*sizeof(int)); offset[0] = 0; for (i = 1; i < nf; i++) offset[i] = offset[i-1] + nadj[i-1]; /* Compute maximum bandwidth before renumbering */ maxbandwidth1 = 0; avebandwidth1 = 0; for (i = 0; i < nf; i++) { int off = offset[i]; int curid = i; for (j = 0; j < nadj[i]; j++) { int adjid = adj[off+j]; int diff = abs(adjid-curid); maxbandwidth1 = (diff > maxbandwidth1) ? diff : maxbandwidth1; avebandwidth1 += diff; nconn++; } } nconn = offset[nf-1]+nadj[nf-1]; avebandwidth1 /= nconn; if (MESH_Num_Regions(mesh)) { fprintf(stderr, "Ave face ID difference on elements before renumbering: %-lf\n", avebandwidth1); fprintf(stderr, "Max face ID difference on elements before renumbering: %-d\n", maxbandwidth1); } else { fprintf(stderr, "Ave face ID difference before renumbering: %-lf\n", avebandwidth1); fprintf(stderr, "Max face ID difference before renumbering: %-d\n", maxbandwidth1); } fprintf(stderr,"\n"); /* Call Graph Renumbering algorithm */ newmap = (int *) malloc(nf*sizeof(int)); Graph_Renumber_GPS(nf, MF_ID(f0)-1, nadj, adj, newmap, &depth, &maxwidth); /* Compute bandwidth after renumbering */ maxbandwidth2 = 0; avebandwidth2 = 0; for (i = 0; i < nf; i++) { int off = offset[i]; int curid = newmap[i]; for (j = 0; j < nadj[i]; j++) { int adjid = newmap[adj[off+j]]; int diff = abs(adjid-curid); maxbandwidth2 = (diff > maxbandwidth2) ? diff : maxbandwidth2; avebandwidth2 += diff; nconn++; } } nconn = offset[nf-1]+nadj[nf-1]; avebandwidth2 /= nconn; if (maxbandwidth2 < maxbandwidth1 && avebandwidth2 < avebandwidth1) { /* Renumber */ idx = 0; i = 0; while ((mf = MESH_Next_Face(mesh,&idx))) { MF_Set_ID(mf,newmap[i]+1); i++; } if (MESH_Num_Regions(mesh)) { fprintf(stderr, "Ave face ID difference on elements after renumbering: %-lf\n", avebandwidth2); fprintf(stderr, "Max face ID difference on elements after renumbering: %-d\n", maxbandwidth2); } else { fprintf(stderr, "Ave face ID difference after renumbering: %-lf\n", avebandwidth2); fprintf(stderr, "Max face ID difference after renumbering: %-d\n", maxbandwidth2); } } else { nf = 0; idx = 0; while ((mf = MESH_Next_Face(mesh,&idx))) MF_Set_ID(mf,++nf); fprintf(stderr,"Bandwidth did not improve. Keeping old numbering with gaps eliminated\n"); } fprintf(stderr,"\n\n\n"); free(nadj); free(adj); free(offset); free(newmap); } if (mtype == MREGION || mtype == MALLTYPE) { int nr = MESH_Num_Regions(mesh); if (nr) { if (mtype == MALLTYPE) { /* Renumbering algorithm already applied on the vertices. Use a region connected to the starting vertex as the first region */ List_ptr vregions = MV_Regions(v0); r0 = List_Entry(vregions,0); List_Delete(vregions); } else { /* Find the region whose center point is a minimum point */ minx = miny = minz = 1.0e+12; r0 = NULL; idx = 0; while ((mr = MESH_Next_Region(mesh,&idx))) { double rxyz[MAXPV3][3]; int nrv; MR_Coords(mr,&nrv,rxyz); xyz[0] = rxyz[0][0]; xyz[1] = rxyz[0][1]; xyz[2] = rxyz[0][2]; for (i = 1; i < nrv; i++) { xyz[0] += rxyz[i][0]; xyz[1] += rxyz[i][1]; xyz[2] += rxyz[i][2]; } xyz[0] /= nrv; xyz[1] /= nrv; xyz[2] /= nrv; if (xyz[0] < minx && xyz[1] < miny && xyz[2] < minz) { minx = xyz[0]; miny = xyz[1]; minz = xyz[2]; r0 = mr; } } } nadj = (int *) malloc(nr*sizeof(int)); nalloc = nr*5; adj = (int *) malloc(nalloc*sizeof(int)); int nentries = 0; i = 0; idx = 0; while ((mr = MESH_Next_Region(mesh,&idx))) { List_ptr vregions, rverts, adjrlist; MRegion_ptr vr, adjr; MVertex_ptr rv; adjrlist = List_New(0); rverts = MR_Vertices(mr); idx2 = 0; while ((rv = List_Next_Entry(rverts,&idx2))) { List_ptr vregions = MV_Regions(rv); idx3 = 0; while ((adjr = List_Next_Entry(vregions,&idx3))) { if (adjr != mr) { int rmarked; #ifdef MSTK_USE_MARKERS rmarked = MEnt_IsMarked(adjr,mkid); #else MEnt_Get_AttVal(adjr, mkatt, &rmarked, &rval, &pval); #endif List_Add(adjrlist,adjr); #ifdef MSTK_USE_MARKERS MEnt_Mark(adjr,mkid); #else MEnt_Set_AttVal(adjr, mkatt, 1, 0.0, NULL); #endif } } List_Delete(vregions); } List_Delete(rverts); #ifdef MSTK_USE_MARKERS List_Unmark(adjrlist,mkid); #endif nadj[i] = List_Num_Entries(adjrlist); if (nentries+nadj[i] > nalloc) { nalloc *= 2; adj = (int *) realloc(adj,nalloc*sizeof(int)); } idx2 = 0; while ((adjr = List_Next_Entry(adjrlist,&idx2))) adj[nentries++] = MR_ID(adjr)-1; List_Delete(adjrlist); i++; } /* Compute offsets into adj array */ offset = (int *) malloc(nr*sizeof(int)); offset[0] = 0; for (i = 1; i < nr; i++) offset[i] = offset[i-1] + nadj[i-1]; /* Compute maximum bandwidth before renumbering */ maxbandwidth1 = 0; avebandwidth1 = 0; for (i = 0; i < nr; i++) { int off = offset[i]; int curid = i; for (j = 0; j < nadj[i]; j++) { int adjid = adj[off+j]; int diff = abs(adjid-curid); maxbandwidth1 = (diff > maxbandwidth1) ? diff : maxbandwidth1; avebandwidth1 += diff; nconn++; } } nconn = offset[nr-1]+nadj[nr-1]; avebandwidth1 /= nconn; fprintf(stderr, "Ave region ID difference before renumbering: %-lf\n", avebandwidth1); fprintf(stderr, "Max region ID difference before renumbering: %-d\n", maxbandwidth1); fprintf(stderr,"\n"); /* Call Graph Renumbering algorithm */ newmap = (int *) malloc(nr*sizeof(int)); Graph_Renumber_GPS(nr, MR_ID(r0)-1, nadj, adj, newmap, &depth, &maxwidth); /* Compute bandwidth after renumbering */ maxbandwidth2 = 0; avebandwidth2 = 0; for (i = 0; i < nr; i++) { int off = offset[i]; int curid = newmap[i]; for (j = 0; j < nadj[i]; j++) { int adjid = newmap[adj[off+j]]; int diff = abs(adjid-curid); maxbandwidth2 = (diff > maxbandwidth2) ? diff : maxbandwidth2; avebandwidth2 += diff; nconn++; } } nconn = offset[nr-1]+nadj[nr-1]; avebandwidth2 /= nconn; if (maxbandwidth2 < maxbandwidth1 && avebandwidth2 < avebandwidth1) { /* Renumber */ idx = 0; i = 0; while ((mr = MESH_Next_Region(mesh,&idx))) { MR_Set_ID(mr,newmap[i]+1); i++; } fprintf(stderr, "Ave region ID difference after renumbering: %-lf\n", avebandwidth2); fprintf(stderr, "Max region ID difference after renumbering: %-d\n", maxbandwidth2); } else { nr = 0; idx = 0; while ((mr = MESH_Next_Region(mesh,&idx))) MR_Set_ID(mr,++nr); fprintf(stderr,"Bandwidth did not improve. Keeping old numbering with gaps eliminated\n"); } fprintf(stderr,"\n\n\n"); free(nadj); free(adj); free(offset); free(newmap); } } #ifdef MSTK_USE_MARKERS MSTK_FreeMarker(mkid); #endif } vidatt = MAttrib_New(mesh,"vidrcm",INT,MVERTEX); idx = 0; while ((mv = MESH_Next_Vertex(mesh,&idx))) { MEnt_Set_AttVal(mv,vidatt,MV_ID(mv),0.0,NULL); } /* We have to reset the max IDs stored in the mesh so that we can correctly assign IDs to new entities */ MESH_Reset_Cached_MaxIDs(mesh); return; }
int MESH_CheckTopo(Mesh_ptr mesh) { int valid = 1; char mesg[256], funcname[32] = "MESH_CheckTopo"; int idx1, idx2, idx3, idx4; MVertex_ptr mv; MEdge_ptr me, ve, fe, re; MFace_ptr mf, vf, ef, rf; MRegion_ptr mr, vr, er, fr; int found, done; int dir; int i, j, k; int nfe; int vid, eid, fid, rid; int gvid, geid, gfid, grid; int gvdim, gedim, gfdim, grdim; int maxiter = 1000; List_ptr vedges, vfaces, vregions; List_ptr efaces; List_ptr fverts, fedges, fregs, fregs1; List_ptr rverts, redges, rfaces; /*****************************************************************/ /* Vertices */ /*****************************************************************/ /* Check that edges connected to vertices reference the vertices */ /* Also check that the classification of the vertex is consistent with respect to the edge */ int first_unknown_classfn = 1; idx1 = 0; while ((mv = MESH_Next_Vertex(mesh,&idx1))) { #ifdef MSTK_HAVE_MPI if (MV_PType(mv) == PGHOST) continue; #endif vid = MV_ID(mv); gvdim = MV_GEntDim(mv); gvid = MV_GEntID(mv); if (gvdim == 4 && first_unknown_classfn) { sprintf(mesg, "Vertex %-d - classification unknown\n", vid); MSTK_Report(funcname, mesg, MSTK_WARN); first_unknown_classfn = 0; } vedges = MV_Edges(mv); if (!vedges) { sprintf(mesg,"Vertex %-d does not have any connected edges\n",vid); MSTK_Report(funcname,mesg,MSTK_WARN); continue; } idx2 = 0; while ((ve = List_Next_Entry(vedges,&idx2))) { eid = ME_ID(ve); if (ME_Vertex(ve,0) != mv && ME_Vertex(ve,1) != mv) { sprintf(mesg,"Vertex %-d connected to edge %-d but edge does not use vertex",vid,eid); MSTK_Report(funcname,mesg,MSTK_ERROR); valid = 0; } } if (gvdim == 1) { /* If vertex is classified on a model edge, then it should be connected to two and only two edges that are classified on the same model edge */ int ne = 0; idx2 = 0; while ((ve = List_Next_Entry(vedges,&idx2))) { gedim = ME_GEntDim(ve); geid = ME_GEntID(ve); if (gedim == 1 && geid == gvid) ne++; } if (ne != 2) { sprintf(mesg,"Vertex %-d classified on model edge %-d but it is not \n connected to two edges classified on this model edge",vid,gvid); MSTK_Report(funcname,mesg,MSTK_WARN); } } List_Delete(vedges); if (gvdim == 2) { MEdge_ptr e0, ecur, enxt; MFace_ptr fcur; int flipped = 0; /* If vertex is classified on a model face, then we should be able to find a ring of faces classified on that model face */ vfaces = MV_Faces(mv); found = 0; idx2 = 0; while ((vf = List_Next_Entry(vfaces,&idx2))) { if (MF_GEntDim(vf) == 2) { found = 1; break; } } List_Delete(vfaces); if (!found) { sprintf(mesg,"Vertex %-d classified on model face %-d but could not \n find connected face classified on this model face",vid,gvid); MSTK_Report(funcname,mesg,MSTK_WARN); valid = 0; } fcur = vf; fedges = MF_Edges(fcur,1,mv); nfe = List_Num_Entries(fedges); e0 = List_Entry(fedges,0); ecur = e0; enxt = List_Entry(fedges,nfe-1); List_Delete(fedges); done = 0; i = 0; while (!done) { ecur = enxt; efaces = ME_Faces(ecur); found = 0; idx3 = 0; while ((ef = List_Next_Entry(efaces,&idx3))) { if (ef != fcur && MF_GEntDim(ef) == 2 && MF_GEntID(ef) == gvid) { fcur = ef; found = 1; break; } } List_Delete(efaces); if (!found) { sprintf(mesg,"Could not find next boundary face connected to vertex %-d",vid); MSTK_Report(funcname,mesg,MSTK_WARN); valid = 0; break; } fedges = MF_Edges(fcur,1,mv); nfe = List_Num_Entries(fedges); if (List_Entry(fedges,0) == ecur) enxt = List_Entry(fedges,nfe-1); else if (List_Entry(fedges,nfe-1) == ecur) { enxt = List_Entry(fedges,0); flipped = 1; } else { sprintf(mesg,"Could not find next edge while traversing around vertex %-d on model face %-d",vid,gvid); MSTK_Report(funcname,mesg,MSTK_ERROR); } List_Delete(fedges); if (enxt == e0) done = 1; if (++i > maxiter) break; } if (!done) { sprintf(mesg,"Vertex %-d classified on model face %-d but could not find ring of faces classified on this model face",vid,gvid); MSTK_Report(funcname,mesg,MSTK_WARN); } if (done && flipped) { List_ptr fregs = MF_Regions(fcur); if (List_Num_Entries(fregs) < 2) { sprintf(mesg,"Inconsistent orientations of boundary faces around vertex %-d",vid); MSTK_Report(funcname,mesg,MSTK_WARN); } if (fregs) List_Delete(fregs); } } } /* while ((mv = MESH_Next_Vertex(mesh,&idx1))) */ /*****************************************************************/ /* Edges */ /*****************************************************************/ first_unknown_classfn = 1; idx1 = 0; while ((me = MESH_Next_Edge(mesh,&idx1))) { #ifdef MSTK_HAVE_MPI if (ME_PType(me) == PGHOST) continue; #endif eid = ME_ID(me); gedim = ME_GEntDim(me); geid = ME_GEntID(me); if (gedim == 4 && first_unknown_classfn) { sprintf(mesg, "Edge %-d - unknown classification", eid); MSTK_Report(funcname, mesg, MSTK_WARN); first_unknown_classfn = 0; } if (ME_Vertex(me,0) == ME_Vertex(me,1)) { sprintf(mesg,"Edge %-d has repeated vertices",eid); MSTK_Report(funcname,mesg,MSTK_ERROR); } for (i = 0; i < 2; i++) { MVertex_ptr ev = ME_Vertex(me,i); vid = MV_ID(ev); gvid = MV_GEntID(ev); gvdim = MV_GEntDim(ev); if (gvdim != 4 && gvdim != 4) { /* vertex and edge classifn is known */ if (gedim < gvdim) { sprintf(mesg,"Edge %-d classified on lower dimensional entity than connected vertex %-d",eid,vid); MSTK_Report(funcname,mesg,MSTK_WARN); valid = 0; } else if (gedim == gvdim && geid != gvid) { sprintf(mesg,"Edge %-d and its vertex %-d classified on different entities of the same dimension",eid,vid); MSTK_Report(funcname,mesg,MSTK_WARN); valid = 0; } } vedges = MV_Edges(ev); if (!List_Contains(vedges,me)) { sprintf(mesg,"Edge %-d sees vertex %-d but not vice versa",eid,vid); MSTK_Report(funcname,mesg,MSTK_ERROR); valid = 0; } List_Delete(vedges); if (gedim == 2) { MFace_ptr ebf[2], fcur, fnxt; MRegion_ptr rcur; int nf, nfr; List_ptr eregs; /* Edge is classified on model face - it should be connected to two and only two faces also classified on this model face */ ebf[0] = ebf[1] = NULL; nf = 0; efaces = ME_Faces(me); idx2 = 0; while ((ef = List_Next_Entry(efaces,&idx2))) { fid = MF_ID(ef); if (MF_GEntDim(ef) == 2) { nf++; if (gedim == 2 && MF_GEntID(ef) != geid) { sprintf(mesg,"Face %-d connected to edge %-d classified on different model face",fid,eid); MSTK_Report(funcname,mesg,MSTK_WARN); valid = 0; } if (ebf[0] == NULL) ebf[0] = ef; else ebf[1] = ef; } } List_Delete(efaces); if (nf != 2) { sprintf(mesg,"Boundary edge %-d is not connected to exactly two\n faces classified on the boundary",eid); MSTK_Report(funcname,mesg,MSTK_ERROR); valid = 0; } eregs = ME_Regions(me); if (!eregs) continue; else List_Delete(eregs); /* Can we go from f0 to f1 in one or two dirs? */ fcur = ebf[0]; fnxt = NULL; fregs = MF_Regions(fcur); if (!fregs) { fid = MF_ID(fcur); sprintf(mesg,"Edge %-d connected to regions but face %-d is not",eid,fid); MSTK_Report(funcname,mesg,MSTK_ERROR); valid = 0; } nfr = List_Num_Entries(fregs); for (i = 0; i < nfr; i++) { rcur = List_Entry(fregs,i); rfaces = MR_Faces(rcur); idx3 = 0; found = 0; while ((rf = List_Next_Entry(rfaces,&idx3))) { if (rf != fcur && MF_UsesEntity(rf,me,1)) { found = 1; fnxt = rf; break; } } List_Delete(rfaces); if (!found) { rid = MR_ID(rcur); sprintf(mesg,"Could not find second face in region %-d using edge %-d",rid,eid); } done = 0; j = 0; while (!done) { fcur = fnxt; fid = MF_ID(fcur); if (fnxt == ebf[1]) { done = 1; break; } fregs1 = MF_Regions(fcur); idx3 = 0; while ((fr = List_Next_Entry(fregs1,&idx3))) { if (fr != rcur) { rcur = fr; found = 1; break; } } List_Delete(fregs1); if (!found) { sprintf(mesg,"Could not find next region around edge %-d",eid); MSTK_Report(funcname,mesg,MSTK_ERROR); valid = 0; break; } rfaces = MR_Faces(rcur); idx3 = 0; found = 0; while ((rf = List_Next_Entry(rfaces,&idx3))) { if (rf != fcur && MF_UsesEntity(rf,me,1)) { found = 1; fnxt = rf; break; } } List_Delete(rfaces); if (!found) { rid = MR_ID(rcur); sprintf(mesg,"Could not find second face in region %-d using edge %-d",rid,eid); } if (++j > maxiter) break; } /* while (!done) */ if (!done) { sprintf(mesg,"Could not traverse around edge %-d from face %-d to face %-d",eid,MF_ID(ebf[0]),MF_ID(ebf[1])); MSTK_Report(funcname,mesg,MSTK_ERROR); valid = 0; } } /* for (i = 0; i < nfr; i++) */ List_Delete(fregs); } /* if (geid == 2) */ } /* for (i = 0; i < 2; i++) */ } /* while ((me = MESH_Next_Edge(mesh,&idx1))) */ /*****************************************************************/ /* Faces */ /*****************************************************************/ first_unknown_classfn = 1; idx1 = 0; while ((mf = MESH_Next_Face(mesh,&idx1))) { #ifdef MSTK_HAVE_MPI if (MF_PType(mf) == PGHOST) continue; #endif fid = MF_ID(mf); gfid = MF_GEntID(mf); gfdim = MF_GEntDim(mf); if (gfdim == 4 && first_unknown_classfn) { sprintf(mesg, "Face %-d - unknown classification", fid); MSTK_Report(funcname, mesg, MSTK_WARN); first_unknown_classfn = 0; } fedges = MF_Edges(mf,1,0); if (List_Num_Entries(fedges) < 3) { sprintf(mesg,"Face %-d has less than 3 edges",fid); MSTK_Report(funcname,mesg,MSTK_ERROR); } idx2 = 0; while ((fe = List_Next_Entry(fedges,&idx2))) { eid = ME_ID(fe); geid = ME_GEntID(fe); gedim = ME_GEntDim(fe); if (gedim != 4 && gfdim != 4) { /* Edge, Face classfn is known */ if (gfdim < gedim) { sprintf(mesg,"Face %-d classified on lower order entity than edge %-d",fid,ME_ID(fe)); MSTK_Report(funcname,mesg,MSTK_WARN); valid = 0; } else if (gedim == gfdim && geid != gfid) { sprintf(mesg,"Face %-d and edge %-d classified on different\n entities of the same dimension",fid,eid); MSTK_Report(funcname,mesg,MSTK_WARN); } } efaces = ME_Faces(fe); if (!List_Contains(efaces,mf)) { sprintf(mesg,"Face %-d refers to edge %-d but not vice versa",fid,ME_ID(fe)); MSTK_Report(funcname,mesg,MSTK_ERROR); valid = 0; } List_Delete(efaces); } List_Delete(fedges); fregs = MF_Regions(mf); if (gfdim == 3) { if (!fregs || List_Num_Entries(fregs) != 2) { sprintf(mesg,"Interior face %-d does not have two connected regions",fid); MSTK_Report(funcname,mesg,MSTK_ERROR); valid = 0; } } if (fregs) { if (List_Num_Entries(fregs) == 2) { if (MR_FaceDir(List_Entry(fregs,0),mf) == MR_FaceDir(List_Entry(fregs,1),mf)) { sprintf(mesg,"Both regions using face %-d in the same sense",fid); MSTK_Report(funcname,mesg,MSTK_ERROR); valid = 0; } } List_Delete(fregs); } } /* while ((mf = MESH_Next_Face(mesh,&idx1))) */ /*****************************************************************/ /* Regions */ /*****************************************************************/ idx1 = 0; while ((mr = MESH_Next_Region(mesh,&idx1))) { #ifdef MSTK_HAVE_MPI if (MR_PType(mr) == PGHOST) continue; #endif rid = MR_ID(mr); grid = MR_GEntID(mr); rfaces = MR_Faces(mr); int nrf = List_Num_Entries(rfaces); if (nrf < 4) { sprintf(mesg,"Region %-d has less than 4 faces",rid); MSTK_Report(funcname,mesg,MSTK_ERROR); } /* Check that face to region and region to face links are consistent with each other */ int *rfdirs = (int *) malloc(nrf*sizeof(int)); i = 0; idx2 = 0; while ((rf = List_Next_Entry(rfaces,&idx2))) { rfdirs[i] = MR_FaceDir_i(mr,i); if (mr != MF_Region(rf,!rfdirs[i])) { sprintf(mesg,"Region %-d to face %-d dir inconsistent with \n face to region dir",rid,MF_ID(rf)); MSTK_Report(funcname,mesg,MSTK_ERROR); valid = 0; } i++; } /* Check that faces of a region have consistent orientation in the region with respect to each other */ for (i = 0; i < nrf; i++) { MFace_ptr rf, rf2; rf = List_Entry(rfaces,i); fedges = MF_Edges(rf,1,0); nfe = List_Num_Entries(fedges); for (j = 0; j < nfe; j++) { fe = List_Entry(fedges,j); int fedir = MF_EdgeDir_i(rf,j); /* Find adjacent face in the region */ found = 0; for (k = 0; k < nrf; k++) { rf2 = List_Entry(rfaces,k); if (rf != rf2 && MF_UsesEntity(rf2,fe,MEDGE)) { found = 1; break; } } if (!found) { sprintf(mesg,"Cannot find another face in region %-d sharing edge %-d (ID = %-d) of face with ID = %-d",MR_ID(mr),j,ME_ID(fe),MF_ID(rf)); MSTK_Report(funcname,mesg,MSTK_ERROR); valid = 0; } int fedir_adj = MF_EdgeDir(rf2,fe); /* If the two faces use the edge in opposite directions then the region should use the faces in the same direction and vice versa */ if (((fedir_adj == fedir) && (rfdirs[i] == rfdirs[k])) || ((fedir_adj != fedir) && (rfdirs[i] != rfdirs[k]))) { sprintf(mesg,"Region %-d faces are inconsistently oriented",MR_ID(mr)); MSTK_Report(funcname,mesg,MSTK_ERROR); valid = 0; } } List_Delete(fedges); } List_Delete(rfaces); free(rfdirs); } return valid; } /* int MESH_CheckTopo */
int main(int argc, char **argv) { int len, ok; int i, j, id, idx, ncells, file_found; int num_parts, imethod; int nproc, rank; Mesh_ptr mesh; char basename[256], meshname[256], gmvfilename[256], method[256]; if (argc == 1) { fprintf(stderr,"Usage: %s --num=num_parts --method=Metis/Zoltan (default=Metis) --file=name.mstk\n",argv[0]); exit(-1); } file_found = 0; num_parts = 4; strcpy(method,"Metis"); for (i = 1; i < argc; i++) { if (strncmp(argv[i],"--file=",7) == 0) { sscanf(argv[i]+7,"%s",meshname); file_found = 1; } else if (strncmp(argv[i],"--num=",6) == 0) sscanf(argv[i]+6,"%d",&num_parts); else if (strncmp(argv[i],"--method=",9) == 0) sscanf(argv[i]+9,"%s",method); else if (strncmp(argv[i],"--help",6) == 0) { fprintf(stderr,"Usage: %s --num=num_parts --method=Metis/Zoltan (default=Metis) --file=name.mstk\n",argv[0]); exit(-1); } } if (strncasecmp(method,"metis",5) == 0) imethod = 0; else if (strncasecmp(method,"zoltan",6) == 0) imethod = 1; else { fprintf(stderr,"vizpart: Partitioning method not recognized\n"); exit(-1); } if (!file_found) { fprintf(stderr,"Must specify input filename using --file argument\n"); exit(-1); } MPI_Init(&argc,&argv); MSTK_Init(); MSTK_Comm comm = MPI_COMM_WORLD; MPI_Comm_size(comm,&nproc); MPI_Comm_rank(comm,&rank); if (imethod == 1 && nproc != num_parts) { fprintf(stderr,"Zoltan based partitioner: Number of processors must equal requested number of partition\n"); exit(-1); } strcpy(basename,meshname); len = strlen(meshname); if (len > 5 && strncmp(&(meshname[len-5]),".mstk",5) == 0) basename[len-5] = '\0'; else strcat(meshname,".mstk"); mesh = MESH_New(UNKNOWN_REP); if (rank == 0) { ok = MESH_InitFromFile(mesh,meshname,comm); if (!ok) { fprintf(stderr,"Cannot file input file %s\n\n\n",meshname); exit(-1); } if ( (ncells = MESH_Num_Regions(mesh)) > 0) { printf("3D mesh with %d regions\n", ncells); } else if ( (ncells = MESH_Num_Faces(mesh)) > 0) printf("2D mesh with %d faces\n", ncells); else { fprintf(stderr,"Mesh is neither solid nor surface mesh. Exiting...\n"); exit(-1); } } Mesh_ptr *submeshes = (Mesh_ptr*) malloc(num_parts*sizeof(Mesh_ptr)); int *part; MESH_Get_Partitioning(mesh, imethod, &part, comm); if(rank == 0) { int del_inmesh = 0; MESH_Partition(mesh, num_parts, part, submeshes); idx = 0; if(MESH_Num_Regions(mesh)) { MAttrib_ptr region_part = MAttrib_New(mesh, "part", INT, MREGION); MRegion_ptr mr; while ((mr = MESH_Next_Region(mesh, &idx))) { id = MR_ID(mr); MEnt_Set_AttVal(mr,region_part,part[id-1],0,NULL); } } else { MAttrib_ptr face_part = MAttrib_New(mesh, "part", INT, MFACE); MFace_ptr mf; while ((mf = MESH_Next_Face(mesh, &idx))) { id = MF_ID(mf); MEnt_Set_AttVal(mf,face_part,part[idx-1],0,NULL); } } sprintf(gmvfilename,"%s_part.gmv",basename); MESH_ExportToGMV(mesh,gmvfilename,0,NULL,NULL,comm); for( i = 0; i < num_parts; i++) { sprintf(gmvfilename,"%s_part.%04d.%04d.gmv",basename,i,0); MESH_ExportToGMV(submeshes[i],gmvfilename,0,NULL,NULL,comm); } for (i = 0; i < num_parts; i++) MESH_Delete(submeshes[i]); } free(submeshes); MESH_Delete(mesh); free(part); MPI_Finalize(); return 1; }
int MESH_PartitionWithZoltan(Mesh_ptr mesh, int nparts, int **part, int noptions, char **options, MSTK_Comm comm) { MEdge_ptr fedge; MFace_ptr mf, oppf, rface; MRegion_ptr mr, oppr; List_ptr fedges, efaces, rfaces, fregions; int i, j, k, id; int nv, ne, nf, nr=0, nfe, nef, nfr, nrf, idx, idx2; int numflag, nedgecut, ipos; int wtflag; int rc; float ver; struct Zoltan_Struct *zz; GRAPH_DATA graph; int changes, numGidEntries, numLidEntries, numImport, numExport; ZOLTAN_ID_PTR importGlobalGids, importLocalGids, exportGlobalGids, exportLocalGids; int *importProcs, *importToPart, *exportProcs, *exportToPart; int rank; MPI_Comm_rank(comm,&rank); rc = Zoltan_Initialize(0, NULL, &ver); if (rc != ZOLTAN_OK){ MSTK_Report("MESH_PartitionWithZoltan","Could not initialize Zoltan",MSTK_FATAL); MPI_Finalize(); exit(0); } /****************************************************************** ** Create a Zoltan library structure for this instance of partition ********************************************************************/ zz = Zoltan_Create(comm); /***************************************************************** ** Figure out partitioning method *****************************************************************/ char partition_method_str[32]; strcpy(partition_method_str,"RCB"); /* Default - Recursive Coordinate Bisection */ if (noptions) { for (i = 0; i < noptions; i++) { if (strncmp(options[i],"LB_PARTITION",12) == 0) { char *result = NULL, instring[256]; strcpy(instring,options[i]); result = strtok(instring,"="); result = strtok(NULL," "); strcpy(partition_method_str,result); } } } if (rank == 0) { char mesg[256]; sprintf(mesg,"Using partitioning method %s for ZOLTAN\n",partition_method_str); MSTK_Report("MESH_PartitionWithZoltan",mesg,MSTK_MESG); } /* General parameters for Zoltan */ Zoltan_Set_Param(zz, "DEBUG_LEVEL", "0"); Zoltan_Set_Param(zz, "LB_METHOD", partition_method_str); Zoltan_Set_Param(zz, "LB_APPROACH", "PARTITION"); Zoltan_Set_Param(zz, "NUM_GID_ENTRIES", "1"); Zoltan_Set_Param(zz, "NUM_LID_ENTRIES", "1"); Zoltan_Set_Param(zz, "RETURN_LISTS", "ALL"); graph.numMyNodes = 0; graph.numAllNbors = 0; graph.nodeGID = NULL; graph.nodeCoords = NULL; graph.nborIndex = NULL; graph.nborGID = NULL; graph.nborProc = NULL; if (strcmp(partition_method_str,"RCB") == 0) { if (rank == 0) { nr = MESH_Num_Regions(mesh); nf = MESH_Num_Faces(mesh); if (!nf && !nr) MSTK_Report("MESH_PartitionWithZoltan","Cannot partition wire meshes", MSTK_FATAL); if (nr == 0) { /* Surface or planar mesh */ int ndim = 2; /* assume mesh is planar */ idx = 0; MVertex_ptr mv; while ((mv = MESH_Next_Vertex(mesh,&idx))) { double vxyz[3]; MV_Coords(mv,vxyz); if (vxyz[2] != 0.0) { ndim = 3; /* non-planar or planar with non-zero z */ break; } } NDIM_4_ZOLTAN = ndim-1; /* ignore last dimension to avoid partitioning in that dimension */ graph.numMyNodes = nf; graph.nodeGID = (ZOLTAN_ID_TYPE *) malloc(sizeof(ZOLTAN_ID_TYPE) * nf); graph.nodeCoords = (double *) malloc(sizeof(double) * NDIM_4_ZOLTAN * nf); idx = 0; while ((mf = MESH_Next_Face(mesh,&idx))) { double fxyz[MAXPV2][3], cen[3]; int nfv; MF_Coords(mf,&nfv,fxyz); cen[0] = cen[1] = cen[2] = 0.0; for (j = 0; j < nfv; j++) for (k = 0; k < NDIM_4_ZOLTAN; k++) cen[k] += fxyz[j][k]; for (k = 0; k < NDIM_4_ZOLTAN; k++) cen[k] /= nfv; id = MF_ID(mf); graph.nodeGID[id-1] = id; memcpy(&(graph.nodeCoords[NDIM_4_ZOLTAN*(id-1)]),cen,NDIM_4_ZOLTAN*sizeof(double)); } } else { /* Volume mesh */ int ndim = 3; NDIM_4_ZOLTAN = ndim-1; /* ignore last dimension to avoid partitioning in that dimension */ graph.numMyNodes = nr; graph.nodeGID = (ZOLTAN_ID_TYPE *) malloc(sizeof(ZOLTAN_ID_TYPE) * nr); graph.nodeCoords = (double *) malloc(sizeof(double) * NDIM_4_ZOLTAN * nr); idx = 0; while ((mr = MESH_Next_Region(mesh,&idx))) { double rxyz[MAXPV3][3], cen[3]; int nrv; MR_Coords(mr,&nrv,rxyz); cen[0] = cen[1] = cen[2] = 0.0; for (j = 0; j < nrv; j++) for (k = 0; k < NDIM_4_ZOLTAN; k++) cen[k] += rxyz[j][k]; for (k = 0; k < NDIM_4_ZOLTAN; k++) cen[k] /= nrv; for (k = 0; k < NDIM_4_ZOLTAN; k++) if (fabs(cen[k]) < 1.0e-10) cen[k] = 0.0; id = MR_ID(mr); graph.nodeGID[id-1] = id; memcpy(&(graph.nodeCoords[NDIM_4_ZOLTAN*(id-1)]),cen,NDIM_4_ZOLTAN*sizeof(double)); } } } MPI_Bcast(&NDIM_4_ZOLTAN,1,MPI_INT,0,comm); /* Set some default values */ Zoltan_Set_Param(zz, "RCB_RECTILINEAR_BLOCKS","1"); // Zoltan_Set_Param(zz, "AVERAGE_CUTS", "1"); if (noptions > 1) { for (i = 1; i < noptions; i++) { char *paramstr = NULL, *valuestr = NULL, instring[256]; strcpy(instring,options[i]); paramstr = strtok(instring,"="); valuestr = strtok(NULL," "); Zoltan_Set_Param(zz,paramstr,valuestr); } } /* Query functions - defined in simpleQueries.h */ Zoltan_Set_Num_Obj_Fn(zz, get_number_of_nodes, &graph); Zoltan_Set_Obj_List_Fn(zz, get_node_list, &graph); Zoltan_Set_Num_Geom_Fn(zz, get_num_dimensions_reduced, &graph); /* reduced dimensions */ Zoltan_Set_Geom_Multi_Fn(zz, get_element_centers_reduced, &graph); /* reduced dimension centers */ } else if (strcmp(partition_method_str,"GRAPH") == 0) { if(rank == 0) { nv = MESH_Num_Vertices(mesh); ne = MESH_Num_Edges(mesh); nf = MESH_Num_Faces(mesh); nr = MESH_Num_Regions(mesh); ipos = 0; /* build nodes and neighbors list, similar as in partition with metis Assign processor 0 the whole mesh, assign other processors a NULL mesh */ if (nr == 0) { if (nf == 0) { MSTK_Report("MESH_PartitionWithZoltan", "Cannot partition wire meshes with Zoltan",MSTK_FATAL); exit(-1); } graph.nodeGID = (ZOLTAN_ID_TYPE *)malloc(sizeof(ZOLTAN_ID_TYPE) * nf); graph.nborIndex = (int *)malloc(sizeof(int) * (nf + 1)); graph.nborGID = (ZOLTAN_ID_TYPE *)malloc(sizeof(ZOLTAN_ID_TYPE) * 2*ne); graph.nborProc = (int *)malloc(sizeof(int) * 2*ne); graph.nborIndex[0] = 0; /* Surface mesh */ idx = 0; i = 0; while ((mf = MESH_Next_Face(mesh,&idx))) { graph.nodeGID[i] = MF_ID(mf); fedges = MF_Edges(mf,1,0); nfe = List_Num_Entries(fedges); idx2 = 0; while ((fedge = List_Next_Entry(fedges,&idx2))) { efaces = ME_Faces(fedge); nef = List_Num_Entries(efaces); if (nef == 1) { continue; /* boundary edge; nothing to do */ } else { int j; for (j = 0; j < nef; j++) { oppf = List_Entry(efaces,j); if (oppf == mf) { graph.nborGID[ipos] = MF_ID(oppf); /* initially set all nodes on processor 0 */ graph.nborProc[ipos] = 0; ipos++; } } } List_Delete(efaces); } List_Delete(fedges); i++; graph.nborIndex[i] = ipos; } graph.numMyNodes = i; graph.numAllNbors = ipos; } else { graph.nodeGID = (ZOLTAN_ID_TYPE *)malloc(sizeof(ZOLTAN_ID_TYPE) * nr); graph.nborIndex = (int *)malloc(sizeof(int) * (nr + 1)); graph.nborGID = (ZOLTAN_ID_TYPE *)malloc(sizeof(ZOLTAN_ID_TYPE) * 2*nf); graph.nborProc = (int *)malloc(sizeof(int) * 2*nf); graph.nborIndex[0] = 0; /* Volume mesh */ idx = 0; i = 0; while ((mr = MESH_Next_Region(mesh,&idx))) { graph.nodeGID[i] = MR_ID(mr); rfaces = MR_Faces(mr); nrf = List_Num_Entries(rfaces); idx2 = 0; while ((rface = List_Next_Entry(rfaces,&idx2))) { fregions = MF_Regions(rface); nfr = List_Num_Entries(fregions); if (nfr > 1) { oppr = List_Entry(fregions,0); if (oppr == mr) oppr = List_Entry(fregions,1); graph.nborGID[ipos] = MR_ID(oppr); /* initially set all nodes on processor 0 */ graph.nborProc[ipos] = 0; ipos++; } List_Delete(fregions); } List_Delete(rfaces); i++; graph.nborIndex[i] = ipos; } graph.numMyNodes = i; graph.numAllNbors = ipos; } } /* Graph parameters */ /* Zoltan_Set_Param(zz, "CHECK_GRAPH", "2"); */ Zoltan_Set_Param(zz, "PHG_EDGE_SIZE_THRESHOLD", ".35"); /* 0-remove all, 1-remove none */ /* Query functions - defined in simpleQueries.h */ Zoltan_Set_Num_Obj_Fn(zz, get_number_of_nodes, &graph); Zoltan_Set_Obj_List_Fn(zz, get_node_list, &graph); Zoltan_Set_Num_Edges_Multi_Fn(zz, get_num_edges_list, &graph); Zoltan_Set_Edge_List_Multi_Fn(zz, get_edge_list, &graph); } /* Partition the graph */ /****************************************************************** ** Zoltan can now partition the graph. ** We assume the number of partitions is ** equal to the number of processes. Process rank 0 will own ** partition 0, process rank 1 will own partition 1, and so on. ******************************************************************/ rc = Zoltan_LB_Partition(zz, /* input (all remaining fields are output) */ &changes, /* 1 if partitioning was changed, 0 otherwise */ &numGidEntries, /* Number of integers used for a global ID */ &numLidEntries, /* Number of integers used for a local ID */ &numImport, /* Number of nodes to be sent to me */ &importGlobalGids, /* Global IDs of nodes to be sent to me */ &importLocalGids, /* Local IDs of nodes to be sent to me */ &importProcs, /* Process rank for source of each incoming node */ &importToPart, /* New partition for each incoming node */ &numExport, /* Number of nodes I must send to other processes*/ &exportGlobalGids, /* Global IDs of the nodes I must send */ &exportLocalGids, /* Local IDs of the nodes I must send */ &exportProcs, /* Process to which I send each of the nodes */ &exportToPart); /* Partition to which each node will belong */ if (rc != ZOLTAN_OK){ if (rank == 0) MSTK_Report("MESH_PartitionWithZoltan","Could not partition mesh with ZOLTAN", MSTK_ERROR); Zoltan_Destroy(&zz); MPI_Finalize(); return 0; } if(rank == 0) { *part = (int *) calloc(graph.numMyNodes,sizeof(int)); for ( i = 0; i < numExport; i++ ) { (*part)[exportGlobalGids[i]-1] = exportToPart[i]; } if (graph.nodeGID) free(graph.nodeGID); if (graph.nodeCoords) free(graph.nodeCoords); if (graph.nborIndex) free(graph.nborIndex); if (graph.nborGID) free(graph.nborGID); if (graph.nborProc) free(graph.nborProc); } else { *part = NULL; } Zoltan_LB_Free_Part(&exportGlobalGids, &exportLocalGids, &exportProcs, &exportToPart); Zoltan_LB_Free_Part(&importGlobalGids, &importLocalGids, &importProcs, &importToPart); Zoltan_Destroy(&zz); return 1; }
int MESH_PartitionWithMetis(Mesh_ptr mesh, int nparts, int **part) { MEdge_ptr fedge; MFace_ptr mf, oppf, rface; MRegion_ptr mr, oppr; List_ptr fedges, efaces, rfaces, fregions; int i, ncells, ipos; int nv, ne, nf, nr, nfe, nef, nfr, nrf, idx, idx2; #ifdef METIS_5 idx_t ngraphvtx, numflag, nedgecut, numparts, ncons; idx_t wtflag, metisopts[METIS_NOPTIONS]; idx_t *vsize, *idxpart; idx_t *xadj, *adjncy, *vwgt, *adjwgt; real_t *tpwgts, *ubvec; #else idxtype ngraphvtx, numflag, nedgecut, numparts; idxtype wtflag, metisopts[5] = {0,0,0,0,0}; idxtype *xadj, *adjncy, *vwgt, *adjwgt, *idxpart; #endif /* First build a nodal graph of the mesh in the format required by metis */ nv = MESH_Num_Vertices(mesh); ne = MESH_Num_Edges(mesh); nf = MESH_Num_Faces(mesh); nr = MESH_Num_Regions(mesh); ipos = 0; if (nr == 0) { if (nf == 0) { fprintf(stderr,"Cannot partition wire meshes\n"); exit(-1); } #ifdef METIS_5 xadj = (idx_t *) malloc((nf+1)*sizeof(idx_t)); adjncy = (idx_t *) malloc(2*ne*sizeof(idx_t)); #else xadj = (idxtype *) malloc((nf+1)*sizeof(idxtype)); adjncy = (idxtype *) malloc(2*ne*sizeof(idxtype)); #endif ncells = nf; /* Surface mesh */ idx = 0; i = 0; xadj[i] = ipos; while ((mf = MESH_Next_Face(mesh,&idx))) { fedges = MF_Edges(mf,1,0); nfe = List_Num_Entries(fedges); idx2 = 0; while ((fedge = List_Next_Entry(fedges,&idx2))) { efaces = ME_Faces(fedge); nef = List_Num_Entries(efaces); if (nef == 1) { continue; /* boundary edge; nothing to do */ } else { int j; for (j = 0; j < nef; j++) { oppf = List_Entry(efaces,j); if (oppf != mf) { adjncy[ipos] = MF_ID(oppf)-1; ipos++; } } } List_Delete(efaces); } List_Delete(fedges); i++; xadj[i] = ipos; } } else { #ifdef METIS_5 xadj = (idx_t *) malloc((nr+1)*sizeof(idx_t)); adjncy = (idx_t *) malloc(2*nf*sizeof(idx_t)); #else xadj = (idxtype *) malloc((nr+1)*sizeof(idxtype)); adjncy = (idxtype *) malloc(2*nf*sizeof(idxtype)); #endif ncells = nr; /* Volume mesh */ idx = 0; i = 0; xadj[i] = ipos; while ((mr = MESH_Next_Region(mesh,&idx))) { rfaces = MR_Faces(mr); nrf = List_Num_Entries(rfaces); idx2 = 0; while ((rface = List_Next_Entry(rfaces,&idx2))) { fregions = MF_Regions(rface); nfr = List_Num_Entries(fregions); if (nfr > 1) { oppr = List_Entry(fregions,0); if (oppr == mr) oppr = List_Entry(fregions,1); adjncy[ipos] = MR_ID(oppr)-1; ipos++; } List_Delete(fregions); } List_Delete(rfaces); i++; xadj[i] = ipos; } } /* Partition the graph */ wtflag = 0; /* No weights are specified */ vwgt = adjwgt = NULL; numflag = 0; /* C style numbering of elements (nodes of the dual graph) */ ngraphvtx = ncells; /* we want the variable to be of type idxtype or idx_t */ numparts = nparts; /* we want the variable to be of type idxtype or idx_t */ #ifdef METIS_5 idxpart = (idx_t *) malloc(ncells*sizeof(idx_t)); ncons = 1; /* Number of constraints */ vsize = NULL; tpwgts = NULL; ubvec = NULL; METIS_SetDefaultOptions(metisopts); metisopts[METIS_OPTION_NUMBERING] = 0; if (nparts <= 8) METIS_PartGraphRecursive(&ngraphvtx,&ncons,xadj,adjncy,vwgt,vsize,adjwgt, &numparts,tpwgts,ubvec,metisopts,&nedgecut, idxpart); else METIS_PartGraphKway(&ngraphvtx,&ncons,xadj,adjncy,vwgt,vsize,adjwgt, &numparts,tpwgts,ubvec,metisopts,&nedgecut,idxpart); #else idxpart = (idxtype *) malloc(ncells*sizeof(idxtype)); if (nparts <= 8) METIS_PartGraphRecursive(&ngraphvtx,xadj,adjncy,vwgt,adjwgt,&wtflag, &numflag,&numparts,metisopts,&nedgecut,idxpart); else METIS_PartGraphKway(&ngraphvtx,xadj,adjncy,vwgt,adjwgt,&wtflag,&numflag, &numparts,metisopts,&nedgecut,idxpart); #endif free(xadj); free(adjncy); *part = (int *) malloc(ncells*sizeof(int)); for (i = 0; i < ncells; i++) (*part)[i] = (int) idxpart[i]; free(idxpart); return 1; }
int MESH_Send_NonVertexEntities_FN(Mesh_ptr mesh, int torank, MSTK_Comm comm, int *numreq, int *maxreq, MPI_Request **requests, int *numptrs2free, int *maxptrs2free, void ***ptrs2free) { int i, j, nv, ne, nf, nr; int nevs, nfes, nrfs, nfe, nrv, nrf, dir; int maxnfe, maxnrf; int *mesh_info; int *list_edge=NULL, *list_face=NULL, *list_region=NULL; MVertex_ptr mv; MEdge_ptr me; MFace_ptr mf; MRegion_ptr mr; List_ptr mfedges, mrfaces, mrverts; RepType rtype; double coor[3]; MPI_Request mpirequest; if (requests == NULL) MSTK_Report("MESH_Surf_SendMesh","MPI requests array is NULL",MSTK_FATAL); if (*maxreq == 0) { *maxreq = 25; *requests = (MPI_Request *) malloc(*maxreq*sizeof(MPI_Request)); *numreq = 0; } else if (*maxreq < (*numreq) + 13) { *maxreq = 2*(*maxreq) + 11; *requests = (MPI_Request *) realloc(*requests,*maxreq*sizeof(MPI_Request)); } ne = MESH_Num_Edges(mesh); nf = MESH_Num_Faces(mesh); nr = MESH_Num_Regions(mesh); /* some other known quantitites - 5 items per edge (2 for verts and 3 for extra data), maxnfe+4 items per face (1 for number of edges, maxnfe for edge indices, anad 3 for extra data), maxnrf+4 items per region (1 for number of faces, maxnrf for face indices and 3 for extra data */ maxnfe = 0; for (i = 0; i < nf; i++) { mf = MESH_Face(mesh,i); nfe = MF_Num_Edges(mf); if (nfe > maxnfe) maxnfe = nfe; } maxnrf = 0; for (i = 0; i < nr; i++) { mr = MESH_Region(mesh,i); nrf = MR_Num_Faces(mr); if (nrf > maxnrf) maxnrf = nrf; } // The amount of extra info we are sending and their meaning is obviously // known on the receiving side too. So nevs, nfes and nrfs can be // calculated without us sending it nevs = (2+3)*ne; nfes = (1 + maxnfe + 3)*nf; nrfs = (1 + maxnrf + 3)*nr; /* Reserve nevs spots for each edge */ list_edge = (int *) malloc(5*ne*sizeof(int)); nevs = 0; /* Store the vertex ids, then the 3 auxilliary data fields */ for(i = 0; i < ne; i++) { me = MESH_Edge(mesh,i); list_edge[nevs] = MV_ID(ME_Vertex(me,0)); list_edge[nevs+1] = MV_ID(ME_Vertex(me,1)); list_edge[nevs+2] = (ME_GEntID(me)<<3) | (ME_GEntDim(me)); list_edge[nevs+3] = (ME_MasterParID(me) <<3) | (ME_OnParBoundary(me)<<2) | (ME_PType(me)); list_edge[nevs+4] = ME_GlobalID(me); nevs += 5; } /* send detailed edge info */ MPI_Isend(list_edge,nevs,MPI_INT,torank,torank,comm,&mpirequest); (*requests)[*numreq] = mpirequest; (*numreq)++; /* Reserve nfes spots for each face */ list_face = (int *) malloc(nfes*sizeof(int)); nfes = 0; /* first store nfe, then the edge ids, then the 3 auxilliary data fields */ for(i = 0; i < nf; i++) { mf = MESH_Face(mesh,i); mfedges = MF_Edges(mf,1,0); nfe = List_Num_Entries(mfedges); list_face[nfes] = nfe; for(j = 0; j < nfe; j++) { dir = MF_EdgeDir_i(mf,j) ? 1 : -1; list_face[nfes+j+1] = dir*ME_ID(List_Entry(mfedges,j)); } list_face[nfes+nfe+1] = (MF_GEntID(mf)<<3) | (MF_GEntDim(mf)); list_face[nfes+nfe+2] = (MF_MasterParID(mf)<<3) | (MF_OnParBoundary(mf)<<2) | (MF_PType(mf)); list_face[nfes+nfe+3] = MF_GlobalID(mf); nfes += (nfe + 4); List_Delete(mfedges); } /* send detailed face info */ MPI_Isend(list_face,nfes,MPI_INT,torank,torank,comm,&mpirequest); (*requests)[*numreq] = mpirequest; (*numreq)++; if (nr) { list_region = (int *) malloc(nrfs*sizeof(int)); nrfs = 0; /* first store nrf, then the face ids, then the 3 auxilliary data fields */ for(i = 0; i < nr; i++) { mr = MESH_Region(mesh,i); mrfaces = MR_Faces(mr); nrf = List_Num_Entries(mrfaces); list_region[nrfs] = nrf; for(j = 0; j < nrf; j++) { dir = MR_FaceDir_i(mr,j) == 1 ? 1 : -1; list_region[nrfs+j+1] = dir*MF_ID(List_Entry(mrfaces,j)); } list_region[nrfs+nrf+1] = (MR_GEntID(mr)<<3) | (MR_GEntDim(mr)); list_region[nrfs+nrf+2] = (MR_MasterParID(mr)<<3) | (MR_PType(mr)); /* MR_PType is 2 bits; 3 bit is 0 */ list_region[nrfs+nrf+3] = MR_GlobalID(mr); nrfs += (nrf + 4); List_Delete(mrfaces); } /* send detailed region info */ MPI_Isend(list_region,nrfs,MPI_INT,torank,torank,comm,&mpirequest); (*requests)[*numreq] = mpirequest; (*numreq)++; } /* collect allocated memory so it can be freed in a higher level routine after MPI_Waitall or MPI_Test has ensured that the send has been completed */ if (ptrs2free == NULL) MSTK_Report("MESH_Surf_SendMesh_FN","ptrs2free array is NULL",MSTK_FATAL); int nptrs = 3; if (*maxptrs2free == 0) { *maxptrs2free = 25; *ptrs2free = (void **) malloc(*maxptrs2free*sizeof(void *)); *numptrs2free = 0; } else if (*maxptrs2free < (*numptrs2free) + nptrs) { *maxptrs2free = 2*(*maxptrs2free) + nptrs; *ptrs2free = (void **) realloc(*ptrs2free,(*maxptrs2free)*sizeof(void *)); } if (ne) (*ptrs2free)[(*numptrs2free)++] = list_edge; if (nf) (*ptrs2free)[(*numptrs2free)++] = list_face; if (nr) (*ptrs2free)[(*numptrs2free)++] = list_region; return 1; }