int MESH_AssignGlobalIDs_Face(Mesh_ptr submesh, MSTK_Comm comm) {
  int i, nf, global_id, mesh_info[10];
  MFace_ptr mf;
  RepType rtype;
  int *global_mesh_info;

  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);
  nf = MESH_Num_Faces(submesh);

  mesh_info[0] = rtype;
  mesh_info[3] = nf;

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

  /* calculate starting global id number for faces*/
  global_id = 1;
  for(i = 0; i < rank; i++) 
    global_id = global_id + global_mesh_info[10*i+3];
  for(i = 0; i < nf; i++) {
    mf = MESH_Face(submesh,i);
    MF_Set_PType(mf,PINTERIOR);
    if (!MF_GlobalID(mf)) MF_Set_GlobalID(mf,global_id++);
    MF_Set_MasterParID(mf,rank);
  }

  free(global_mesh_info);
  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;
}
Beispiel #3
0
  int MESH_Send_Attribute(Mesh_ptr mesh, MAttrib_ptr attrib, int torank, 
                          MSTK_Comm comm, int *numreq, int *maxreq, 
                          MPI_Request **requests,
                          int *numptrs2free, int *maxptrs2free,
                          void ***ptrs2free) {
    int j, k;
    int num, ncomp, ival;
    double rval;
    void *pval;
    double *rval_arr;
    int *list_info;

    MType mtype;
    MAttType att_type;
    MEntity_ptr ment;
    MPI_Request mpirequest;

    if (requests == NULL)
      MSTK_Report("MSTK_SendMSet","Invalid MPI request buffer",MSTK_FATAL);
  
    if (*maxreq == 0) {
      *maxreq = 25;
      *requests = (MPI_Request *) malloc(*maxreq*sizeof(MPI_Request));
      *numreq = 0;
    }
    else if (*maxreq < (*numreq) + 2) {
      *maxreq *= 2;
      *requests = (MPI_Request *) realloc(*requests,*maxreq*sizeof(MPI_Request));
    }
  

    /* get attribute properties */
    att_type = MAttrib_Get_Type(attrib);
    ncomp = MAttrib_Get_NumComps(attrib);
    mtype = MAttrib_Get_EntDim(attrib);
  
    /* attribute entity type, used before ghost list established, so no ghost */
    switch (mtype) {
    case MVERTEX:
      num = MESH_Num_Vertices(mesh);
      break;
    case MEDGE:
      num = MESH_Num_Edges(mesh);
      break;
    case MFACE:
      num = MESH_Num_Faces(mesh);
      break;
    case MREGION:
      num = MESH_Num_Regions(mesh);
      break;
    default:
      num = 0;
#ifdef DEBUG2
      MSTK_Report("MESH_SendAttr()","Cannot send attributes on entity type MALLTYPE",MSTK_WARN);
#endif
      return 0;
    }
    
    /* attribute index and global id */ 

    list_info = (int *) malloc(num*sizeof(int));

    /* attribute values */

    int *list_value_int = NULL;
    double *list_value_double = NULL;
    if (att_type == INT)
      list_value_int = (int *) malloc(num*ncomp*sizeof(int));
    else
      list_value_double = (double *) malloc(num*ncomp*sizeof(double));
  
    /* collect data */
    for(j = 0; j < num; j++) {
      switch (mtype) {
      case MVERTEX:
        ment = MESH_Vertex(mesh,j);
        break;
      case MEDGE:
        ment = MESH_Edge(mesh,j);
        break;
      case MFACE:
        ment = MESH_Face(mesh,j);
        break;
      case MREGION:
        ment = MESH_Region(mesh,j);
        break;
      default:
        MSTK_Report("MESH_SendAttr()","Invalid entity type",MSTK_WARN);
        return 0;
      }
    
      MEnt_Get_AttVal(ment,attrib,&ival,&rval,&pval);

      list_info[j] = MEnt_GlobalID(ment);
      if (att_type == INT)
        list_value_int[j] = ival;
      else {
        if(ncomp == 1)
          list_value_double[j] = rval;
        if(ncomp > 1) {
          rval_arr = (double *)pval;
          for(k = 0; k < ncomp; k++)
            list_value_double[ncomp*j+k] = rval_arr[k];
        }
      }
    }
  
    /* send entity global IDs */

    MPI_Isend(list_info,num,MPI_INT,torank,torank,comm,&mpirequest);
    (*requests)[*numreq] = mpirequest;
    (*numreq)++;

    /* send values */

    if (att_type == INT) {
      MPI_Isend(list_value_int,num*ncomp,MPI_INT,torank,torank,comm,&mpirequest);
      (*requests)[*numreq] = mpirequest;
    }
    else {
      MPI_Isend(list_value_double,num*ncomp,MPI_DOUBLE,torank,torank,comm,
                &mpirequest);
      (*requests)[*numreq] = mpirequest;
    }
    (*numreq)++;

    /* track the buffers used for sending so that they can be released later */

    if (*maxptrs2free == 0) {
      *maxptrs2free = 25;
      *ptrs2free = (void **) malloc(*maxptrs2free*sizeof(void *));
      *numptrs2free = 0;
    }
    else if (*maxptrs2free < (*numptrs2free) + 2) {
      *maxptrs2free = 2*(*maxptrs2free) + 2;
      *ptrs2free = (void **) realloc(*ptrs2free,(*maxptrs2free)*sizeof(void *));
    }

    (*ptrs2free)[(*numptrs2free)++] = list_info;
    if (att_type == INT) 
      (*ptrs2free)[(*numptrs2free)++] = list_value_int;  
    else
      (*ptrs2free)[(*numptrs2free)++] = list_value_double;

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