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