int MESH_ConcatSubMesh_Face(Mesh_ptr mesh, int num, Mesh_ptr *submeshes) {
    int nfv, nfe, i, j, k, ival;
    MVertex_ptr mv, new_mv, sub_mv;
    MEdge_ptr me, new_me, sub_me;
    MFace_ptr new_mf, sub_mf;
    List_ptr mfverts, mfedges;
    int add_face, idx, global_id, iloc, *loc;
    double coor[3], rval;
    void *pval;
    Mesh_ptr submesh;

    List_ptr parbndry_verts = List_New(10);
    List_ptr parbndry_edges = List_New(10);

    MEdge_ptr *fedges = (MEdge_ptr *) malloc(MAXPV2*sizeof(MEdge_ptr));
    int *fedirs = (int *) malloc(MAXPV2*sizeof(int));

    MAttrib_ptr parbndryatt = MAttrib_New(mesh, "on_parbndry", INT, MVERTEX);
    
    /* collect edges and vertices on the partition boundary */
    int num_parbndry_edges = 0;
    idx = 0;
    while ((me = MESH_Next_Edge(mesh,&idx))) 
      if (ME_PType(me) != PINTERIOR) {
        List_Add(parbndry_edges,me);
        num_parbndry_edges++;
      }
    int num_parbndry_verts = 0;
    idx = 0;
    while ((mv = MESH_Next_Vertex(mesh,&idx)))
      if (MV_PType(mv) != PINTERIOR) {
        List_Add(parbndry_verts,mv);
        MEnt_Set_AttVal(mv, parbndryatt, 1, 0.0, NULL);
        num_parbndry_verts++;
      }
    /* sort based on global ID */
    List_Sort(parbndry_edges,num_parbndry_edges,sizeof(MEdge_ptr),compareGlobalID);
    List_Sort(parbndry_verts,num_parbndry_verts,sizeof(MVertex_ptr),compareGlobalID);

    int *parbndry_vert_gids = (int *) malloc(num_parbndry_verts*sizeof(int));
    int *parbndry_edge_gids = (int *) malloc(num_parbndry_edges*sizeof(int));

    /* store them in array for binary search */
    for (i = 0; i < num_parbndry_edges; i++) {
      me = List_Entry(parbndry_edges,i);
      parbndry_edge_gids[i] = ME_GlobalID(me);
    }
    for (i = 0; i < num_parbndry_verts; i++) {
      mv = List_Entry(parbndry_verts,i);
      parbndry_vert_gids[i] = MV_GlobalID(mv);
    }

    
    /* Make list of new edges and vertices which will be updated
       with each mesh that is concatenated */
    int max_vnew = 0, max_enew = 0;
    for (i = 0; i < num; i++) {
      max_vnew += MESH_Num_Vertices(submeshes[i]);
      max_enew += MESH_Num_Edges(submeshes[i]);
    }

    int num_new_verts = 0, num_new_edges = 0; 
    int *new_vert_gids = (int *) malloc(max_vnew*sizeof(int));
    int *new_edge_gids = (int *) malloc(max_enew*sizeof(int));

    List_ptr new_verts = List_New(max_vnew);
    List_ptr new_edges = List_New(max_enew);


    /* Now process each mesh and add a layer of ghost elements from
       each of them to the main partition */
    
    for (i = 0; i < num; i++) {
      submesh = submeshes[i];
    
      MAttrib_ptr vidatt = MAttrib_New(submesh, "tempvid", POINTER, MVERTEX);
      MAttrib_ptr eidatt = MAttrib_New(submesh, "tempeid", POINTER, MEDGE);

      idx = 0;
      while ((sub_mf = MESH_Next_Face(submesh, &idx))) {
        add_face = 0;
      
        /* Find matching vertices between the submesh and main mesh */
      
        mfverts = MF_Vertices(sub_mf,1,0);
        nfv = List_Num_Entries(mfverts);
        for (j = 0; j < nfv; j++) {
          sub_mv = List_Entry(mfverts,j);
        
          /* Does the vertex have a known counterpart on the partition
           * boundary of the main mesh? */
          MEnt_Get_AttVal(sub_mv, vidatt, &ival, &rval, &mv);

          if (mv) {
            int on_parbndry=0;
            MEnt_Get_AttVal(mv, parbndryatt, &on_parbndry, &rval, &pval);
            if (on_parbndry)
              add_face = 1; 
          } else {
        
            /* Does the global ID of this vertex of the sub mesh face
             * match the global ID of a partition boundary vertex in
             * the main mesh? */
            
            global_id = MV_GlobalID(sub_mv);
            loc = (int *) bsearch(&global_id, parbndry_vert_gids, num_parbndry_verts, sizeof(int),
                                  compareINT);
            if (loc) {  /* found a match */
              add_face = 1; 
              iloc = loc - parbndry_vert_gids;
              mv = List_Entry(parbndry_verts,iloc); 
              /* here set the ghost vertex property, only necessary when the input submeshes are not consistent */
              if (MV_PType(mv) == PGHOST && MV_PType(sub_mv) != PGHOST) {
                MV_Set_GEntDim(mv,MV_GEntDim(sub_mv));
                MV_Set_GEntID(mv,MV_GEntID(sub_mv));
              }
              
              MEnt_Set_AttVal(sub_mv, vidatt, 0, 0.0, mv);
            }
          }
        }
        List_Delete(mfverts);
      
        /* Find matching edges between the submesh and main mesh */
      
        mfedges = MF_Edges(sub_mf,1,0);
        nfe = List_Num_Entries(mfedges);
        for (j = 0; j < nfe; j++) {
          sub_me = List_Entry(mfedges,j);

          /* Does the edge have a known counterpart on the partition
           * boundary of the main mesh */
          MEnt_Get_AttVal(sub_me, eidatt, &ival, &rval, &me);

          if (!me) {
            /* Does the global ID of this edge of the sub mesh face
             * match the global ID of a partition boundary edge in the
             * main mesh? */
            
            global_id = ME_GlobalID(sub_me);
            loc = (int *) bsearch(&global_id, parbndry_edge_gids, num_parbndry_edges, sizeof(int),
                                  compareINT);
            if (loc) {
              iloc = loc - parbndry_edge_gids;
              me = List_Entry(parbndry_edges,iloc); 
              /* here set the ghost edge property, only necessary when the input submeshes are not consistent */
              if (ME_PType(me) == PGHOST && ME_PType(sub_me) != PGHOST) {
                ME_Set_GEntDim(me,ME_GEntDim(sub_me));
                ME_Set_GEntID(me,ME_GEntID(sub_me));
              }
              
	      MEnt_Set_AttVal(sub_me, eidatt, 0, 0.0, me);
            }
          }
        }
        
        if (!add_face) {
          List_Delete(mfedges);
          continue;
        }

        new_mf = MF_New(mesh); /* add face */
        MF_Set_GEntDim(new_mf,MF_GEntDim(sub_mf));
        MF_Set_GEntID(new_mf,MF_GEntID(sub_mf));
        MF_Set_PType(new_mf,PGHOST);
        MF_Set_MasterParID(new_mf,MF_MasterParID(sub_mf));
        MF_Set_GlobalID(new_mf,MF_GlobalID(sub_mf));
      
        nfe = List_Num_Entries(mfedges);
        for (j = 0; j < nfe; j++) {
          sub_me = List_Entry(mfedges,j);
          global_id = ME_GlobalID(sub_me);
          fedirs[j] = MF_EdgeDir_i(sub_mf,j) == 1 ? 1 : 0;

          new_me = NULL;	  
	  MEnt_Get_AttVal(sub_me, eidatt, &ival, &rval, &new_me);

          if (!new_me) {
            /* search in the ghost layer if another edge with
             * this global ID has been added */
            loc = (int *) bsearch(&global_id, new_edge_gids, num_new_edges,
                                  sizeof(int), compareINT);
            if (loc) {
              iloc = loc - new_edge_gids;
              new_me = List_Entry(new_edges, iloc);
              MEnt_Set_AttVal(sub_me, eidatt, 0, 0.0, new_me);
            }
          }

          if (new_me) {
            if (MV_GlobalID(ME_Vertex(new_me,0)) != MV_GlobalID(ME_Vertex(sub_me,0)))
              fedirs[j] = 1 - fedirs[j];  /* if the edge dir is not the same, reverse the edge dir */
          } else  {  /* add a new edge to main mesh */
            
            new_me = ME_New(mesh);
            ME_Set_GEntDim(new_me,ME_GEntDim(sub_me));
            ME_Set_GEntID(new_me,ME_GEntID(sub_me));
            ME_Set_PType(new_me,PGHOST);
            ME_Set_MasterParID(new_me,ME_MasterParID(sub_me));
            ME_Set_GlobalID(new_me,ME_GlobalID(sub_me));
	  
            MEnt_Set_AttVal(sub_me, eidatt, 0, 0.0, new_me);
	    List_Add(new_edges, new_me);
          
            for (k = 0; k < 2; k++) {
              sub_mv = ME_Vertex(sub_me,k);
              global_id = MV_GlobalID(sub_mv);

              new_mv = NULL;
              MEnt_Get_AttVal(sub_mv, vidatt, &ival, &rval, &new_mv);
	      if (!new_mv) {
		/* search in the ghost layer if another vertex with
                 * this global ID has been added */
                loc = (int *) bsearch(&global_id, new_vert_gids, num_new_verts,
                                      sizeof(int), compareINT);
                if (loc) {
                  iloc = loc - new_vert_gids;
                  new_mv = List_Entry(new_verts, iloc);
                  MEnt_Set_AttVal(sub_mv, vidatt, 0, 0.0, new_mv);
                }
              }

              if (!new_mv) {  /* add a new vertex to main mesh */
                new_mv = MV_New(mesh);
                MV_Set_GEntDim(new_mv,MV_GEntDim(sub_mv));
                MV_Set_GEntID(new_mv,MV_GEntID(sub_mv));
                MV_Set_PType(new_mv,PGHOST);
                MV_Set_MasterParID(new_mv,MV_MasterParID(sub_mv));
                MV_Set_GlobalID(new_mv,MV_GlobalID(sub_mv));
                MV_Coords(sub_mv,coor);
                MV_Set_Coords(new_mv,coor);
	      
                MEnt_Set_AttVal(sub_mv, vidatt, 0, 0.0, new_mv);
		List_Add(new_verts, new_mv);
              }
              ME_Set_Vertex(new_me,k,new_mv);  /* set edge-vertex */
            }
          }								
          fedges[j] = new_me;
        }
        MF_Set_Edges(new_mf,nfe,fedges,fedirs); /* set face-edge */

        List_Delete(mfedges);
      }

      idx = 0;
      while ((sub_mv = MESH_Next_Vertex(submesh, &idx)))
	MEnt_Rem_AttVal(sub_mv, vidatt);
      MAttrib_Delete(vidatt);
      idx = 0;
      while ((sub_me = MESH_Next_Edge(submesh, &idx)))
	MEnt_Rem_AttVal(sub_me, eidatt);
      MAttrib_Delete(eidatt);

      /* Sort the added entity lists by GlobalID */
      num_new_edges = List_Num_Entries(new_edges);
      List_Sort(new_edges, num_new_edges, sizeof(MEdge_ptr), compareGlobalID);
      for (j = 0; j < num_new_edges; j++)
        new_edge_gids[j] = ME_GlobalID(List_Entry(new_edges, j));

      num_new_verts = List_Num_Entries(new_verts);
      List_Sort(new_verts, num_new_verts, sizeof(MVertex_ptr), compareGlobalID);
      for (j = 0; j < num_new_verts; j++)
        new_vert_gids[j] = MV_GlobalID(List_Entry(new_verts, j));
    }

    idx = 0;
    while ((mv = List_Next_Entry(parbndry_verts, &idx)))
      MEnt_Rem_AttVal(mv, parbndryatt);
    MAttrib_Delete(parbndryatt);
    
    List_Delete(parbndry_edges);
    List_Delete(parbndry_verts);
    List_Delete(new_edges);
    List_Delete(new_verts);

    free(parbndry_vert_gids);
    free(parbndry_edge_gids);
    free(new_vert_gids);
    free(new_edge_gids);

    free(fedges);
    free(fedirs);

    return 1;
  }
  int MESH_AssignGlobalIDs_Region(Mesh_ptr submesh, MSTK_Comm comm) {
  int i, j, k, nfv, nbf, nof, ngf, nf, nr, mesh_info[10], global_id;
  MVertex_ptr mv;
  MFace_ptr mf;
  MRegion_ptr mr;
  List_ptr boundary_faces, mfverts;
  int *loc, face_id[MAXPV2+3],index_nbf, max_nbf, iloc, is_boundary;
  int *global_mesh_info, *list_face, *recv_list_face, *face_ov_label, *id_on_ov_list;

  int rank, num;
  MPI_Comm_rank(comm,&rank);
  MPI_Comm_size(comm,&num);

  for (i = 0; i < 10; i++) mesh_info[i] = 0;
  nf = MESH_Num_Faces(submesh);
  nr = MESH_Num_Regions(submesh);
  mesh_info[3] = nf;
  mesh_info[4] = nr;
  /* 
     collect 'boundary' faces
     if endpoints are either GHOST or OVERLAP, then it is a boundary face 
  */
  nbf = 0;  boundary_faces = List_New(10);
  for(i = 0; i < nf; i++) {
    is_boundary = 1;       
    mf = MESH_Face(submesh,i);
    mfverts = MF_Vertices(mf,1,0);
    nfv = List_Num_Entries(mfverts);
    for(j = 0; j < nfv; j++) {
      mv = List_Entry(mfverts,j);
      if(MV_PType(mv)!=PGHOST && MV_PType(mv)!=POVERLAP)
	is_boundary = 0;
    }
    if(is_boundary) {
      MF_Flag_OnParBoundary(mf);
      List_Add(boundary_faces,mf);
      nbf++;
    }
    List_Delete(mfverts);
  }
  /* printf("num of boundary faces %d, on rank %d\n", nbf,rank); */
  mesh_info[6] = nbf;

  List_Sort(boundary_faces,nbf,sizeof(MFace_ptr),compareFaceID);

  global_mesh_info = (int *)malloc(10*num*sizeof(int));
  MPI_Allgather(mesh_info,10,MPI_INT,global_mesh_info,10,MPI_INT,comm);

  max_nbf = 0;
  for(i = 0; i < num; i++)
    if(max_nbf < global_mesh_info[10*i+6])
      max_nbf = global_mesh_info[10*i+6];

  list_face = (int *)malloc(max_nbf*(MAXPV2+1)*sizeof(int));
  recv_list_face = (int *)malloc(num*max_nbf*(MAXPV2+1)*sizeof(int));

  /* indicate if a face is overlapped */
  face_ov_label = (int *)malloc(num*max_nbf*sizeof(int));
  for (i = 0; i < num*max_nbf; i++)
    face_ov_label[i] = 0;
  id_on_ov_list = (int *)malloc(max_nbf*sizeof(int));
  /* pack face information to send  */
  index_nbf = 0;
  for(i = 0; i < nbf; i++) {
    mf = List_Entry(boundary_faces,i);
    mfverts = MF_Vertices(mf,1,0);
    nfv = List_Num_Entries(mfverts);
    list_face[index_nbf] = nfv;
    for(j = 0; j < nfv; j++) 
      list_face[index_nbf+j+1] = MV_GlobalID(List_Entry(mfverts,j));
    index_nbf += MAXPV2+1;
    List_Delete(mfverts);
  }

  MPI_Allgather(list_face,(MAXPV2+1)*max_nbf,MPI_INT,recv_list_face,(MAXPV2+1)*max_nbf,MPI_INT,comm);

  ngf = 0;
  /* for processor other than 0 */
  if(rank > 0) {
    for(i = 0; i < nbf; i++) {
      mf = List_Entry(boundary_faces,i);
      if(MF_GlobalID(mf) > 0) continue;  /* if already assigned */
      mfverts = MF_Vertices(mf,1,0);
      nfv = List_Num_Entries(mfverts);
      face_id[0] = nfv;                           
      for(k = 0; k < nfv; k++) 
	face_id[k+1] = MV_GlobalID(List_Entry(mfverts,k));
      List_Delete(mfverts);
      for(j = 0; j < rank; j++) {
	loc = (int *)bsearch(&face_id,
			     &recv_list_face[(MAXPV2+1)*max_nbf*j],
			     global_mesh_info[10*j+6],
			     (MAXPV2+1)*sizeof(int),
			     compareFaceINT);
	if(loc) {
	  iloc = (int)(loc - &recv_list_face[(MAXPV2+1)*max_nbf*j])/(MAXPV2+1);
	  MF_Set_PType(mf,PGHOST);
	  MF_Set_MasterParID(mf,j);
	  face_ov_label[max_nbf*j+iloc] |= 1;
	  id_on_ov_list[i] = iloc;
	  ngf++;
	  break;
	}
      }
    }
  }

 
  /* num of ghost verts */
  mesh_info[9] = ngf;
  MPI_Allgather(mesh_info,10,MPI_INT,global_mesh_info,10,MPI_INT,comm);
  /* since this is a OR reduction, we can use MPI_IN_PLACE, send buffer same as recv buffer */
  MPI_Allreduce(MPI_IN_PLACE,face_ov_label,num*max_nbf,MPI_INT,MPI_LOR,comm);    

  /* Assign global ID for non ghost face */
  global_id = 1;
  for(i = 0; i < rank; i++) 
    global_id = global_id + global_mesh_info[10*i+3] - global_mesh_info[10*i+9];
  for(i = 0; i < nf; i++) {
    mf = MESH_Face(submesh,i);
    if (MF_PType(mf) == PGHOST) continue;
    if (!MF_GlobalID(mf)) MF_Set_GlobalID(mf,global_id++);
    MF_Set_MasterParID(mf,rank);
  }

  /* label OVERLAP face */
  nof = 0;
  for(i = 0; i < nbf; i++) 
    if(face_ov_label[rank*max_nbf+i]) {
      mf = List_Entry(boundary_faces,i);
      MF_Set_PType(mf,POVERLAP);
      nof++;
  }
  /* printf("num of ghost faces %d, overlap faces %d on rank %d\n", ngf, nof, rank);  */

  /* this time only global id are sent */
  for(i = 0; i < nbf; i++) {
    mf = List_Entry(boundary_faces,i);
    list_face[i] = MF_GlobalID(mf);
  }

  MPI_Allgather(list_face,max_nbf,MPI_INT,recv_list_face,max_nbf,MPI_INT,comm);

  for(i = 0; i < nbf; i++) {
    mf = List_Entry(boundary_faces,i);
    if(MF_PType(mf)==PGHOST) {
      int gid = recv_list_face[MF_MasterParID(mf)*max_nbf+id_on_ov_list[i]];
#ifdef DEBUG
      if (MF_GlobalID(mf) && MF_GlobalID(mf) != gid)
	MSTK_Report("MESH_AssignGlobalIDs_region",
		    "Ghost face already has different global ID", MSTK_WARN);
#endif
      MF_Set_GlobalID(mf, gid);
    }
  }

  /* assign region global id */
  global_id = 1;
  for(i = 0; i < rank; i++) 
    global_id = global_id + global_mesh_info[10*i+4];
  for(i = 0; i < nr; i++) {
    mr = MESH_Region(submesh,i);
    MR_Set_PType(mr,PINTERIOR);
    if (!MR_GlobalID(mr)) MR_Set_GlobalID(mr,global_id++);
    MR_Set_MasterParID(mr,rank);
  }


  List_Delete(boundary_faces);
  free(global_mesh_info);
  free(face_ov_label);
  free(id_on_ov_list);
  free(list_face);
  free(recv_list_face);


  return 1;
}
  int MESH_AssignGlobalIDs_Vertex(Mesh_ptr submesh, int have_GIDs, MSTK_Comm comm) {
  int i, j, nv, nbv, ne, nf, nr, mesh_info[10];
  MVertex_ptr mv;
  List_ptr boundary_verts;
  RepType rtype;
  int index_nbv, max_nbv, iloc, num_ghost_verts, global_id;
  int *global_mesh_info, *vertex_ov_label, *vertex_ov_global_id, *id_on_ov_list;

  int rank, num;
  MPI_Comm_rank(comm,&rank);
  MPI_Comm_size(comm,&num);

  for (i = 0; i < 10; i++) mesh_info[i] = 0;

  rtype = MESH_RepType(submesh);
  nv = MESH_Num_Vertices(submesh);
  ne = MESH_Num_Edges(submesh);
  nf = MESH_Num_Faces(submesh);
  nr = MESH_Num_Regions(submesh);

  mesh_info[0] = rtype;
  mesh_info[1] = nv;
  mesh_info[2] = ne;
  mesh_info[3] = nf;
  mesh_info[4] = nr;

  /* calculate number of boundary vertices */ 
  nbv = 0;  boundary_verts = List_New(10);
  if (nr) {
    for(i = 0; i < nv; i++) {
      mv = MESH_Vertex(submesh,i);
      if (vertex_on_boundary3D(mv)) {
        MV_Flag_OnParBoundary(mv);
        List_Add(boundary_verts,mv);
        nbv++;
      }
    }
  }
  else {
    for(i = 0; i < nv; i++) {
      mv = MESH_Vertex(submesh,i);
      if (vertex_on_boundary2D(mv)) {
        MV_Flag_OnParBoundary(mv);
        List_Add(boundary_verts,mv);
        nbv++;
      }
    }
  }
  mesh_info[5] = nbv;
  
  /* 
     gather submeshes information
     right now we only need nv and nbv, and later num_ghost_verts, but we gather all mesh_info
  */
  global_mesh_info = (int *)malloc(10*num*sizeof(int));
  MPI_Allgather(mesh_info,10,MPI_INT,global_mesh_info,10,MPI_INT,comm);

  /* get largest number of boundary vertices of all the processors */
  max_nbv = 0;
  for(i = 0; i < num; i++)
    if(max_nbv < global_mesh_info[10*i+5])
      max_nbv = global_mesh_info[10*i+5];

  if (have_GIDs) {
    int *list_boundary_vertex_gid = (int *)malloc(max_nbv*sizeof(int));

    int *recv_list_vertex_gid = (int *)malloc(num*max_nbv*sizeof(int));
    
    /* sort boundary vertices based on Global ID, for binary search */
    List_Sort(boundary_verts,nbv,sizeof(MVertex_ptr),compareGlobalID);

    /* only global ids are sent */
    index_nbv = 0;
    for(i = 0; i < nbv; i++) {
      mv = List_Entry(boundary_verts,i);
      list_boundary_vertex_gid[index_nbv] = MV_GlobalID(mv);
      index_nbv++;
    }
    
    MPI_Allgather(list_boundary_vertex_gid,max_nbv,MPI_INT,recv_list_vertex_gid,max_nbv,MPI_INT,comm);
    
    /* indicate if a vertex is overlapped */
    vertex_ov_label = (int *)malloc(num*max_nbv*sizeof(int));
    
    /* 
       store the local boundary id on ov processor
       it is used to assign global id of local ghost vertices
       no need to store master partition id, MV_MasterParID(mv) is already assigned
    */
    id_on_ov_list = (int *)malloc(max_nbv*sizeof(int));

    for (i = 0; i < num*max_nbv; i++)
      vertex_ov_label[i] = 0;
    num_ghost_verts = 0;
    /* for processor other than 0 */
    if(rank > 0) {
      for(i = 0; i < nbv; i++) {
        mv = List_Entry(boundary_verts,i);
        int gid = MV_GlobalID(mv);
        /* check which previous processor has a vertex with same Global ID*/
        for(j = 0; j < rank; j++) {
          /* since each processor has sorted the boundary vertices, use binary search */
          int *loc = (int *)bsearch(&gid,
                                    &recv_list_vertex_gid[max_nbv*j],
                                    global_mesh_info[10*j+5],
                                    sizeof(int),
                                    compareINT);
          /* if found the vertex on previous processors */
          if(loc) {
            /* here the location iloc is relative to the beginning of the jth processor */
            iloc = (int)(loc - &recv_list_vertex_gid[max_nbv*j]);
            MV_Set_PType(mv,PGHOST);
            MV_Set_MasterParID(mv,j);
            num_ghost_verts++;
            /* label the original vertex as overlapped */
            vertex_ov_label[max_nbv*j+iloc] |= 1;
            id_on_ov_list[i] = iloc;
            /* if found on processor j, no need to test for j+1,j+2...*/
            break;
          }
        }
      }
    }

    free(list_boundary_vertex_gid);
    free(recv_list_vertex_gid);
  }
  else {
    double coor[3];

    int *list_boundary_vertex = (int *)malloc(max_nbv*sizeof(int));
    double *list_boundary_coor = (double *)malloc(3*max_nbv*sizeof(double));

    int *recv_list_vertex = (int *)malloc(num*max_nbv*sizeof(int));
    double *recv_list_coor = (double *)malloc(3*num*max_nbv*sizeof(double));
    
    /* sort boundary vertices based on coordinate value, for binary search */
    List_Sort(boundary_verts,nbv,sizeof(MVertex_ptr),compareVertexCoor);

    /* only local id and coordinate values are sent */
    index_nbv = 0;
    for(i = 0; i < nbv; i++) {
      mv = List_Entry(boundary_verts,i);
      list_boundary_vertex[index_nbv] = MV_ID(mv);
      MV_Coords(mv,coor);
      list_boundary_coor[index_nbv*3] = coor[0];
      list_boundary_coor[index_nbv*3+1] = coor[1];
      list_boundary_coor[index_nbv*3+2] = coor[2];
      index_nbv++;
    }
    
    MPI_Allgather(list_boundary_vertex,max_nbv,MPI_INT,recv_list_vertex,max_nbv,MPI_INT,comm);
    MPI_Allgather(list_boundary_coor,3*max_nbv,MPI_DOUBLE,recv_list_coor,3*max_nbv,MPI_DOUBLE,comm);
    
    /* indicate if a vertex is overlapped */
    vertex_ov_label = (int *)malloc(num*max_nbv*sizeof(int));
    
    /* 
       store the local boundary id on ov processor
       it is used to assign global id of local ghost vertices
       no need to store master partition id, MV_MasterParID(mv) is already assigned
    */
    id_on_ov_list = (int *)malloc(max_nbv*sizeof(int));

    for (i = 0; i < num*max_nbv; i++)
      vertex_ov_label[i] = 0;
    num_ghost_verts = 0;
    /* for processor other than 0 */
    if(rank > 0) {
      for(i = 0; i < nbv; i++) {
        mv = List_Entry(boundary_verts,i);
        MV_Coords(mv,coor);
        /* check which previous processor has the same coordinate vertex */
        for(j = 0; j < rank; j++) {
          /* since each processor has sorted the boundary vertices, use binary search */
          double *loc = (double *)bsearch(&coor,
                                          &recv_list_coor[3*max_nbv*j],
                                          global_mesh_info[10*j+5],
                                          3*sizeof(double),
                                          compareCoorDouble);
          /* if found the vertex on previous processors */
          if(loc) {
            /* here the location iloc is relative to the beginning of the jth processor */
            iloc = (int)(loc - &recv_list_coor[3*max_nbv*j])/3;
            MV_Set_PType(mv,PGHOST);
            MV_Set_MasterParID(mv,j);
            num_ghost_verts++;
            /* label the original vertex as overlapped */
            vertex_ov_label[max_nbv*j+iloc] |= 1;
            id_on_ov_list[i] = iloc;
            /* if found on processor j, no need to test for j+1,j+2...*/
            break;
          }
        }
      }
    }

    free(list_boundary_coor);
    free(recv_list_coor);
    free(list_boundary_vertex);
    free(recv_list_vertex);
  }


  /* num of ghost verts */
  mesh_info[9] = num_ghost_verts;
  /* update ghost verts number */
  MPI_Allgather(mesh_info,10,MPI_INT,global_mesh_info,10,MPI_INT,comm);
  /* since this is a OR reduction, we can use MPI_IN_PLACE, send buffer same as recv buffer */
  MPI_Allreduce(MPI_IN_PLACE,vertex_ov_label,num*max_nbv,MPI_INT,MPI_LOR,comm);    

  /* calculate starting global id number for vertices*/
  if (!have_GIDs) {
    global_id = 1;
    for(i = 0; i < rank; i++) 
      global_id = global_id + global_mesh_info[10*i+1] - global_mesh_info[10*i+9];
    for(i = 0; i < nv; i++) {
      mv = MESH_Vertex(submesh,i);
      if (MV_PType(mv) == PGHOST)
        continue;
      MV_Set_GlobalID(mv,global_id++);
      MV_Set_MasterParID(mv,rank);
    }
  }

      

  /* store overlapped vertices IDs and broadast */
  vertex_ov_global_id = (int *)malloc(num*max_nbv*sizeof(int));
  for(i = 0; i < num*max_nbv; i++) 
    vertex_ov_global_id[i] = 0;
  for(i = 0; i < nbv; i++) {
    if(vertex_ov_label[rank*max_nbv+i]) {
      mv = List_Entry(boundary_verts,i);
      MV_Set_PType(mv,POVERLAP);
      vertex_ov_global_id[rank*max_nbv+i] = MV_GlobalID(mv);
    }
  }

  MPI_Allreduce(MPI_IN_PLACE,vertex_ov_global_id,num*max_nbv,MPI_INT,MPI_MAX,comm);    

  for(i = 0; i < nbv; i++) {
    mv = List_Entry(boundary_verts,i);
    if(MV_PType(mv) == PGHOST) 
      MV_Set_GlobalID(mv,vertex_ov_global_id[MV_MasterParID(mv)*max_nbv+id_on_ov_list[i]]);
  }



  List_Delete(boundary_verts);
  free(global_mesh_info);
  free(vertex_ov_label);
  free(vertex_ov_global_id);
  free(id_on_ov_list);
  
  return 1;
}
  int MESH_AssignGlobalIDs_Edge(Mesh_ptr submesh, MSTK_Comm comm) {
  int i, j, k, nbe, noe, nge, ne, mesh_info[10], global_id;
  MVertex_ptr mv;
  MEdge_ptr me;
  List_ptr boundary_edges;
  int *loc, edge_id[2], max_nbe, index_nbe, iloc, is_boundary;
  int *global_mesh_info, *list_edge, *recv_list_edge, *edge_ov_label, *id_on_ov_list;
  
  int rank, num;
  MPI_Comm_rank(comm,&rank);
  MPI_Comm_size(comm,&num);

  for (i = 0; i < 10; i++) mesh_info[i] = 0;
  ne = MESH_Num_Edges(submesh);
  mesh_info[2] = ne;

  /* 
     collect 'boundary' edges
     if endpoints are either GHOST or OVERLAP, then it is a boundary edge 
  */
  nbe = 0;  boundary_edges = List_New(10);
  for(i = 0; i < ne; i++) {
    is_boundary = 1;       
    me = MESH_Edge(submesh,i);
    for( k = 0; k < 2; k++) {
      mv = ME_Vertex(me,k);
      if(MV_PType(mv)!=PGHOST && MV_PType(mv)!=POVERLAP)
	is_boundary = 0;
    }
    if(is_boundary) {
      ME_Flag_OnParBoundary(me);
      List_Add(boundary_edges,me);
      nbe++;
    }
  }
  mesh_info[6] = nbe;

  List_Sort(boundary_edges,nbe,sizeof(MEdge_ptr),compareEdgeID);

  global_mesh_info = (int *)malloc(10*num*sizeof(int));
  MPI_Allgather(mesh_info,10,MPI_INT,global_mesh_info,10,MPI_INT,comm);

  max_nbe = 0;
  for(i = 0; i < num; i++)
    if(max_nbe < global_mesh_info[10*i+6])
      max_nbe = global_mesh_info[10*i+6];

  list_edge = (int *)malloc(max_nbe*2*sizeof(int));
  recv_list_edge = (int *)malloc(num*max_nbe*2*sizeof(int));

  /* indicate if a edge is overlapped */
  edge_ov_label = (int *)malloc(num*max_nbe*sizeof(int));
  for (i = 0; i < num*max_nbe; i++)
    edge_ov_label[i] = 0;
  id_on_ov_list = (int *)malloc(max_nbe*sizeof(int));
  /* pack edge information to send  */
  index_nbe = 0;
  for(i = 0; i < nbe; i++) {
    me = List_Entry(boundary_edges,i);
    list_edge[index_nbe] = MV_GlobalID(ME_Vertex(me,0));
    list_edge[index_nbe+1] = MV_GlobalID(ME_Vertex(me,1));
    index_nbe += 2;
  }

  MPI_Allgather(list_edge,2*max_nbe,MPI_INT,recv_list_edge,2*max_nbe,MPI_INT,comm);

  nge = 0;
  /* for processor other than 0 */
  if(rank > 0) {
    for(i = 0; i < nbe; i++) {
      me = List_Entry(boundary_edges,i);
      if(ME_GlobalID(me) > 0) continue;  /* if already assigned */
      edge_id[0] = MV_GlobalID(ME_Vertex(me,0));
      edge_id[1] = MV_GlobalID(ME_Vertex(me,1));
      for(j = 0; j < rank; j++) {
	loc = (int *)bsearch(&edge_id,
			     &recv_list_edge[2*max_nbe*j],
			     global_mesh_info[10*j+6],
			     2*sizeof(int),
			     compareEdgeINT);
	if(loc) {
	  iloc = (int)(loc - &recv_list_edge[2*max_nbe*j])/2;
	  ME_Set_PType(me,PGHOST);
	  ME_Set_MasterParID(me,j);
	  edge_ov_label[max_nbe*j+iloc] |= 1;
	  id_on_ov_list[i] = iloc;
	  nge++;
	  break;
	}
      }
    }
  }

 
  /* num of ghost verts */
  mesh_info[9] = nge;
  MPI_Allgather(mesh_info,10,MPI_INT,global_mesh_info,10,MPI_INT,comm);
  /* since this is a OR reduction, we can use MPI_IN_PLACE, send buffer same as recv buffer */
  MPI_Allreduce(MPI_IN_PLACE,edge_ov_label,num*max_nbe,MPI_INT,MPI_LOR,comm);    

  /* Assign global ID for non ghost edge */
  global_id = 1;
  for(i = 0; i < rank; i++) 
    global_id = global_id + global_mesh_info[10*i+2] - global_mesh_info[10*i+9];
  for(i = 0; i < ne; i++) {
    me = MESH_Edge(submesh,i);
    if (ME_PType(me) == PGHOST) continue;
    if (!ME_GlobalID(me)) ME_Set_GlobalID(me,global_id++);
    ME_Set_MasterParID(me,rank);
  }

  /* label OVERLAP edge */
  noe = 0;
  for(i = 0; i < nbe; i++) 
    if(edge_ov_label[rank*max_nbe+i]) {
      me = List_Entry(boundary_edges,i);
      ME_Set_PType(me,POVERLAP);
      noe++;
  }

  /* this time only global id are sent */
  for(i = 0; i < nbe; i++) {
    me = List_Entry(boundary_edges,i);
    list_edge[i] = ME_GlobalID(me);
  }

  MPI_Allgather(list_edge,max_nbe,MPI_INT,recv_list_edge,max_nbe,MPI_INT,comm);

  for(i = 0; i < nbe; i++) {
    me = List_Entry(boundary_edges,i);
    if(ME_PType(me)==PGHOST) {
      int gid = recv_list_edge[ME_MasterParID(me)*max_nbe+id_on_ov_list[i]];
#ifdef DEBUG
      if (ME_GlobalID(me) && ME_GlobalID(me) != gid)
	MSTK_Report("MESH_Assign_GlobalIDs_Edge",
		    "Ghost edge already has different global ID",
		    MSTK_WARN);
#endif
	ME_Set_GlobalID(me, gid);
    }
  }



  List_Delete(boundary_edges);
  free(global_mesh_info);
  free(edge_ov_label);
  free(id_on_ov_list);
  free(list_edge);
  free(recv_list_edge);
  return 1;
}
Exemple #5
0
  int MESH_Send_Vertices(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;
    MVertex_ptr mv;
    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) + 11) {
      *maxreq = 2*(*maxreq) + 11;
      *requests = (MPI_Request *) realloc(*requests,*maxreq*sizeof(MPI_Request));
    }

  
    /* Now send out detailed vertex info */

    nv = MESH_Num_Vertices(mesh);
    int *list_vertex = (int *) malloc(3*nv*sizeof(int));

    /* Store the 3 auxilliary data fields - would be nice if we didn't
     * have to send the data in such an error prone way (with bit
     * shifting) that requires knowledge of the internal structure of
     * MEntity */
    
    for(i = 0; i < nv; i++) {
      mv = MESH_Vertex(mesh,i);
      list_vertex[3*i] = (MV_GEntID(mv)<<3) | (MV_GEntDim(mv));
      list_vertex[3*i+1] = (MV_MasterParID(mv) <<3) | MV_OnParBoundary(mv)<<2 | (MV_PType(mv));
      list_vertex[3*i+2] = MV_GlobalID(mv);
    }

    /* send vertices */
    MPI_Isend(list_vertex,3*nv,MPI_INT,torank,torank,comm,&mpirequest);
    (*requests)[*numreq] = mpirequest;
    (*numreq)++;;

    int nptrs = 1;

    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 *));
    }

    (*ptrs2free)[(*numptrs2free)++] = list_vertex;

    return 1;
  }