[694] | 1 | #include "mpi.hpp" |
---|
[688] | 2 | #include <map> |
---|
| 3 | #include "cputime.hpp" /* time */ |
---|
| 4 | #include "meshutil.hpp" |
---|
| 5 | #include "polyg.hpp" |
---|
| 6 | #include "circle.hpp" |
---|
| 7 | #include "intersect.hpp" |
---|
| 8 | #include "intersection_ym.hpp" |
---|
| 9 | #include "errhandle.hpp" |
---|
| 10 | #include "mpi_routing.hpp" |
---|
| 11 | #include "gridRemap.hpp" |
---|
| 12 | |
---|
[2269] | 13 | #include <fstream> |
---|
| 14 | #include "client.hpp" |
---|
| 15 | |
---|
[688] | 16 | #include "mapper.hpp" |
---|
| 17 | |
---|
[2269] | 18 | namespace MemTrack |
---|
| 19 | { |
---|
| 20 | void TrackDumpBlocks(std::ofstream& memReport, size_t memtrack_blocks, size_t memtrack_size); |
---|
| 21 | } |
---|
| 22 | |
---|
[688] | 23 | namespace sphereRemap { |
---|
| 24 | |
---|
| 25 | /* A subdivition of an array into N sub-arays |
---|
| 26 | can be represented by the length of the N arrays |
---|
| 27 | or by the offsets when each of the N arrays starts. |
---|
| 28 | This function convertes from the former to the latter. */ |
---|
| 29 | void cptOffsetsFromLengths(const int *lengths, int *offsets, int sz) |
---|
| 30 | { |
---|
[1614] | 31 | offsets[0] = 0; |
---|
| 32 | for (int i = 1; i < sz; i++) |
---|
| 33 | offsets[i] = offsets[i-1] + lengths[i-1]; |
---|
[688] | 34 | } |
---|
| 35 | |
---|
| 36 | |
---|
[1614] | 37 | void Mapper::setSourceMesh(const double* boundsLon, const double* boundsLat, const double* area, int nVertex, int nbCells, const double* pole, const long int* globalId) |
---|
[688] | 38 | { |
---|
| 39 | srcGrid.pole = Coord(pole[0], pole[1], pole[2]); |
---|
| 40 | |
---|
[1614] | 41 | int mpiRank, mpiSize; |
---|
[1639] | 42 | MPI_Comm_rank(communicator, &mpiRank); |
---|
| 43 | MPI_Comm_size(communicator, &mpiSize); |
---|
[688] | 44 | |
---|
[1614] | 45 | sourceElements.reserve(nbCells); |
---|
| 46 | sourceMesh.reserve(nbCells); |
---|
[688] | 47 | sourceGlobalId.resize(nbCells) ; |
---|
| 48 | |
---|
| 49 | if (globalId==NULL) |
---|
| 50 | { |
---|
| 51 | long int offset ; |
---|
| 52 | long int nb=nbCells ; |
---|
[1639] | 53 | MPI_Scan(&nb,&offset,1,MPI_LONG,MPI_SUM,communicator) ; |
---|
[688] | 54 | offset=offset-nb ; |
---|
| 55 | for(int i=0;i<nbCells;i++) sourceGlobalId[i]=offset+i ; |
---|
| 56 | } |
---|
| 57 | else sourceGlobalId.assign(globalId,globalId+nbCells); |
---|
| 58 | |
---|
[1614] | 59 | for (int i = 0; i < nbCells; i++) |
---|
| 60 | { |
---|
| 61 | int offs = i*nVertex; |
---|
| 62 | Elt elt(boundsLon + offs, boundsLat + offs, nVertex); |
---|
| 63 | elt.src_id.rank = mpiRank; |
---|
| 64 | elt.src_id.ind = i; |
---|
[688] | 65 | elt.src_id.globalId = sourceGlobalId[i]; |
---|
[1614] | 66 | sourceElements.push_back(elt); |
---|
[2269] | 67 | sourceMesh.push_back(make_shared<Node>(elt.x, cptRadius(elt), &sourceElements.back())); |
---|
[1614] | 68 | cptEltGeom(sourceElements[i], Coord(pole[0], pole[1], pole[2])); |
---|
| 69 | if (area!=NULL) sourceElements[i].given_area=area[i] ; |
---|
| 70 | else sourceElements[i].given_area=sourceElements[i].area ; |
---|
| 71 | } |
---|
[688] | 72 | |
---|
| 73 | } |
---|
| 74 | |
---|
[1614] | 75 | void Mapper::setTargetMesh(const double* boundsLon, const double* boundsLat, const double* area, int nVertex, int nbCells, const double* pole, const long int* globalId) |
---|
[688] | 76 | { |
---|
| 77 | tgtGrid.pole = Coord(pole[0], pole[1], pole[2]); |
---|
| 78 | |
---|
[1614] | 79 | int mpiRank, mpiSize; |
---|
[1639] | 80 | MPI_Comm_rank(communicator, &mpiRank); |
---|
| 81 | MPI_Comm_size(communicator, &mpiSize); |
---|
[688] | 82 | |
---|
[1614] | 83 | targetElements.reserve(nbCells); |
---|
| 84 | targetMesh.reserve(nbCells); |
---|
[688] | 85 | |
---|
| 86 | targetGlobalId.resize(nbCells) ; |
---|
| 87 | if (globalId==NULL) |
---|
| 88 | { |
---|
| 89 | long int offset ; |
---|
| 90 | long int nb=nbCells ; |
---|
[1639] | 91 | MPI_Scan(&nb,&offset,1,MPI_LONG,MPI_SUM,communicator) ; |
---|
[688] | 92 | offset=offset-nb ; |
---|
| 93 | for(int i=0;i<nbCells;i++) targetGlobalId[i]=offset+i ; |
---|
| 94 | } |
---|
| 95 | else targetGlobalId.assign(globalId,globalId+nbCells); |
---|
| 96 | |
---|
[1614] | 97 | for (int i = 0; i < nbCells; i++) |
---|
| 98 | { |
---|
| 99 | int offs = i*nVertex; |
---|
| 100 | Elt elt(boundsLon + offs, boundsLat + offs, nVertex); |
---|
| 101 | targetElements.push_back(elt); |
---|
[2269] | 102 | targetMesh.push_back(make_shared<Node>(elt.x, cptRadius(elt), &sourceElements.back())); |
---|
[1614] | 103 | cptEltGeom(targetElements[i], Coord(pole[0], pole[1], pole[2])); |
---|
| 104 | if (area!=NULL) targetElements[i].given_area=area[i] ; |
---|
| 105 | else targetElements[i].given_area=targetElements[i].area ; |
---|
| 106 | } |
---|
[688] | 107 | |
---|
| 108 | |
---|
| 109 | } |
---|
| 110 | |
---|
| 111 | void Mapper::setSourceValue(const double* val) |
---|
| 112 | { |
---|
| 113 | int size=sourceElements.size() ; |
---|
| 114 | for(int i=0;i<size;++i) sourceElements[i].val=val[i] ; |
---|
| 115 | } |
---|
| 116 | |
---|
| 117 | void Mapper::getTargetValue(double* val) |
---|
| 118 | { |
---|
| 119 | int size=targetElements.size() ; |
---|
| 120 | for(int i=0;i<size;++i) val[i]=targetElements[i].val ; |
---|
| 121 | } |
---|
| 122 | |
---|
[1158] | 123 | vector<double> Mapper::computeWeights(int interpOrder, bool renormalize, bool quantity) |
---|
[688] | 124 | { |
---|
[1614] | 125 | vector<double> timings; |
---|
| 126 | int mpiSize, mpiRank; |
---|
[1639] | 127 | MPI_Comm_size(communicator, &mpiSize); |
---|
| 128 | MPI_Comm_rank(communicator, &mpiRank); |
---|
[688] | 129 | |
---|
| 130 | this->buildSSTree(sourceMesh, targetMesh); |
---|
[2269] | 131 | |
---|
[688] | 132 | |
---|
[1614] | 133 | if (mpiRank == 0 && verbose) cout << "Computing intersections ..." << endl; |
---|
| 134 | double tic = cputime(); |
---|
| 135 | computeIntersection(&targetElements[0], targetElements.size()); |
---|
| 136 | timings.push_back(cputime() - tic); |
---|
[688] | 137 | |
---|
[1614] | 138 | tic = cputime(); |
---|
| 139 | if (interpOrder == 2) { |
---|
| 140 | if (mpiRank == 0 && verbose) cout << "Computing grads ..." << endl; |
---|
| 141 | buildMeshTopology(); |
---|
| 142 | computeGrads(); |
---|
| 143 | } |
---|
| 144 | timings.push_back(cputime() - tic); |
---|
[688] | 145 | |
---|
[1614] | 146 | /* Prepare computation of weights */ |
---|
| 147 | /* compute number of intersections which for the first order case |
---|
[688] | 148 | corresponds to the number of edges in the remap matrix */ |
---|
[1614] | 149 | int nIntersections = 0; |
---|
| 150 | for (int j = 0; j < targetElements.size(); j++) |
---|
| 151 | { |
---|
| 152 | Elt &elt = targetElements[j]; |
---|
| 153 | for (list<Polyg*>::iterator it = elt.is.begin(); it != elt.is.end(); it++) |
---|
| 154 | nIntersections++; |
---|
| 155 | } |
---|
[2538] | 156 | |
---|
[1614] | 157 | /* overallocate for NMAX neighbours for each elements */ |
---|
[2538] | 158 | /* |
---|
[1614] | 159 | remapMatrix = new double[nIntersections*NMAX]; |
---|
| 160 | srcAddress = new int[nIntersections*NMAX]; |
---|
| 161 | srcRank = new int[nIntersections*NMAX]; |
---|
| 162 | dstAddress = new int[nIntersections*NMAX]; |
---|
[688] | 163 | sourceWeightId =new long[nIntersections*NMAX]; |
---|
| 164 | targetWeightId =new long[nIntersections*NMAX]; |
---|
[2538] | 165 | */ |
---|
[688] | 166 | |
---|
[1614] | 167 | if (mpiRank == 0 && verbose) cout << "Remapping..." << endl; |
---|
| 168 | tic = cputime(); |
---|
| 169 | nWeights = remap(&targetElements[0], targetElements.size(), interpOrder, renormalize, quantity); |
---|
| 170 | timings.push_back(cputime() - tic); |
---|
[688] | 171 | |
---|
| 172 | for (int i = 0; i < targetElements.size(); i++) targetElements[i].delete_intersections(); |
---|
| 173 | |
---|
[1614] | 174 | return timings; |
---|
[688] | 175 | } |
---|
| 176 | |
---|
| 177 | /** |
---|
| 178 | @param elements are cells of the target grid that are distributed over CPUs |
---|
| 179 | indepentently of the distribution of the SS-tree. |
---|
| 180 | @param nbElements is the size of the elements array. |
---|
| 181 | @param order is the order of interpolaton (must be 1 or 2). |
---|
| 182 | */ |
---|
[1158] | 183 | int Mapper::remap(Elt *elements, int nbElements, int order, bool renormalize, bool quantity) |
---|
[688] | 184 | { |
---|
[1614] | 185 | int mpiSize, mpiRank; |
---|
[1639] | 186 | MPI_Comm_size(communicator, &mpiSize); |
---|
| 187 | MPI_Comm_rank(communicator, &mpiRank); |
---|
[688] | 188 | |
---|
[1614] | 189 | /* create list of intersections (super mesh elements) for each rank */ |
---|
| 190 | multimap<int, Polyg *> *elementList = new multimap<int, Polyg *>[mpiSize]; |
---|
| 191 | for (int j = 0; j < nbElements; j++) |
---|
| 192 | { |
---|
| 193 | Elt& e = elements[j]; |
---|
| 194 | for (list<Polyg *>::iterator it = e.is.begin(); it != e.is.end(); it++) |
---|
| 195 | elementList[(*it)->id.rank].insert(pair<int, Polyg *>((*it)->id.ind, *it)); |
---|
| 196 | } |
---|
[688] | 197 | |
---|
[1614] | 198 | int *nbSendElement = new int[mpiSize]; |
---|
[2538] | 199 | int *nbSendVertex = new int[mpiSize]; |
---|
[1614] | 200 | int **sendElement = new int*[mpiSize]; /* indices of elements required from other rank */ |
---|
| 201 | double **recvValue = new double*[mpiSize]; |
---|
| 202 | double **recvArea = new double*[mpiSize]; |
---|
| 203 | double **recvGivenArea = new double*[mpiSize]; |
---|
[2538] | 204 | int **recvVertexOffset = new int*[mpiSize]; |
---|
[1614] | 205 | Coord **recvGrad = new Coord*[mpiSize]; |
---|
| 206 | GloId **recvNeighIds = new GloId*[mpiSize]; /* ids of the of the source neighbours which also contribute through gradient */ |
---|
| 207 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 208 | { |
---|
| 209 | /* get size for allocation */ |
---|
| 210 | int last = -1; /* compares unequal to any index */ |
---|
| 211 | int index = -1; /* increased to starting index 0 in first iteration */ |
---|
| 212 | for (multimap<int, Polyg *>::iterator it = elementList[rank].begin(); it != elementList[rank].end(); ++it) |
---|
| 213 | { |
---|
| 214 | if (last != it->first) |
---|
| 215 | index++; |
---|
| 216 | (it->second)->id.ind = index; |
---|
| 217 | last = it->first; |
---|
| 218 | } |
---|
| 219 | nbSendElement[rank] = index + 1; |
---|
[688] | 220 | |
---|
[1614] | 221 | /* if size is non-zero allocate and collect indices of elements on other ranks that we intersect */ |
---|
| 222 | if (nbSendElement[rank] > 0) |
---|
| 223 | { |
---|
| 224 | sendElement[rank] = new int[nbSendElement[rank]]; |
---|
| 225 | last = -1; |
---|
| 226 | index = -1; |
---|
| 227 | for (multimap<int, Polyg *>::iterator it = elementList[rank].begin(); it != elementList[rank].end(); ++it) |
---|
| 228 | { |
---|
| 229 | if (last != it->first) |
---|
| 230 | index++; |
---|
| 231 | sendElement[rank][index] = it->first; |
---|
| 232 | last = it->first; |
---|
| 233 | } |
---|
| 234 | } |
---|
| 235 | } |
---|
[688] | 236 | |
---|
[1614] | 237 | /* communicate sizes of source elements to be sent (index lists and later values and gradients) */ |
---|
| 238 | int *nbRecvElement = new int[mpiSize]; |
---|
[2538] | 239 | int *nbRecvVertex = new int[mpiSize]; |
---|
[1639] | 240 | MPI_Alltoall(nbSendElement, 1, MPI_INT, nbRecvElement, 1, MPI_INT, communicator); |
---|
[688] | 241 | |
---|
[1614] | 242 | /* communicate indices of source elements on other ranks whoes value and gradient we need (since intersection) */ |
---|
| 243 | int nbSendRequest = 0; |
---|
| 244 | int nbRecvRequest = 0; |
---|
| 245 | int **recvElement = new int*[mpiSize]; |
---|
| 246 | double **sendValue = new double*[mpiSize]; |
---|
| 247 | double **sendArea = new double*[mpiSize]; |
---|
| 248 | double **sendGivenArea = new double*[mpiSize]; |
---|
[2538] | 249 | int **sendVertexOffset = new int*[mpiSize]; |
---|
[1614] | 250 | Coord **sendGrad = new Coord*[mpiSize]; |
---|
| 251 | GloId **sendNeighIds = new GloId*[mpiSize]; |
---|
[2538] | 252 | MPI_Request *sendRequest = new MPI_Request[6*mpiSize]; |
---|
| 253 | MPI_Request *recvRequest = new MPI_Request[6*mpiSize]; |
---|
[1614] | 254 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 255 | { |
---|
| 256 | if (nbSendElement[rank] > 0) |
---|
| 257 | { |
---|
[1639] | 258 | MPI_Issend(sendElement[rank], nbSendElement[rank], MPI_INT, rank, 0, communicator, &sendRequest[nbSendRequest]); |
---|
[1614] | 259 | nbSendRequest++; |
---|
| 260 | } |
---|
[688] | 261 | |
---|
[1614] | 262 | if (nbRecvElement[rank] > 0) |
---|
| 263 | { |
---|
| 264 | recvElement[rank] = new int[nbRecvElement[rank]]; |
---|
[1639] | 265 | MPI_Irecv(recvElement[rank], nbRecvElement[rank], MPI_INT, rank, 0, communicator, &recvRequest[nbRecvRequest]); |
---|
[1614] | 266 | nbRecvRequest++; |
---|
| 267 | } |
---|
| 268 | } |
---|
[2538] | 269 | MPI_Status *status = new MPI_Status[6*mpiSize]; |
---|
[1614] | 270 | |
---|
[1639] | 271 | MPI_Waitall(nbSendRequest, sendRequest, status); |
---|
[2538] | 272 | MPI_Waitall(nbRecvRequest, recvRequest, status); |
---|
[688] | 273 | |
---|
[2538] | 274 | nbSendRequest = 0; |
---|
| 275 | nbRecvRequest = 0; |
---|
| 276 | |
---|
| 277 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 278 | { |
---|
| 279 | if (nbRecvElement[rank] > 0) |
---|
| 280 | { |
---|
| 281 | nbRecvVertex[rank]=0 ; |
---|
| 282 | for (int j = 0; j < nbRecvElement[rank]; j++) nbRecvVertex[rank]+=sstree.localElements[recvElement[rank][j]].n+1; |
---|
| 283 | MPI_Issend(&nbRecvVertex[rank], 1, MPI_INT, rank, 0, communicator, &sendRequest[nbSendRequest]); |
---|
| 284 | nbSendRequest++; |
---|
| 285 | } |
---|
| 286 | |
---|
| 287 | if (nbSendElement[rank] > 0) |
---|
| 288 | { |
---|
| 289 | MPI_Irecv(&nbSendVertex[rank], 1, MPI_INT, rank, 0, communicator, &recvRequest[nbRecvRequest]); |
---|
| 290 | nbRecvRequest++; |
---|
| 291 | } |
---|
| 292 | } |
---|
| 293 | |
---|
| 294 | MPI_Waitall(nbSendRequest, sendRequest, status); |
---|
| 295 | MPI_Waitall(nbRecvRequest, recvRequest, status); |
---|
| 296 | |
---|
| 297 | |
---|
[1614] | 298 | /* for all indices that have been received from requesting ranks: pack values and gradients, then send */ |
---|
| 299 | nbSendRequest = 0; |
---|
| 300 | nbRecvRequest = 0; |
---|
| 301 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 302 | { |
---|
| 303 | if (nbRecvElement[rank] > 0) |
---|
| 304 | { |
---|
[2538] | 305 | sendValue[rank] = new double[nbRecvElement[rank]]; |
---|
| 306 | sendArea[rank] = new double[nbRecvElement[rank]]; |
---|
| 307 | sendGivenArea[rank] = new double[nbRecvElement[rank]]; |
---|
| 308 | sendVertexOffset[rank] = new int[nbRecvElement[rank]]; |
---|
| 309 | |
---|
| 310 | if (order == 2) |
---|
| 311 | { |
---|
| 312 | sendNeighIds[rank] = new GloId[nbRecvVertex[rank]]; |
---|
| 313 | sendGrad[rank] = new Coord[nbRecvVertex[rank]]; |
---|
| 314 | } |
---|
| 315 | else sendNeighIds[rank] = new GloId[nbRecvElement[rank]]; |
---|
| 316 | |
---|
[1614] | 317 | int jj = 0; // jj == j if no weight writing |
---|
[2538] | 318 | sendVertexOffset[rank][0]=0 ; |
---|
[1614] | 319 | for (int j = 0; j < nbRecvElement[rank]; j++) |
---|
| 320 | { |
---|
| 321 | sendValue[rank][j] = sstree.localElements[recvElement[rank][j]].val; |
---|
| 322 | sendArea[rank][j] = sstree.localElements[recvElement[rank][j]].area; |
---|
| 323 | sendGivenArea[rank][j] = sstree.localElements[recvElement[rank][j]].given_area; |
---|
[2538] | 324 | |
---|
| 325 | if (j==0) |
---|
| 326 | { |
---|
| 327 | if (order==2) sendVertexOffset[rank][j]=sstree.localElements[recvElement[rank][j]].n+1 ; |
---|
| 328 | else sendVertexOffset[rank][j]=1 ; |
---|
| 329 | } |
---|
| 330 | else |
---|
| 331 | { |
---|
| 332 | if (order == 2) sendVertexOffset[rank][j] = sendVertexOffset[rank][j-1]+sstree.localElements[recvElement[rank][j]].n+1 ; |
---|
| 333 | else sendVertexOffset[rank][j] = sendVertexOffset[rank][j-1]+1 ; |
---|
| 334 | } |
---|
| 335 | |
---|
[1614] | 336 | if (order == 2) |
---|
| 337 | { |
---|
| 338 | sendGrad[rank][jj] = sstree.localElements[recvElement[rank][j]].grad; |
---|
[844] | 339 | // cout<<"grad "<<jj<<" "<<recvElement[rank][j]<<" "<<sendGrad[rank][jj]<<" "<<sstree.localElements[recvElement[rank][j]].grad<<endl ; |
---|
[1614] | 340 | sendNeighIds[rank][jj] = sstree.localElements[recvElement[rank][j]].src_id; |
---|
| 341 | jj++; |
---|
[2538] | 342 | for (int i = 0; i < sstree.localElements[recvElement[rank][j]].n ; i++) |
---|
[1614] | 343 | { |
---|
| 344 | sendGrad[rank][jj] = sstree.localElements[recvElement[rank][j]].gradNeigh[i]; |
---|
[844] | 345 | // cout<<"grad "<<jj<<" "<<sendGrad[rank][jj]<<" "<<sstree.localElements[recvElement[rank][j]].grad<<endl ; |
---|
| 346 | sendNeighIds[rank][jj] = sstree.localElements[recvElement[rank][j]].neighId[i]; |
---|
[1614] | 347 | jj++; |
---|
| 348 | } |
---|
| 349 | } |
---|
| 350 | else |
---|
| 351 | sendNeighIds[rank][j] = sstree.localElements[recvElement[rank][j]].src_id; |
---|
| 352 | } |
---|
[1639] | 353 | MPI_Issend(sendValue[rank], nbRecvElement[rank], MPI_DOUBLE, rank, 0, communicator, &sendRequest[nbSendRequest]); |
---|
[1614] | 354 | nbSendRequest++; |
---|
[1639] | 355 | MPI_Issend(sendArea[rank], nbRecvElement[rank], MPI_DOUBLE, rank, 0, communicator, &sendRequest[nbSendRequest]); |
---|
[1614] | 356 | nbSendRequest++; |
---|
[1639] | 357 | MPI_Issend(sendGivenArea[rank], nbRecvElement[rank], MPI_DOUBLE, rank, 0, communicator, &sendRequest[nbSendRequest]); |
---|
[1614] | 358 | nbSendRequest++; |
---|
[2538] | 359 | MPI_Issend(sendVertexOffset[rank], nbRecvElement[rank], MPI_INT, rank, 0, communicator, &sendRequest[nbSendRequest]); |
---|
| 360 | nbSendRequest++; |
---|
| 361 | |
---|
[1614] | 362 | if (order == 2) |
---|
| 363 | { |
---|
[2538] | 364 | MPI_Issend(sendGrad[rank], 3*nbRecvVertex[rank], MPI_DOUBLE, rank, 0, communicator, &sendRequest[nbSendRequest]); |
---|
[1614] | 365 | nbSendRequest++; |
---|
[2538] | 366 | MPI_Issend(sendNeighIds[rank], 4*nbRecvVertex[rank], MPI_INT, rank, 0, communicator, &sendRequest[nbSendRequest]); |
---|
[688] | 367 | //ym --> attention taille GloId |
---|
[1614] | 368 | nbSendRequest++; |
---|
| 369 | } |
---|
| 370 | else |
---|
| 371 | { |
---|
[1639] | 372 | MPI_Issend(sendNeighIds[rank], 4*nbRecvElement[rank], MPI_INT, rank, 0, communicator, &sendRequest[nbSendRequest]); |
---|
[688] | 373 | //ym --> attention taille GloId |
---|
[1614] | 374 | nbSendRequest++; |
---|
| 375 | } |
---|
| 376 | } |
---|
| 377 | if (nbSendElement[rank] > 0) |
---|
| 378 | { |
---|
[2538] | 379 | recvValue[rank] = new double[nbSendElement[rank]]; |
---|
| 380 | recvArea[rank] = new double[nbSendElement[rank]]; |
---|
| 381 | recvGivenArea[rank] = new double[nbSendElement[rank]]; |
---|
| 382 | recvVertexOffset[rank] = new int[nbSendElement[rank]]; |
---|
| 383 | |
---|
| 384 | if (order == 2) |
---|
| 385 | { |
---|
| 386 | recvNeighIds[rank] = new GloId[nbSendVertex[rank]]; |
---|
| 387 | recvGrad[rank] = new Coord[nbSendVertex[rank]]; |
---|
| 388 | } |
---|
| 389 | else recvNeighIds[rank] = new GloId[nbSendElement[rank]]; |
---|
| 390 | |
---|
| 391 | |
---|
[1639] | 392 | MPI_Irecv(recvValue[rank], nbSendElement[rank], MPI_DOUBLE, rank, 0, communicator, &recvRequest[nbRecvRequest]); |
---|
[1614] | 393 | nbRecvRequest++; |
---|
[1639] | 394 | MPI_Irecv(recvArea[rank], nbSendElement[rank], MPI_DOUBLE, rank, 0, communicator, &recvRequest[nbRecvRequest]); |
---|
[1614] | 395 | nbRecvRequest++; |
---|
[1639] | 396 | MPI_Irecv(recvGivenArea[rank], nbSendElement[rank], MPI_DOUBLE, rank, 0, communicator, &recvRequest[nbRecvRequest]); |
---|
[1614] | 397 | nbRecvRequest++; |
---|
[2538] | 398 | MPI_Irecv(recvVertexOffset[rank], nbSendElement[rank], MPI_INT, rank, 0, communicator, &recvRequest[nbRecvRequest]); |
---|
| 399 | nbRecvRequest++; |
---|
| 400 | |
---|
[1614] | 401 | if (order == 2) |
---|
| 402 | { |
---|
[2538] | 403 | MPI_Irecv(recvGrad[rank], 3*nbSendVertex[rank], MPI_DOUBLE, rank, 0, communicator, &recvRequest[nbRecvRequest]); |
---|
[1614] | 404 | nbRecvRequest++; |
---|
[2538] | 405 | MPI_Irecv(recvNeighIds[rank], 4*nbSendVertex[rank], MPI_INT, rank, 0, communicator, &recvRequest[nbRecvRequest]); |
---|
[688] | 406 | //ym --> attention taille GloId |
---|
[1614] | 407 | nbRecvRequest++; |
---|
| 408 | } |
---|
| 409 | else |
---|
| 410 | { |
---|
[1639] | 411 | MPI_Irecv(recvNeighIds[rank], 4*nbSendElement[rank], MPI_INT, rank, 0, communicator, &recvRequest[nbRecvRequest]); |
---|
[688] | 412 | //ym --> attention taille GloId |
---|
[1614] | 413 | nbRecvRequest++; |
---|
| 414 | } |
---|
| 415 | } |
---|
| 416 | } |
---|
[1158] | 417 | |
---|
[2538] | 418 | MPI_Waitall(nbSendRequest, sendRequest, status); |
---|
[1639] | 419 | MPI_Waitall(nbRecvRequest, recvRequest, status); |
---|
[1614] | 420 | |
---|
[688] | 421 | |
---|
[1614] | 422 | /* now that all values and gradients are available use them to computed interpolated values on target |
---|
| 423 | and also to compute weights */ |
---|
| 424 | int i = 0; |
---|
| 425 | for (int j = 0; j < nbElements; j++) |
---|
| 426 | { |
---|
| 427 | Elt& e = elements[j]; |
---|
[688] | 428 | |
---|
[1614] | 429 | /* since for the 2nd order case source grid elements can contribute to a destination grid element over several "paths" |
---|
| 430 | (step1: gradient is computed using neighbours on same grid, step2: intersection uses several elements on other grid) |
---|
| 431 | accumulate them so that there is only one final weight between two elements */ |
---|
| 432 | map<GloId,double> wgt_map; |
---|
[688] | 433 | |
---|
[1614] | 434 | /* for destination element `e` loop over all intersetions/the corresponding source elements */ |
---|
| 435 | for (list<Polyg *>::iterator it = e.is.begin(); it != e.is.end(); it++) |
---|
| 436 | { |
---|
| 437 | /* it is the intersection element, so it->x and it->area are barycentre and area of intersection element (super mesh) |
---|
| 438 | but it->id is id of the source element that it intersects */ |
---|
| 439 | int n1 = (*it)->id.ind; |
---|
| 440 | int rank = (*it)->id.rank; |
---|
| 441 | double fk = recvValue[rank][n1]; |
---|
| 442 | double srcArea = recvArea[rank][n1]; |
---|
| 443 | double srcGivenArea = recvGivenArea[rank][n1]; |
---|
| 444 | double w = (*it)->area; |
---|
[1158] | 445 | if (quantity) w/=srcArea ; |
---|
[1614] | 446 | else w=w*srcGivenArea/srcArea*e.area/e.given_area ; |
---|
[688] | 447 | |
---|
[1614] | 448 | /* first order: src value times weight (weight = supermesh area), later divide by target area */ |
---|
[2538] | 449 | // int kk = (order == 2) ? n1 * (NMAX + 1) : n1; |
---|
| 450 | int offset ; |
---|
| 451 | int nvertex ; |
---|
| 452 | if (n1==0) |
---|
| 453 | { |
---|
| 454 | offset=0 ; |
---|
| 455 | nvertex=recvVertexOffset[rank][0] ; |
---|
| 456 | } |
---|
| 457 | else |
---|
| 458 | { |
---|
| 459 | offset = recvVertexOffset[rank][n1-1]; |
---|
| 460 | nvertex = recvVertexOffset[rank][n1]-recvVertexOffset[rank][n1-1]; |
---|
| 461 | } |
---|
| 462 | |
---|
| 463 | GloId neighID = recvNeighIds[rank][offset]; |
---|
[1614] | 464 | wgt_map[neighID] += w; |
---|
[688] | 465 | |
---|
[1614] | 466 | if (order == 2) |
---|
| 467 | { |
---|
[2538] | 468 | for (int k = 0; k < nvertex; k++) |
---|
[1614] | 469 | { |
---|
[2538] | 470 | int kk = offset + k; |
---|
[1614] | 471 | GloId neighID = recvNeighIds[rank][kk]; |
---|
| 472 | if (neighID.ind != -1) wgt_map[neighID] += w * scalarprod(recvGrad[rank][kk], (*it)->x); |
---|
| 473 | } |
---|
[688] | 474 | |
---|
[1614] | 475 | } |
---|
| 476 | } |
---|
[844] | 477 | |
---|
| 478 | double renorm=0; |
---|
| 479 | if (renormalize) |
---|
[1614] | 480 | { |
---|
| 481 | if (quantity) for (map<GloId,double>::iterator it = wgt_map.begin(); it != wgt_map.end(); it++) renorm+=it->second ; |
---|
| 482 | else for (map<GloId,double>::iterator it = wgt_map.begin(); it != wgt_map.end(); it++) renorm+=it->second / e.area; |
---|
| 483 | } |
---|
[844] | 484 | else renorm=1. ; |
---|
| 485 | |
---|
| 486 | for (map<GloId,double>::iterator it = wgt_map.begin(); it != wgt_map.end(); it++) |
---|
[1614] | 487 | { |
---|
[2538] | 488 | if (quantity) this->remapMatrix.push_back((it->second ) / renorm); |
---|
| 489 | else this->remapMatrix.push_back((it->second / e.area) / renorm); |
---|
| 490 | this->srcAddress.push_back(it->first.ind); |
---|
| 491 | this->srcRank.push_back(it->first.rank); |
---|
| 492 | this->dstAddress.push_back(j); |
---|
| 493 | this->sourceWeightId.push_back(it->first.globalId) ; |
---|
| 494 | this->targetWeightId.push_back(targetGlobalId[j]) ; |
---|
[1614] | 495 | i++; |
---|
| 496 | } |
---|
| 497 | } |
---|
[688] | 498 | |
---|
[1614] | 499 | /* free all memory allocated in this function */ |
---|
| 500 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 501 | { |
---|
| 502 | if (nbSendElement[rank] > 0) |
---|
| 503 | { |
---|
| 504 | delete[] sendElement[rank]; |
---|
| 505 | delete[] recvValue[rank]; |
---|
| 506 | delete[] recvArea[rank]; |
---|
| 507 | delete[] recvGivenArea[rank]; |
---|
[2538] | 508 | delete[] recvVertexOffset[rank]; |
---|
[1614] | 509 | if (order == 2) |
---|
| 510 | { |
---|
| 511 | delete[] recvGrad[rank]; |
---|
| 512 | } |
---|
| 513 | delete[] recvNeighIds[rank]; |
---|
| 514 | } |
---|
| 515 | if (nbRecvElement[rank] > 0) |
---|
| 516 | { |
---|
| 517 | delete[] recvElement[rank]; |
---|
| 518 | delete[] sendValue[rank]; |
---|
| 519 | delete[] sendArea[rank]; |
---|
| 520 | delete[] sendGivenArea[rank]; |
---|
[2538] | 521 | delete[] sendVertexOffset[rank]; |
---|
[1614] | 522 | if (order == 2) |
---|
| 523 | delete[] sendGrad[rank]; |
---|
| 524 | delete[] sendNeighIds[rank]; |
---|
| 525 | } |
---|
| 526 | } |
---|
| 527 | delete[] status; |
---|
| 528 | delete[] sendRequest; |
---|
| 529 | delete[] recvRequest; |
---|
| 530 | delete[] elementList; |
---|
| 531 | delete[] nbSendElement; |
---|
| 532 | delete[] nbRecvElement; |
---|
[2538] | 533 | delete[] nbSendVertex; |
---|
| 534 | delete[] nbRecvVertex; |
---|
[1614] | 535 | delete[] sendElement; |
---|
| 536 | delete[] recvElement; |
---|
| 537 | delete[] sendValue; |
---|
[2269] | 538 | delete[] sendArea ; |
---|
| 539 | delete[] sendGivenArea ; |
---|
[1614] | 540 | delete[] recvValue; |
---|
[2269] | 541 | delete[] recvArea ; |
---|
| 542 | delete[] recvGivenArea ; |
---|
[1614] | 543 | delete[] sendGrad; |
---|
| 544 | delete[] recvGrad; |
---|
| 545 | delete[] sendNeighIds; |
---|
| 546 | delete[] recvNeighIds; |
---|
| 547 | return i; |
---|
[688] | 548 | } |
---|
| 549 | |
---|
| 550 | void Mapper::computeGrads() |
---|
| 551 | { |
---|
[1614] | 552 | /* array of pointers to collect local elements and elements received from other cpu */ |
---|
| 553 | vector<Elt*> globalElements(sstree.nbLocalElements + nbNeighbourElements); |
---|
| 554 | int index = 0; |
---|
| 555 | for (int i = 0; i < sstree.nbLocalElements; i++, index++) |
---|
| 556 | globalElements[index] = &(sstree.localElements[i]); |
---|
| 557 | for (int i = 0; i < nbNeighbourElements; i++, index++) |
---|
| 558 | globalElements[index] = &neighbourElements[i]; |
---|
[688] | 559 | |
---|
[1614] | 560 | update_baryc(sstree.localElements, sstree.nbLocalElements); |
---|
| 561 | computeGradients(&globalElements[0], sstree.nbLocalElements); |
---|
[688] | 562 | } |
---|
| 563 | |
---|
| 564 | /** for each element of the source grid, finds all the neighbouring elements that share an edge |
---|
| 565 | (filling array neighbourElements). This is used later to compute gradients */ |
---|
| 566 | void Mapper::buildMeshTopology() |
---|
| 567 | { |
---|
[1614] | 568 | int mpiSize, mpiRank; |
---|
[1639] | 569 | MPI_Comm_size(communicator, &mpiSize); |
---|
| 570 | MPI_Comm_rank(communicator, &mpiRank); |
---|
[688] | 571 | |
---|
[2269] | 572 | vector<NodePtr> *routingList = new vector<NodePtr>[mpiSize]; |
---|
[1614] | 573 | vector<vector<int> > routes(sstree.localTree.leafs.size()); |
---|
[688] | 574 | |
---|
[1614] | 575 | sstree.routeIntersections(routes, sstree.localTree.leafs); |
---|
[688] | 576 | |
---|
[1614] | 577 | for (int i = 0; i < routes.size(); ++i) |
---|
| 578 | for (int k = 0; k < routes[i].size(); ++k) |
---|
| 579 | routingList[routes[i][k]].push_back(sstree.localTree.leafs[i]); |
---|
| 580 | routingList[mpiRank].clear(); |
---|
[688] | 581 | |
---|
| 582 | |
---|
[1614] | 583 | CMPIRouting mpiRoute(communicator); |
---|
| 584 | mpiRoute.init(routes); |
---|
| 585 | int nRecv = mpiRoute.getTotalSourceElement(); |
---|
[923] | 586 | // cout << mpiRank << " NRECV " << nRecv << "(" << routes.size() << ")"<< endl; |
---|
[688] | 587 | |
---|
[1614] | 588 | int *nbSendNode = new int[mpiSize]; |
---|
| 589 | int *nbRecvNode = new int[mpiSize]; |
---|
| 590 | int *sendMessageSize = new int[mpiSize]; |
---|
| 591 | int *recvMessageSize = new int[mpiSize]; |
---|
[688] | 592 | |
---|
[1614] | 593 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 594 | { |
---|
| 595 | nbSendNode[rank] = routingList[rank].size(); |
---|
| 596 | sendMessageSize[rank] = 0; |
---|
| 597 | for (size_t j = 0; j < routingList[rank].size(); j++) |
---|
| 598 | { |
---|
[2269] | 599 | Elt *elt = (Elt *) (routingList[rank][j]->data); |
---|
[1614] | 600 | sendMessageSize[rank] += packedPolygonSize(*elt); |
---|
| 601 | } |
---|
| 602 | } |
---|
[688] | 603 | |
---|
[1639] | 604 | MPI_Alltoall(nbSendNode, 1, MPI_INT, nbRecvNode, 1, MPI_INT, communicator); |
---|
| 605 | MPI_Alltoall(sendMessageSize, 1, MPI_INT, recvMessageSize, 1, MPI_INT, communicator); |
---|
[688] | 606 | |
---|
[1614] | 607 | char **sendBuffer = new char*[mpiSize]; |
---|
| 608 | char **recvBuffer = new char*[mpiSize]; |
---|
| 609 | int *pos = new int[mpiSize]; |
---|
[688] | 610 | |
---|
[1614] | 611 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 612 | { |
---|
| 613 | if (nbSendNode[rank] > 0) sendBuffer[rank] = new char[sendMessageSize[rank]]; |
---|
| 614 | if (nbRecvNode[rank] > 0) recvBuffer[rank] = new char[recvMessageSize[rank]]; |
---|
| 615 | } |
---|
[688] | 616 | |
---|
[1614] | 617 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 618 | { |
---|
| 619 | pos[rank] = 0; |
---|
| 620 | for (size_t j = 0; j < routingList[rank].size(); j++) |
---|
| 621 | { |
---|
[2269] | 622 | Elt *elt = (Elt *) (routingList[rank][j]->data); |
---|
[1614] | 623 | packPolygon(*elt, sendBuffer[rank], pos[rank]); |
---|
| 624 | } |
---|
| 625 | } |
---|
| 626 | delete [] routingList; |
---|
[688] | 627 | |
---|
| 628 | |
---|
[1614] | 629 | int nbSendRequest = 0; |
---|
| 630 | int nbRecvRequest = 0; |
---|
[1639] | 631 | MPI_Request *sendRequest = new MPI_Request[mpiSize]; |
---|
| 632 | MPI_Request *recvRequest = new MPI_Request[mpiSize]; |
---|
| 633 | MPI_Status *status = new MPI_Status[mpiSize]; |
---|
[688] | 634 | |
---|
[1614] | 635 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 636 | { |
---|
| 637 | if (nbSendNode[rank] > 0) |
---|
| 638 | { |
---|
[1639] | 639 | MPI_Issend(sendBuffer[rank], sendMessageSize[rank], MPI_CHAR, rank, 0, communicator, &sendRequest[nbSendRequest]); |
---|
[1614] | 640 | nbSendRequest++; |
---|
| 641 | } |
---|
| 642 | if (nbRecvNode[rank] > 0) |
---|
| 643 | { |
---|
[1639] | 644 | MPI_Irecv(recvBuffer[rank], recvMessageSize[rank], MPI_CHAR, rank, 0, communicator, &recvRequest[nbRecvRequest]); |
---|
[1614] | 645 | nbRecvRequest++; |
---|
| 646 | } |
---|
| 647 | } |
---|
[688] | 648 | |
---|
[1639] | 649 | MPI_Waitall(nbRecvRequest, recvRequest, status); |
---|
| 650 | MPI_Waitall(nbSendRequest, sendRequest, status); |
---|
[688] | 651 | |
---|
[1614] | 652 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 653 | if (nbSendNode[rank] > 0) delete [] sendBuffer[rank]; |
---|
| 654 | delete [] sendBuffer; |
---|
[688] | 655 | |
---|
[1614] | 656 | char **sendBuffer2 = new char*[mpiSize]; |
---|
| 657 | char **recvBuffer2 = new char*[mpiSize]; |
---|
[688] | 658 | |
---|
[1614] | 659 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 660 | { |
---|
| 661 | nbSendNode[rank] = 0; |
---|
| 662 | sendMessageSize[rank] = 0; |
---|
[688] | 663 | |
---|
[1614] | 664 | if (nbRecvNode[rank] > 0) |
---|
| 665 | { |
---|
| 666 | set<NodePtr> neighbourList; |
---|
| 667 | pos[rank] = 0; |
---|
| 668 | for (int j = 0; j < nbRecvNode[rank]; j++) |
---|
| 669 | { |
---|
| 670 | Elt elt; |
---|
| 671 | unpackPolygon(elt, recvBuffer[rank], pos[rank]); |
---|
[2269] | 672 | NodePtr node=make_shared<Node>(elt.x, cptRadius(elt), &elt); |
---|
| 673 | findNeighbour(sstree.localTree.root, node, neighbourList); |
---|
[1614] | 674 | } |
---|
| 675 | nbSendNode[rank] = neighbourList.size(); |
---|
| 676 | for (set<NodePtr>::iterator it = neighbourList.begin(); it != neighbourList.end(); it++) |
---|
| 677 | { |
---|
| 678 | Elt *elt = (Elt *) ((*it)->data); |
---|
| 679 | sendMessageSize[rank] += packedPolygonSize(*elt); |
---|
| 680 | } |
---|
[688] | 681 | |
---|
[1614] | 682 | sendBuffer2[rank] = new char[sendMessageSize[rank]]; |
---|
| 683 | pos[rank] = 0; |
---|
[688] | 684 | |
---|
[1614] | 685 | for (set<NodePtr>::iterator it = neighbourList.begin(); it != neighbourList.end(); it++) |
---|
| 686 | { |
---|
| 687 | Elt *elt = (Elt *) ((*it)->data); |
---|
| 688 | packPolygon(*elt, sendBuffer2[rank], pos[rank]); |
---|
| 689 | } |
---|
| 690 | } |
---|
| 691 | } |
---|
| 692 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 693 | if (nbRecvNode[rank] > 0) delete [] recvBuffer[rank]; |
---|
| 694 | delete [] recvBuffer; |
---|
[688] | 695 | |
---|
| 696 | |
---|
[1639] | 697 | MPI_Barrier(communicator); |
---|
| 698 | MPI_Alltoall(nbSendNode, 1, MPI_INT, nbRecvNode, 1, MPI_INT, communicator); |
---|
| 699 | MPI_Alltoall(sendMessageSize, 1, MPI_INT, recvMessageSize, 1, MPI_INT, communicator); |
---|
[688] | 700 | |
---|
[1614] | 701 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 702 | if (nbRecvNode[rank] > 0) recvBuffer2[rank] = new char[recvMessageSize[rank]]; |
---|
[688] | 703 | |
---|
[1614] | 704 | nbSendRequest = 0; |
---|
| 705 | nbRecvRequest = 0; |
---|
[688] | 706 | |
---|
[1614] | 707 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 708 | { |
---|
| 709 | if (nbSendNode[rank] > 0) |
---|
| 710 | { |
---|
[1639] | 711 | MPI_Issend(sendBuffer2[rank], sendMessageSize[rank], MPI_CHAR, rank, 0, communicator, &sendRequest[nbSendRequest]); |
---|
[1614] | 712 | nbSendRequest++; |
---|
| 713 | } |
---|
| 714 | if (nbRecvNode[rank] > 0) |
---|
| 715 | { |
---|
[1639] | 716 | MPI_Irecv(recvBuffer2[rank], recvMessageSize[rank], MPI_CHAR, rank, 0, communicator, &recvRequest[nbRecvRequest]); |
---|
[1614] | 717 | nbRecvRequest++; |
---|
| 718 | } |
---|
| 719 | } |
---|
[688] | 720 | |
---|
[1639] | 721 | MPI_Waitall(nbRecvRequest, recvRequest, status); |
---|
| 722 | MPI_Waitall(nbSendRequest, sendRequest, status); |
---|
[688] | 723 | |
---|
[1614] | 724 | int nbNeighbourNodes = 0; |
---|
| 725 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 726 | nbNeighbourNodes += nbRecvNode[rank]; |
---|
[688] | 727 | |
---|
[1614] | 728 | neighbourElements = new Elt[nbNeighbourNodes]; |
---|
| 729 | nbNeighbourElements = nbNeighbourNodes; |
---|
[688] | 730 | |
---|
[1614] | 731 | int index = 0; |
---|
| 732 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 733 | { |
---|
| 734 | pos[rank] = 0; |
---|
| 735 | for (int j = 0; j < nbRecvNode[rank]; j++) |
---|
| 736 | { |
---|
| 737 | unpackPolygon(neighbourElements[index], recvBuffer2[rank], pos[rank]); |
---|
| 738 | neighbourElements[index].id.ind = sstree.localTree.leafs.size() + index; |
---|
| 739 | index++; |
---|
| 740 | } |
---|
| 741 | } |
---|
| 742 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 743 | { |
---|
| 744 | if (nbRecvNode[rank] > 0) delete [] recvBuffer2[rank]; |
---|
| 745 | if (nbSendNode[rank] > 0) delete [] sendBuffer2[rank]; |
---|
| 746 | } |
---|
| 747 | delete [] recvBuffer2; |
---|
| 748 | delete [] sendBuffer2; |
---|
| 749 | delete [] sendMessageSize; |
---|
| 750 | delete [] recvMessageSize; |
---|
| 751 | delete [] nbSendNode; |
---|
| 752 | delete [] nbRecvNode; |
---|
| 753 | delete [] sendRequest; |
---|
| 754 | delete [] recvRequest; |
---|
| 755 | delete [] status; |
---|
| 756 | delete [] pos; |
---|
[688] | 757 | |
---|
[1614] | 758 | /* re-compute on received elements to avoid having to send this information */ |
---|
| 759 | neighbourNodes.resize(nbNeighbourNodes); |
---|
[2269] | 760 | for(int i=0;i<neighbourNodes.size();i++) neighbourNodes[i]=make_shared<Node>() ; |
---|
[1614] | 761 | setCirclesAndLinks(neighbourElements, neighbourNodes); |
---|
| 762 | cptAllEltsGeom(neighbourElements, nbNeighbourNodes, srcGrid.pole); |
---|
[688] | 763 | |
---|
[1614] | 764 | /* the local SS tree must include nodes from other cpus if they are potential |
---|
[688] | 765 | intersector of a local node */ |
---|
[1614] | 766 | sstree.localTree.insertNodes(neighbourNodes); |
---|
[688] | 767 | |
---|
[1614] | 768 | /* for every local element, |
---|
[688] | 769 | use the SS-tree to find all elements (including neighbourElements) |
---|
| 770 | who are potential neighbours because their circles intersect, |
---|
[1614] | 771 | then check all canditates for common edges to build up connectivity information |
---|
| 772 | */ |
---|
| 773 | for (int j = 0; j < sstree.localTree.leafs.size(); j++) |
---|
| 774 | { |
---|
[2269] | 775 | NodePtr node = sstree.localTree.leafs[j]; |
---|
[688] | 776 | |
---|
[1614] | 777 | /* find all leafs whoes circles that intersect node's circle and save into node->intersectors */ |
---|
[2269] | 778 | node->search(sstree.localTree.root); |
---|
[688] | 779 | |
---|
[2269] | 780 | Elt *elt = (Elt *)(node->data); |
---|
[688] | 781 | |
---|
[1614] | 782 | for (int i = 0; i < elt->n; i++) elt->neighbour[i] = NOT_FOUND; |
---|
[688] | 783 | |
---|
[1614] | 784 | /* for element `elt` loop through all nodes in the SS-tree |
---|
[688] | 785 | whoes circles intersect with the circle around `elt` (the SS intersectors) |
---|
| 786 | and check if they are neighbours in the sense that the two elements share an edge. |
---|
| 787 | If they do, save this information for elt */ |
---|
[2269] | 788 | for (list<NodePtr>::iterator it = (node->intersectors).begin(); it != (node->intersectors).end(); ++it) |
---|
[1614] | 789 | { |
---|
| 790 | Elt *elt2 = (Elt *)((*it)->data); |
---|
| 791 | set_neighbour(*elt, *elt2); |
---|
| 792 | } |
---|
[688] | 793 | |
---|
[844] | 794 | /* |
---|
[1614] | 795 | for (int i = 0; i < elt->n; i++) |
---|
| 796 | { |
---|
| 797 | if (elt->neighbour[i] == NOT_FOUND) |
---|
| 798 | error_exit("neighbour not found"); |
---|
| 799 | } |
---|
[844] | 800 | */ |
---|
[1614] | 801 | } |
---|
[688] | 802 | } |
---|
| 803 | |
---|
| 804 | /** @param elements are the target grid elements */ |
---|
| 805 | void Mapper::computeIntersection(Elt *elements, int nbElements) |
---|
| 806 | { |
---|
[1614] | 807 | int mpiSize, mpiRank; |
---|
[1639] | 808 | MPI_Comm_size(communicator, &mpiSize); |
---|
| 809 | MPI_Comm_rank(communicator, &mpiRank); |
---|
[688] | 810 | |
---|
[1639] | 811 | MPI_Barrier(communicator); |
---|
[688] | 812 | |
---|
[2269] | 813 | vector<NodePtr> *routingList = new vector<NodePtr>[mpiSize]; |
---|
[688] | 814 | |
---|
[2269] | 815 | vector<NodePtr> routeNodes; routeNodes.reserve(nbElements); |
---|
[1614] | 816 | for (int j = 0; j < nbElements; j++) |
---|
| 817 | { |
---|
| 818 | elements[j].id.ind = j; |
---|
| 819 | elements[j].id.rank = mpiRank; |
---|
[2269] | 820 | routeNodes.push_back(make_shared<Node>(elements[j].x, cptRadius(elements[j]), &elements[j])); |
---|
[1614] | 821 | } |
---|
[688] | 822 | |
---|
[1614] | 823 | vector<vector<int> > routes(routeNodes.size()); |
---|
| 824 | sstree.routeIntersections(routes, routeNodes); |
---|
| 825 | for (int i = 0; i < routes.size(); ++i) |
---|
| 826 | for (int k = 0; k < routes[i].size(); ++k) |
---|
| 827 | routingList[routes[i][k]].push_back(routeNodes[i]); |
---|
[688] | 828 | |
---|
[1614] | 829 | if (verbose >= 2) |
---|
| 830 | { |
---|
| 831 | cout << " --> rank " << mpiRank << " nbElements " << nbElements << " : "; |
---|
| 832 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 833 | cout << routingList[rank].size() << " "; |
---|
| 834 | cout << endl; |
---|
| 835 | } |
---|
[1639] | 836 | MPI_Barrier(communicator); |
---|
[688] | 837 | |
---|
[1614] | 838 | int *nbSendNode = new int[mpiSize]; |
---|
| 839 | int *nbRecvNode = new int[mpiSize]; |
---|
| 840 | int *sentMessageSize = new int[mpiSize]; |
---|
| 841 | int *recvMessageSize = new int[mpiSize]; |
---|
[688] | 842 | |
---|
[1614] | 843 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 844 | { |
---|
| 845 | nbSendNode[rank] = routingList[rank].size(); |
---|
| 846 | sentMessageSize[rank] = 0; |
---|
| 847 | for (size_t j = 0; j < routingList[rank].size(); j++) |
---|
| 848 | { |
---|
[2269] | 849 | Elt *elt = (Elt *) (routingList[rank][j]->data); |
---|
[1614] | 850 | sentMessageSize[rank] += packedPolygonSize(*elt); |
---|
| 851 | } |
---|
| 852 | } |
---|
[688] | 853 | |
---|
[1639] | 854 | MPI_Alltoall(nbSendNode, 1, MPI_INT, nbRecvNode, 1, MPI_INT, communicator); |
---|
| 855 | MPI_Alltoall(sentMessageSize, 1, MPI_INT, recvMessageSize, 1, MPI_INT, communicator); |
---|
[688] | 856 | |
---|
[1614] | 857 | int total = 0; |
---|
[688] | 858 | |
---|
[1614] | 859 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 860 | { |
---|
| 861 | total = total + nbRecvNode[rank]; |
---|
| 862 | } |
---|
[688] | 863 | |
---|
[1614] | 864 | if (verbose >= 2) cout << "---> rank " << mpiRank << " : compute intersection : total received nodes " << total << endl; |
---|
| 865 | char **sendBuffer = new char*[mpiSize]; |
---|
| 866 | char **recvBuffer = new char*[mpiSize]; |
---|
| 867 | int *pos = new int[mpiSize]; |
---|
[688] | 868 | |
---|
[1614] | 869 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 870 | { |
---|
| 871 | if (nbSendNode[rank] > 0) sendBuffer[rank] = new char[sentMessageSize[rank]]; |
---|
| 872 | if (nbRecvNode[rank] > 0) recvBuffer[rank] = new char[recvMessageSize[rank]]; |
---|
| 873 | } |
---|
[688] | 874 | |
---|
[1614] | 875 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 876 | { |
---|
| 877 | pos[rank] = 0; |
---|
| 878 | for (size_t j = 0; j < routingList[rank].size(); j++) |
---|
| 879 | { |
---|
[2269] | 880 | Elt* elt = (Elt *) (routingList[rank][j]->data); |
---|
[1614] | 881 | packPolygon(*elt, sendBuffer[rank], pos[rank]); |
---|
| 882 | } |
---|
| 883 | } |
---|
| 884 | delete [] routingList; |
---|
[688] | 885 | |
---|
[1614] | 886 | int nbSendRequest = 0; |
---|
| 887 | int nbRecvRequest = 0; |
---|
[1639] | 888 | MPI_Request *sendRequest = new MPI_Request[mpiSize]; |
---|
| 889 | MPI_Request *recvRequest = new MPI_Request[mpiSize]; |
---|
| 890 | MPI_Status *status = new MPI_Status[mpiSize]; |
---|
[688] | 891 | |
---|
[1614] | 892 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 893 | { |
---|
| 894 | if (nbSendNode[rank] > 0) |
---|
| 895 | { |
---|
[1639] | 896 | MPI_Issend(sendBuffer[rank], sentMessageSize[rank], MPI_CHAR, rank, 0, communicator, &sendRequest[nbSendRequest]); |
---|
[1614] | 897 | nbSendRequest++; |
---|
| 898 | } |
---|
| 899 | if (nbRecvNode[rank] > 0) |
---|
| 900 | { |
---|
[1639] | 901 | MPI_Irecv(recvBuffer[rank], recvMessageSize[rank], MPI_CHAR, rank, 0, communicator, &recvRequest[nbRecvRequest]); |
---|
[1614] | 902 | nbRecvRequest++; |
---|
| 903 | } |
---|
| 904 | } |
---|
[688] | 905 | |
---|
[1639] | 906 | MPI_Waitall(nbRecvRequest, recvRequest, status); |
---|
| 907 | MPI_Waitall(nbSendRequest, sendRequest, status); |
---|
[1614] | 908 | char **sendBuffer2 = new char*[mpiSize]; |
---|
| 909 | char **recvBuffer2 = new char*[mpiSize]; |
---|
[688] | 910 | |
---|
[1614] | 911 | double tic = cputime(); |
---|
| 912 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 913 | { |
---|
| 914 | sentMessageSize[rank] = 0; |
---|
[688] | 915 | |
---|
[1614] | 916 | if (nbRecvNode[rank] > 0) |
---|
| 917 | { |
---|
| 918 | Elt *recvElt = new Elt[nbRecvNode[rank]]; |
---|
| 919 | pos[rank] = 0; |
---|
| 920 | for (int j = 0; j < nbRecvNode[rank]; j++) |
---|
| 921 | { |
---|
| 922 | unpackPolygon(recvElt[j], recvBuffer[rank], pos[rank]); |
---|
| 923 | cptEltGeom(recvElt[j], tgtGrid.pole); |
---|
[2269] | 924 | NodePtr recvNode = make_shared<Node>(recvElt[j].x, cptRadius(recvElt[j]), &recvElt[j]); |
---|
| 925 | recvNode->search(sstree.localTree.root); |
---|
[1614] | 926 | /* for a node holding an element of the target, loop throught candidates for intersecting source */ |
---|
[2269] | 927 | for (list<NodePtr>::iterator it = (recvNode->intersectors).begin(); it != (recvNode->intersectors).end(); ++it) |
---|
[1614] | 928 | { |
---|
| 929 | Elt *elt2 = (Elt *) ((*it)->data); |
---|
| 930 | /* recvElt is target, elt2 is source */ |
---|
| 931 | // intersect(&recvElt[j], elt2); |
---|
| 932 | intersect_ym(&recvElt[j], elt2); |
---|
| 933 | } |
---|
[688] | 934 | |
---|
[1614] | 935 | if (recvElt[j].is.size() > 0) sentMessageSize[rank] += packIntersectionSize(recvElt[j]); |
---|
[688] | 936 | |
---|
[1614] | 937 | // here recvNode goes out of scope |
---|
| 938 | } |
---|
[688] | 939 | |
---|
[1614] | 940 | if (sentMessageSize[rank] > 0) |
---|
| 941 | { |
---|
| 942 | sentMessageSize[rank] += sizeof(int); |
---|
| 943 | sendBuffer2[rank] = new char[sentMessageSize[rank]]; |
---|
| 944 | *((int *) sendBuffer2[rank]) = 0; |
---|
| 945 | pos[rank] = sizeof(int); |
---|
| 946 | for (int j = 0; j < nbRecvNode[rank]; j++) |
---|
| 947 | { |
---|
| 948 | packIntersection(recvElt[j], sendBuffer2[rank], pos[rank]); |
---|
| 949 | //FIXME should be deleted: recvElt[j].delete_intersections(); // intersection areas have been packed to buffer and won't be used any more |
---|
| 950 | } |
---|
| 951 | } |
---|
[2269] | 952 | |
---|
| 953 | for (int j = 0; j < nbRecvNode[rank]; j++) recvElt[j].delete_intersections(); |
---|
| 954 | |
---|
[1614] | 955 | delete [] recvElt; |
---|
[688] | 956 | |
---|
[1614] | 957 | } |
---|
| 958 | } |
---|
| 959 | delete [] pos; |
---|
[688] | 960 | |
---|
[1614] | 961 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 962 | { |
---|
| 963 | if (nbSendNode[rank] > 0) delete [] sendBuffer[rank]; |
---|
| 964 | if (nbRecvNode[rank] > 0) delete [] recvBuffer[rank]; |
---|
| 965 | nbSendNode[rank] = 0; |
---|
| 966 | } |
---|
[688] | 967 | |
---|
[1614] | 968 | if (verbose >= 2) cout << "Rank " << mpiRank << " Compute (internal) intersection " << cputime() - tic << " s" << endl; |
---|
[1639] | 969 | MPI_Alltoall(sentMessageSize, 1, MPI_INT, recvMessageSize, 1, MPI_INT, communicator); |
---|
[688] | 970 | |
---|
[1614] | 971 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 972 | if (recvMessageSize[rank] > 0) |
---|
| 973 | recvBuffer2[rank] = new char[recvMessageSize[rank]]; |
---|
[688] | 974 | |
---|
[1614] | 975 | nbSendRequest = 0; |
---|
| 976 | nbRecvRequest = 0; |
---|
[688] | 977 | |
---|
[1614] | 978 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 979 | { |
---|
| 980 | if (sentMessageSize[rank] > 0) |
---|
| 981 | { |
---|
[1639] | 982 | MPI_Issend(sendBuffer2[rank], sentMessageSize[rank], MPI_CHAR, rank, 0, communicator, &sendRequest[nbSendRequest]); |
---|
[1614] | 983 | nbSendRequest++; |
---|
| 984 | } |
---|
| 985 | if (recvMessageSize[rank] > 0) |
---|
| 986 | { |
---|
[1639] | 987 | MPI_Irecv(recvBuffer2[rank], recvMessageSize[rank], MPI_CHAR, rank, 0, communicator, &recvRequest[nbRecvRequest]); |
---|
[1614] | 988 | nbRecvRequest++; |
---|
| 989 | } |
---|
| 990 | } |
---|
[688] | 991 | |
---|
[1639] | 992 | MPI_Waitall(nbRecvRequest, recvRequest, status); |
---|
| 993 | MPI_Waitall(nbSendRequest, sendRequest, status); |
---|
[688] | 994 | |
---|
[1614] | 995 | delete [] sendRequest; |
---|
| 996 | delete [] recvRequest; |
---|
| 997 | delete [] status; |
---|
| 998 | for (int rank = 0; rank < mpiSize; rank++) |
---|
| 999 | { |
---|
| 1000 | if (nbRecvNode[rank] > 0) |
---|
| 1001 | { |
---|
| 1002 | if (sentMessageSize[rank] > 0) |
---|
| 1003 | delete [] sendBuffer2[rank]; |
---|
| 1004 | } |
---|
[688] | 1005 | |
---|
[1614] | 1006 | if (recvMessageSize[rank] > 0) |
---|
| 1007 | { |
---|
| 1008 | unpackIntersection(elements, recvBuffer2[rank]); |
---|
| 1009 | delete [] recvBuffer2[rank]; |
---|
| 1010 | } |
---|
| 1011 | } |
---|
| 1012 | delete [] sendBuffer2; |
---|
| 1013 | delete [] recvBuffer2; |
---|
| 1014 | delete [] sendBuffer; |
---|
| 1015 | delete [] recvBuffer; |
---|
[688] | 1016 | |
---|
[1614] | 1017 | delete [] nbSendNode; |
---|
| 1018 | delete [] nbRecvNode; |
---|
| 1019 | delete [] sentMessageSize; |
---|
| 1020 | delete [] recvMessageSize; |
---|
[688] | 1021 | } |
---|
| 1022 | |
---|
| 1023 | Mapper::~Mapper() |
---|
| 1024 | { |
---|
[2538] | 1025 | //delete [] remapMatrix; |
---|
| 1026 | //delete [] srcAddress; |
---|
| 1027 | //delete [] srcRank; |
---|
| 1028 | //delete [] dstAddress; |
---|
| 1029 | //delete [] sourceWeightId ; |
---|
| 1030 | //delete [] targetWeightId ; |
---|
| 1031 | //if (neighbourElements) delete [] neighbourElements; |
---|
[2269] | 1032 | CTimer::release() ; |
---|
| 1033 | } |
---|
[688] | 1034 | |
---|
| 1035 | } |
---|