00001 #include <stdlib.h>
00002 #include "ParFUM.h"
00003 #include "ParFUM_internals.h"
00004
00005
00006
00007
00008
00009
00010 void addToLists(int *listIndex,CkVec<CkVec<int> > &lists,int chunk,int node);
00011
00012
00013 void ParFUM_generateGlobalNodeNumbers(int fem_mesh, MPI_Comm comm){
00014 int myID;
00015 int numberChunks;
00016
00017 MPI_Comm_rank(comm,&myID);
00018 MPI_Comm_size(comm,&numberChunks);
00019 FEM_Mesh *mesh = (FEM_chunk::get("ParFUM_renumberMesh"))->lookup(fem_mesh,"ParFUM_renumberMesh");
00020
00021 FEM_Node *node = &(mesh->node);
00022 FEM_IndexAttribute *globalAttr = (FEM_IndexAttribute *)node->lookup(FEM_GLOBALNO,"ParFUM_renumberMesh");
00023 AllocTable2d<int> &globalTable = globalAttr->get();
00024
00025 int *primaryChunk = new int[node->size()];
00026 IDXL_Side *shared = &(node->shared);
00027
00028 for(int i=0;i<node->size();i++){
00029
00030 globalTable[i][0]=-1;
00031
00032
00033 primaryChunk[i] = -1;
00034 }
00035
00036
00037
00038
00039 for(int j=0;j<shared->size();j++){
00040
00041 const IDXL_List &list = shared->getLocalList(j);
00042 if(list.getDest() > myID){
00043
00044 for(int k=0;k<list.size();k++){
00045 int sharedNode = list[k];
00046 if(primaryChunk[sharedNode] == -1){
00047
00048 primaryChunk[sharedNode] = myID;
00049 }
00050 }
00051 }else{
00052
00053
00054 for(int k=0;k<list.size();k++){
00055 int sharedNode = list[k];
00056 if(primaryChunk[sharedNode] == -1){
00057 primaryChunk[sharedNode] = list.getDest();
00058 }else{
00059 if(primaryChunk[sharedNode] > list.getDest()){
00060 primaryChunk[sharedNode] = list.getDest();
00061 }
00062 }
00063 }
00064 }
00065 }
00066
00067 int numberOwned = 0;
00068 for(int i=0;i<node->size();i++){
00069 if(primaryChunk[i] == -1 || primaryChunk[i] == myID){
00070 numberOwned++;
00071 }
00072 }
00073
00074 int globalUniqueNodes;
00075 MPI_Scan((void*)&numberOwned , (void* ) &globalUniqueNodes, 1, MPI_INT, MPI_SUM, comm) ;
00076
00077
00078 int firstGlobalNumber = globalUniqueNodes - numberOwned;
00079 int counter = firstGlobalNumber;
00080 for(int i=0;i<node->size();i++){
00081 if(primaryChunk[i] == -1 || primaryChunk[i] == myID){
00082 globalTable[i][0] = counter;
00083 counter++;
00084 }
00085 }
00086 CkAssert(counter == globalUniqueNodes);
00087
00088
00089
00090
00091 int *sendListIndex = new int[numberChunks];
00092 CkVec<CkVec<int> > listsToSend;
00093
00094 int *recvListIndex = new int[numberChunks];
00095 CkVec<CkVec<int> > listsToRecv;
00096
00097 for(int j=0;j<numberChunks;j++){
00098 sendListIndex[j] = recvListIndex[j] = -1;
00099 }
00100
00101
00102
00103
00104 for(int j=0;j<shared->size();j++){
00105 const IDXL_List &list = shared->getLocalList(j);
00106
00107 if(list.getDest() > myID){
00108
00109 for(int k=0;k<list.size();k++){
00110 int sharedNode = list[k];
00111 if(primaryChunk[sharedNode] == myID){
00112
00113 addToLists(sendListIndex,listsToSend,list.getDest(),sharedNode);
00114 }
00115 }
00116 }else{
00117
00118 for(int k=0;k<list.size();k++){
00119 int sharedNode = list[k];
00120 if(primaryChunk[sharedNode] == list.getDest()){
00121
00122
00123 addToLists(recvListIndex,listsToRecv,primaryChunk[sharedNode],sharedNode);
00124 }
00125 }
00126 }
00127 }
00128
00129
00130
00131 int numberRequests=0;
00132 for(int j=0;j<numberChunks;j++){
00133 if(recvListIndex[j] != -1){
00134 numberRequests++;
00135 }
00136 }
00137 MPI_Request *req = new MPI_Request[numberRequests];
00138 MPI_Status *status = new MPI_Status[numberRequests];
00139 int **recvdGlobalNumbers = new int *[numberRequests];
00140
00141 int countRequests = 0;
00142 for(int j=0;j<numberChunks;j++){
00143 if(recvListIndex[j] != -1){
00144
00145
00146 CkVec<int> &recvList = listsToRecv[recvListIndex[j]];
00147 recvdGlobalNumbers[countRequests] = new int[recvList.size()];
00148 MPI_Irecv((void*)recvdGlobalNumbers[countRequests] , recvList.size(), MPI_INT, j, 42, comm, &req[countRequests]) ;
00149 countRequests++;
00150 }
00151 }
00152
00153 for(int j=0;j<numberChunks;j++){
00154 if(sendListIndex[j] != -1){
00155
00156 CkVec<int> &sendList = listsToSend[sendListIndex[j]];
00157 int *sendGlobalNumbers = new int[sendList.size()];
00158 for(int k=0;k<sendList.size();k++){
00159 sendGlobalNumbers[k] = globalTable[sendList[k]][0];
00160 }
00161 MPI_Send((void*)sendGlobalNumbers , sendList.size(), MPI_INT, j, 42, comm);
00162 delete [] sendGlobalNumbers;
00163 }
00164 }
00165
00166 MPI_Waitall(numberRequests, req, status);
00167
00168
00169
00170 countRequests = 0;
00171 for(int j=0;j<numberChunks;j++){
00172 if(recvListIndex[j] != -1){
00173 CkVec<int> &recvList = listsToRecv[recvListIndex[j]];
00174 for(int k=0;k<recvList.size();k++){
00175 int recvNode = recvList[k];
00176 int recvGlobalNumber = recvdGlobalNumbers[countRequests][k] ;
00177 CkAssert(globalTable[recvNode][0] == -1);
00178 globalTable[recvNode][0] = recvGlobalNumber;
00179 }
00180 countRequests++;
00181 }
00182 }
00183
00184
00185
00186
00187
00188
00189
00190
00191
00192
00193
00194
00195 delete []primaryChunk;
00196 delete [] sendListIndex;
00197 delete [] recvListIndex;
00198 }
00199
00200
00201 void addToLists(int *listIndex,CkVec<CkVec<int> > &lists,int chunk,int node){
00202 if(listIndex[chunk] == -1){
00203 CkVec<int> tmpVec;
00204 int index = lists.push_back_v(tmpVec);
00205 listIndex[chunk] = index;
00206 }
00207 CkVec<int> &vec = lists[listIndex[chunk]];
00208 vec.push_back(node);
00209 }