1 #include <petscsf.h> 2 #include <petscdmswarm.h> 3 #include <petscdmda.h> 4 #include <petsc/private/dmswarmimpl.h> /*I "petscdmswarm.h" I*/ 5 #include "../src/dm/impls/swarm/data_bucket.h" 6 #include "../src/dm/impls/swarm/data_ex.h" 7 8 /* 9 User loads desired location (MPI rank) into field DMSwarm_rank 10 */ 11 PetscErrorCode DMSwarmMigrate_Push_Basic(DM dm, PetscBool remove_sent_points) { 12 DM_Swarm *swarm = (DM_Swarm *)dm->data; 13 DMSwarmDataEx de; 14 PetscInt p, npoints, *rankval, n_points_recv; 15 PetscMPIInt rank, nrank; 16 void *point_buffer, *recv_points; 17 size_t sizeof_dmswarm_point; 18 PetscBool debug = PETSC_FALSE; 19 20 PetscFunctionBegin; 21 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 22 23 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL)); 24 PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval)); 25 PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm), 0, &de)); 26 PetscCall(DMSwarmDataExTopologyInitialize(de)); 27 for (p = 0; p < npoints; ++p) { 28 nrank = rankval[p]; 29 if (nrank != rank) { PetscCall(DMSwarmDataExTopologyAddNeighbour(de, nrank)); } 30 } 31 PetscCall(DMSwarmDataExTopologyFinalize(de)); 32 PetscCall(DMSwarmDataExInitializeSendCount(de)); 33 for (p = 0; p < npoints; p++) { 34 nrank = rankval[p]; 35 if (nrank != rank) { PetscCall(DMSwarmDataExAddToSendCount(de, nrank, 1)); } 36 } 37 PetscCall(DMSwarmDataExFinalizeSendCount(de)); 38 PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db, &sizeof_dmswarm_point, &point_buffer)); 39 PetscCall(DMSwarmDataExPackInitialize(de, sizeof_dmswarm_point)); 40 for (p = 0; p < npoints; p++) { 41 nrank = rankval[p]; 42 if (nrank != rank) { 43 /* copy point into buffer */ 44 PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db, p, point_buffer)); 45 /* insert point buffer into DMSwarmDataExchanger */ 46 PetscCall(DMSwarmDataExPackData(de, nrank, 1, point_buffer)); 47 } 48 } 49 PetscCall(DMSwarmDataExPackFinalize(de)); 50 PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval)); 51 52 if (remove_sent_points) { 53 DMSwarmDataField gfield; 54 55 PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, DMSwarmField_rank, &gfield)); 56 PetscCall(DMSwarmDataFieldGetAccess(gfield)); 57 PetscCall(DMSwarmDataFieldGetEntries(gfield, (void **)&rankval)); 58 59 /* remove points which left processor */ 60 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL)); 61 for (p = 0; p < npoints; p++) { 62 nrank = rankval[p]; 63 if (nrank != rank) { 64 /* kill point */ 65 PetscCall(DMSwarmDataFieldRestoreAccess(gfield)); 66 67 PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, p)); 68 69 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL)); /* you need to update npoints as the list size decreases! */ 70 PetscCall(DMSwarmDataFieldGetAccess(gfield)); 71 PetscCall(DMSwarmDataFieldGetEntries(gfield, (void **)&rankval)); 72 p--; /* check replacement point */ 73 } 74 } 75 PetscCall(DMSwarmDataFieldRestoreEntries(gfield, (void **)&rankval)); 76 PetscCall(DMSwarmDataFieldRestoreAccess(gfield)); 77 } 78 PetscCall(DMSwarmDataExBegin(de)); 79 PetscCall(DMSwarmDataExEnd(de)); 80 PetscCall(DMSwarmDataExGetRecvData(de, &n_points_recv, (void **)&recv_points)); 81 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL)); 82 PetscCall(DMSwarmDataBucketSetSizes(swarm->db, npoints + n_points_recv, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 83 for (p = 0; p < n_points_recv; p++) { 84 void *data_p = (void *)((char *)recv_points + p * sizeof_dmswarm_point); 85 86 PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db, npoints + p, data_p)); 87 } 88 if (debug) PetscCall(DMSwarmDataExView(de)); 89 PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db, &point_buffer)); 90 PetscCall(DMSwarmDataExDestroy(de)); 91 PetscFunctionReturn(0); 92 } 93 94 PetscErrorCode DMSwarmMigrate_DMNeighborScatter(DM dm, DM dmcell, PetscBool remove_sent_points, PetscInt *npoints_prior_migration) { 95 DM_Swarm *swarm = (DM_Swarm *)dm->data; 96 DMSwarmDataEx de; 97 PetscInt r, p, npoints, *rankval, n_points_recv; 98 PetscMPIInt rank, _rank; 99 const PetscMPIInt *neighbourranks; 100 void *point_buffer, *recv_points; 101 size_t sizeof_dmswarm_point; 102 PetscInt nneighbors; 103 PetscMPIInt mynneigh, *myneigh; 104 105 PetscFunctionBegin; 106 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 107 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL)); 108 PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval)); 109 PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm), 0, &de)); 110 PetscCall(DMGetNeighbors(dmcell, &nneighbors, &neighbourranks)); 111 PetscCall(DMSwarmDataExTopologyInitialize(de)); 112 for (r = 0; r < nneighbors; r++) { 113 _rank = neighbourranks[r]; 114 if ((_rank != rank) && (_rank > 0)) { PetscCall(DMSwarmDataExTopologyAddNeighbour(de, _rank)); } 115 } 116 PetscCall(DMSwarmDataExTopologyFinalize(de)); 117 PetscCall(DMSwarmDataExTopologyGetNeighbours(de, &mynneigh, &myneigh)); 118 PetscCall(DMSwarmDataExInitializeSendCount(de)); 119 for (p = 0; p < npoints; p++) { 120 if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) { 121 for (r = 0; r < mynneigh; r++) { 122 _rank = myneigh[r]; 123 PetscCall(DMSwarmDataExAddToSendCount(de, _rank, 1)); 124 } 125 } 126 } 127 PetscCall(DMSwarmDataExFinalizeSendCount(de)); 128 PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db, &sizeof_dmswarm_point, &point_buffer)); 129 PetscCall(DMSwarmDataExPackInitialize(de, sizeof_dmswarm_point)); 130 for (p = 0; p < npoints; p++) { 131 if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) { 132 for (r = 0; r < mynneigh; r++) { 133 _rank = myneigh[r]; 134 /* copy point into buffer */ 135 PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db, p, point_buffer)); 136 /* insert point buffer into DMSwarmDataExchanger */ 137 PetscCall(DMSwarmDataExPackData(de, _rank, 1, point_buffer)); 138 } 139 } 140 } 141 PetscCall(DMSwarmDataExPackFinalize(de)); 142 PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval)); 143 if (remove_sent_points) { 144 DMSwarmDataField PField; 145 146 PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, DMSwarmField_rank, &PField)); 147 PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&rankval)); 148 /* remove points which left processor */ 149 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL)); 150 for (p = 0; p < npoints; p++) { 151 if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) { 152 /* kill point */ 153 PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, p)); 154 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL)); /* you need to update npoints as the list size decreases! */ 155 PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&rankval)); /* update date point increase realloc performed */ 156 p--; /* check replacement point */ 157 } 158 } 159 } 160 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, npoints_prior_migration, NULL, NULL)); 161 PetscCall(DMSwarmDataExBegin(de)); 162 PetscCall(DMSwarmDataExEnd(de)); 163 PetscCall(DMSwarmDataExGetRecvData(de, &n_points_recv, (void **)&recv_points)); 164 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL)); 165 PetscCall(DMSwarmDataBucketSetSizes(swarm->db, npoints + n_points_recv, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 166 for (p = 0; p < n_points_recv; p++) { 167 void *data_p = (void *)((char *)recv_points + p * sizeof_dmswarm_point); 168 169 PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db, npoints + p, data_p)); 170 } 171 PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db, &point_buffer)); 172 PetscCall(DMSwarmDataExDestroy(de)); 173 PetscFunctionReturn(0); 174 } 175 176 PetscErrorCode DMSwarmMigrate_CellDMScatter(DM dm, PetscBool remove_sent_points) { 177 DM_Swarm *swarm = (DM_Swarm *)dm->data; 178 PetscInt p, npoints, npointsg = 0, npoints2, npoints2g, *rankval, npoints_prior_migration; 179 PetscSF sfcell = NULL; 180 const PetscSFNode *LA_sfcell; 181 DM dmcell; 182 Vec pos; 183 PetscBool error_check = swarm->migrate_error_on_missing_point; 184 PetscMPIInt size, rank; 185 186 PetscFunctionBegin; 187 PetscCall(DMSwarmGetCellDM(dm, &dmcell)); 188 PetscCheck(dmcell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid if cell DM provided"); 189 190 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 191 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 192 193 #if 1 194 { 195 PetscInt *p_cellid; 196 PetscInt npoints_curr, range = 0; 197 PetscSFNode *sf_cells; 198 199 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints_curr, NULL, NULL)); 200 PetscCall(PetscMalloc1(npoints_curr, &sf_cells)); 201 202 PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval)); 203 PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&p_cellid)); 204 for (p = 0; p < npoints_curr; p++) { 205 sf_cells[p].rank = 0; 206 sf_cells[p].index = p_cellid[p]; 207 if (p_cellid[p] > range) { range = p_cellid[p]; } 208 } 209 PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&p_cellid)); 210 PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval)); 211 212 /* PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm),&sfcell)); */ 213 PetscCall(PetscSFCreate(PETSC_COMM_SELF, &sfcell)); 214 PetscCall(PetscSFSetGraph(sfcell, range, npoints_curr, NULL, PETSC_OWN_POINTER, sf_cells, PETSC_OWN_POINTER)); 215 } 216 #endif 217 218 PetscCall(DMSwarmCreateLocalVectorFromField(dm, DMSwarmPICField_coor, &pos)); 219 PetscCall(DMLocatePoints(dmcell, pos, DM_POINTLOCATION_NONE, &sfcell)); 220 PetscCall(DMSwarmDestroyLocalVectorFromField(dm, DMSwarmPICField_coor, &pos)); 221 222 if (error_check) { PetscCall(DMSwarmGetSize(dm, &npointsg)); } 223 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL)); 224 PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval)); 225 PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell)); 226 for (p = 0; p < npoints; p++) { rankval[p] = LA_sfcell[p].index; } 227 PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval)); 228 PetscCall(PetscSFDestroy(&sfcell)); 229 230 if (size > 1) { 231 PetscCall(DMSwarmMigrate_DMNeighborScatter(dm, dmcell, remove_sent_points, &npoints_prior_migration)); 232 } else { 233 DMSwarmDataField PField; 234 PetscInt npoints_curr; 235 236 /* remove points which the domain */ 237 PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, DMSwarmField_rank, &PField)); 238 PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&rankval)); 239 240 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints_curr, NULL, NULL)); 241 for (p = 0; p < npoints_curr; p++) { 242 if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) { 243 /* kill point */ 244 PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, p)); 245 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints_curr, NULL, NULL)); /* you need to update npoints as the list size decreases! */ 246 PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&rankval)); /* update date point increase realloc performed */ 247 p--; /* check replacement point */ 248 } 249 } 250 PetscCall(DMSwarmGetSize(dm, &npoints_prior_migration)); 251 } 252 253 /* locate points newly received */ 254 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints2, NULL, NULL)); 255 256 #if 0 257 { /* safe alternative - however this performs two point locations on: (i) the initial points set and; (ii) the (initial + received) point set */ 258 PetscScalar *LA_coor; 259 PetscInt bs; 260 DMSwarmDataField PField; 261 262 PetscCall(DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor)); 263 PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*npoints2,(const PetscScalar*)LA_coor,&pos)); 264 PetscCall(DMLocatePoints(dmcell,pos,DM_POINTLOCATION_NONE,&sfcell)); 265 266 PetscCall(VecDestroy(&pos)); 267 PetscCall(DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor)); 268 269 PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell)); 270 PetscCall(DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval)); 271 for (p=0; p<npoints2; p++) { 272 rankval[p] = LA_sfcell[p].index; 273 } 274 PetscCall(PetscSFDestroy(&sfcell)); 275 PetscCall(DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval)); 276 277 /* remove points which left processor */ 278 PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,DMSwarmField_rank,&PField)); 279 PetscCall(DMSwarmDataFieldGetEntries(PField,(void**)&rankval)); 280 281 PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL)); 282 for (p=0; p<npoints2; p++) { 283 if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) { 284 /* kill point */ 285 PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db,p)); 286 PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL)); /* you need to update npoints as the list size decreases! */ 287 PetscCall(DMSwarmDataFieldGetEntries(PField,(void**)&rankval)); /* update date point increase realloc performed */ 288 p--; /* check replacement point */ 289 } 290 } 291 } 292 #endif 293 294 { /* perform two point locations: (i) on the initial points set prior to communication; and (ii) on the new (received) points */ 295 PetscScalar *LA_coor; 296 PetscInt npoints_from_neighbours, bs; 297 DMSwarmDataField PField; 298 299 npoints_from_neighbours = npoints2 - npoints_prior_migration; 300 301 PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, &bs, NULL, (void **)&LA_coor)); 302 PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, bs, bs * npoints_from_neighbours, (const PetscScalar *)&LA_coor[bs * npoints_prior_migration], &pos)); 303 304 PetscCall(DMLocatePoints(dmcell, pos, DM_POINTLOCATION_NONE, &sfcell)); 305 306 PetscCall(VecDestroy(&pos)); 307 PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, &bs, NULL, (void **)&LA_coor)); 308 309 PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell)); 310 PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval)); 311 for (p = 0; p < npoints_from_neighbours; p++) { rankval[npoints_prior_migration + p] = LA_sfcell[p].index; } 312 PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval)); 313 PetscCall(PetscSFDestroy(&sfcell)); 314 315 /* remove points which left processor */ 316 PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, DMSwarmField_rank, &PField)); 317 PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&rankval)); 318 319 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints2, NULL, NULL)); 320 for (p = npoints_prior_migration; p < npoints2; p++) { 321 if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) { 322 /* kill point */ 323 PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, p)); 324 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints2, NULL, NULL)); /* you need to update npoints as the list size decreases! */ 325 PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&rankval)); /* update date point increase realloc performed */ 326 p--; /* check replacement point */ 327 } 328 } 329 } 330 331 { 332 PetscInt *p_cellid; 333 334 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints2, NULL, NULL)); 335 PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval)); 336 PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&p_cellid)); 337 for (p = 0; p < npoints2; p++) { p_cellid[p] = rankval[p]; } 338 PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&p_cellid)); 339 PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval)); 340 } 341 342 /* check for error on removed points */ 343 if (error_check) { 344 PetscCall(DMSwarmGetSize(dm, &npoints2g)); 345 PetscCheck(npointsg == npoints2g, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Points from the DMSwarm must remain constant during migration (initial %" PetscInt_FMT " - final %" PetscInt_FMT ")", npointsg, npoints2g); 346 } 347 PetscFunctionReturn(0); 348 } 349 350 PetscErrorCode DMSwarmMigrate_CellDMExact(DM dm, PetscBool remove_sent_points) { 351 PetscFunctionBegin; 352 PetscFunctionReturn(0); 353 } 354 355 /* 356 Redundant as this assumes points can only be sent to a single rank 357 */ 358 PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm, PetscInt *globalsize) { 359 DM_Swarm *swarm = (DM_Swarm *)dm->data; 360 DMSwarmDataEx de; 361 PetscInt p, npoints, *rankval, n_points_recv; 362 PetscMPIInt rank, nrank, negrank; 363 void *point_buffer, *recv_points; 364 size_t sizeof_dmswarm_point; 365 366 PetscFunctionBegin; 367 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 368 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL)); 369 *globalsize = npoints; 370 PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval)); 371 PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm), 0, &de)); 372 PetscCall(DMSwarmDataExTopologyInitialize(de)); 373 for (p = 0; p < npoints; p++) { 374 negrank = rankval[p]; 375 if (negrank < 0) { 376 nrank = -negrank - 1; 377 PetscCall(DMSwarmDataExTopologyAddNeighbour(de, nrank)); 378 } 379 } 380 PetscCall(DMSwarmDataExTopologyFinalize(de)); 381 PetscCall(DMSwarmDataExInitializeSendCount(de)); 382 for (p = 0; p < npoints; p++) { 383 negrank = rankval[p]; 384 if (negrank < 0) { 385 nrank = -negrank - 1; 386 PetscCall(DMSwarmDataExAddToSendCount(de, nrank, 1)); 387 } 388 } 389 PetscCall(DMSwarmDataExFinalizeSendCount(de)); 390 PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db, &sizeof_dmswarm_point, &point_buffer)); 391 PetscCall(DMSwarmDataExPackInitialize(de, sizeof_dmswarm_point)); 392 for (p = 0; p < npoints; p++) { 393 negrank = rankval[p]; 394 if (negrank < 0) { 395 nrank = -negrank - 1; 396 rankval[p] = nrank; 397 /* copy point into buffer */ 398 PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db, p, point_buffer)); 399 /* insert point buffer into DMSwarmDataExchanger */ 400 PetscCall(DMSwarmDataExPackData(de, nrank, 1, point_buffer)); 401 rankval[p] = negrank; 402 } 403 } 404 PetscCall(DMSwarmDataExPackFinalize(de)); 405 PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval)); 406 PetscCall(DMSwarmDataExBegin(de)); 407 PetscCall(DMSwarmDataExEnd(de)); 408 PetscCall(DMSwarmDataExGetRecvData(de, &n_points_recv, (void **)&recv_points)); 409 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL)); 410 PetscCall(DMSwarmDataBucketSetSizes(swarm->db, npoints + n_points_recv, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 411 for (p = 0; p < n_points_recv; p++) { 412 void *data_p = (void *)((char *)recv_points + p * sizeof_dmswarm_point); 413 414 PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db, npoints + p, data_p)); 415 } 416 PetscCall(DMSwarmDataExView(de)); 417 PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db, &point_buffer)); 418 PetscCall(DMSwarmDataExDestroy(de)); 419 PetscFunctionReturn(0); 420 } 421 422 typedef struct { 423 PetscMPIInt owner_rank; 424 PetscReal min[3], max[3]; 425 } CollectBBox; 426 427 PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM dm, PetscInt *globalsize) { 428 DM_Swarm *swarm = (DM_Swarm *)dm->data; 429 DMSwarmDataEx de; 430 PetscInt p, pk, npoints, *rankval, n_points_recv, n_bbox_recv, dim, neighbour_cells; 431 PetscMPIInt rank, nrank; 432 void *point_buffer, *recv_points; 433 size_t sizeof_dmswarm_point, sizeof_bbox_ctx; 434 PetscBool isdmda; 435 CollectBBox *bbox, *recv_bbox; 436 const PetscMPIInt *dmneighborranks; 437 DM dmcell; 438 439 PetscFunctionBegin; 440 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 441 442 PetscCall(DMSwarmGetCellDM(dm, &dmcell)); 443 PetscCheck(dmcell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid if cell DM provided"); 444 isdmda = PETSC_FALSE; 445 PetscObjectTypeCompare((PetscObject)dmcell, DMDA, &isdmda); 446 PetscCheck(isdmda, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only DMDA support for CollectBoundingBox"); 447 448 PetscCall(DMGetDimension(dm, &dim)); 449 sizeof_bbox_ctx = sizeof(CollectBBox); 450 PetscMalloc1(1, &bbox); 451 bbox->owner_rank = rank; 452 453 /* compute the bounding box based on the overlapping / stenctil size */ 454 { 455 Vec lcoor; 456 457 PetscCall(DMGetCoordinatesLocal(dmcell, &lcoor)); 458 if (dim >= 1) { 459 PetscCall(VecStrideMin(lcoor, 0, NULL, &bbox->min[0])); 460 PetscCall(VecStrideMax(lcoor, 0, NULL, &bbox->max[0])); 461 } 462 if (dim >= 2) { 463 PetscCall(VecStrideMin(lcoor, 1, NULL, &bbox->min[1])); 464 PetscCall(VecStrideMax(lcoor, 1, NULL, &bbox->max[1])); 465 } 466 if (dim == 3) { 467 PetscCall(VecStrideMin(lcoor, 2, NULL, &bbox->min[2])); 468 PetscCall(VecStrideMax(lcoor, 2, NULL, &bbox->max[2])); 469 } 470 } 471 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL)); 472 *globalsize = npoints; 473 PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval)); 474 PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm), 0, &de)); 475 /* use DMDA neighbours */ 476 PetscCall(DMDAGetNeighbors(dmcell, &dmneighborranks)); 477 if (dim == 1) { 478 neighbour_cells = 3; 479 } else if (dim == 2) { 480 neighbour_cells = 9; 481 } else { 482 neighbour_cells = 27; 483 } 484 PetscCall(DMSwarmDataExTopologyInitialize(de)); 485 for (p = 0; p < neighbour_cells; p++) { 486 if ((dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank)) { PetscCall(DMSwarmDataExTopologyAddNeighbour(de, dmneighborranks[p])); } 487 } 488 PetscCall(DMSwarmDataExTopologyFinalize(de)); 489 PetscCall(DMSwarmDataExInitializeSendCount(de)); 490 for (p = 0; p < neighbour_cells; p++) { 491 if ((dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank)) { PetscCall(DMSwarmDataExAddToSendCount(de, dmneighborranks[p], 1)); } 492 } 493 PetscCall(DMSwarmDataExFinalizeSendCount(de)); 494 /* send bounding boxes */ 495 PetscCall(DMSwarmDataExPackInitialize(de, sizeof_bbox_ctx)); 496 for (p = 0; p < neighbour_cells; p++) { 497 nrank = dmneighborranks[p]; 498 if ((nrank >= 0) && (nrank != rank)) { 499 /* insert bbox buffer into DMSwarmDataExchanger */ 500 PetscCall(DMSwarmDataExPackData(de, nrank, 1, bbox)); 501 } 502 } 503 PetscCall(DMSwarmDataExPackFinalize(de)); 504 /* recv bounding boxes */ 505 PetscCall(DMSwarmDataExBegin(de)); 506 PetscCall(DMSwarmDataExEnd(de)); 507 PetscCall(DMSwarmDataExGetRecvData(de, &n_bbox_recv, (void **)&recv_bbox)); 508 /* Wrong, should not be using PETSC_COMM_WORLD */ 509 for (p = 0; p < n_bbox_recv; p++) { 510 PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[rank %d]: box from %d : range[%+1.4e,%+1.4e]x[%+1.4e,%+1.4e]\n", rank, recv_bbox[p].owner_rank, (double)recv_bbox[p].min[0], (double)recv_bbox[p].max[0], (double)recv_bbox[p].min[1], 511 (double)recv_bbox[p].max[1])); 512 } 513 PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, stdout)); 514 /* of course this is stupid as this "generic" function should have a better way to know what the coordinates are called */ 515 PetscCall(DMSwarmDataExInitializeSendCount(de)); 516 for (pk = 0; pk < n_bbox_recv; pk++) { 517 PetscReal *array_x, *array_y; 518 519 PetscCall(DMSwarmGetField(dm, "coorx", NULL, NULL, (void **)&array_x)); 520 PetscCall(DMSwarmGetField(dm, "coory", NULL, NULL, (void **)&array_y)); 521 for (p = 0; p < npoints; p++) { 522 if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0])) { 523 if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1])) { PetscCall(DMSwarmDataExAddToSendCount(de, recv_bbox[pk].owner_rank, 1)); } 524 } 525 } 526 PetscCall(DMSwarmRestoreField(dm, "coory", NULL, NULL, (void **)&array_y)); 527 PetscCall(DMSwarmRestoreField(dm, "coorx", NULL, NULL, (void **)&array_x)); 528 } 529 PetscCall(DMSwarmDataExFinalizeSendCount(de)); 530 PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db, &sizeof_dmswarm_point, &point_buffer)); 531 PetscCall(DMSwarmDataExPackInitialize(de, sizeof_dmswarm_point)); 532 for (pk = 0; pk < n_bbox_recv; pk++) { 533 PetscReal *array_x, *array_y; 534 535 PetscCall(DMSwarmGetField(dm, "coorx", NULL, NULL, (void **)&array_x)); 536 PetscCall(DMSwarmGetField(dm, "coory", NULL, NULL, (void **)&array_y)); 537 for (p = 0; p < npoints; p++) { 538 if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0])) { 539 if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1])) { 540 /* copy point into buffer */ 541 PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db, p, point_buffer)); 542 /* insert point buffer into DMSwarmDataExchanger */ 543 PetscCall(DMSwarmDataExPackData(de, recv_bbox[pk].owner_rank, 1, point_buffer)); 544 } 545 } 546 } 547 PetscCall(DMSwarmRestoreField(dm, "coory", NULL, NULL, (void **)&array_y)); 548 PetscCall(DMSwarmRestoreField(dm, "coorx", NULL, NULL, (void **)&array_x)); 549 } 550 PetscCall(DMSwarmDataExPackFinalize(de)); 551 PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval)); 552 PetscCall(DMSwarmDataExBegin(de)); 553 PetscCall(DMSwarmDataExEnd(de)); 554 PetscCall(DMSwarmDataExGetRecvData(de, &n_points_recv, (void **)&recv_points)); 555 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL)); 556 PetscCall(DMSwarmDataBucketSetSizes(swarm->db, npoints + n_points_recv, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 557 for (p = 0; p < n_points_recv; p++) { 558 void *data_p = (void *)((char *)recv_points + p * sizeof_dmswarm_point); 559 560 PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db, npoints + p, data_p)); 561 } 562 PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db, &point_buffer)); 563 PetscFree(bbox); 564 PetscCall(DMSwarmDataExView(de)); 565 PetscCall(DMSwarmDataExDestroy(de)); 566 PetscFunctionReturn(0); 567 } 568 569 /* General collection when no order, or neighbour information is provided */ 570 /* 571 User provides context and collect() method 572 Broadcast user context 573 574 for each context / rank { 575 collect(swarm,context,n,list) 576 } 577 */ 578 PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM dm, PetscErrorCode (*collect)(DM, void *, PetscInt *, PetscInt **), size_t ctx_size, void *ctx, PetscInt *globalsize) { 579 DM_Swarm *swarm = (DM_Swarm *)dm->data; 580 DMSwarmDataEx de; 581 PetscInt p, r, npoints, n_points_recv; 582 PetscMPIInt size, rank; 583 void *point_buffer, *recv_points; 584 void *ctxlist; 585 PetscInt *n2collect, **collectlist; 586 size_t sizeof_dmswarm_point; 587 588 PetscFunctionBegin; 589 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 590 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 591 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL)); 592 *globalsize = npoints; 593 /* Broadcast user context */ 594 PetscMalloc(ctx_size * size, &ctxlist); 595 PetscCallMPI(MPI_Allgather(ctx, ctx_size, MPI_CHAR, ctxlist, ctx_size, MPI_CHAR, PetscObjectComm((PetscObject)dm))); 596 PetscCall(PetscMalloc1(size, &n2collect)); 597 PetscCall(PetscMalloc1(size, &collectlist)); 598 for (r = 0; r < size; r++) { 599 PetscInt _n2collect; 600 PetscInt *_collectlist; 601 void *_ctx_r; 602 603 _n2collect = 0; 604 _collectlist = NULL; 605 if (r != rank) { /* don't collect data from yourself */ 606 _ctx_r = (void *)((char *)ctxlist + r * ctx_size); 607 PetscCall(collect(dm, _ctx_r, &_n2collect, &_collectlist)); 608 } 609 n2collect[r] = _n2collect; 610 collectlist[r] = _collectlist; 611 } 612 PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm), 0, &de)); 613 /* Define topology */ 614 PetscCall(DMSwarmDataExTopologyInitialize(de)); 615 for (r = 0; r < size; r++) { 616 if (n2collect[r] > 0) { PetscCall(DMSwarmDataExTopologyAddNeighbour(de, (PetscMPIInt)r)); } 617 } 618 PetscCall(DMSwarmDataExTopologyFinalize(de)); 619 /* Define send counts */ 620 PetscCall(DMSwarmDataExInitializeSendCount(de)); 621 for (r = 0; r < size; r++) { 622 if (n2collect[r] > 0) { PetscCall(DMSwarmDataExAddToSendCount(de, r, n2collect[r])); } 623 } 624 PetscCall(DMSwarmDataExFinalizeSendCount(de)); 625 /* Pack data */ 626 PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db, &sizeof_dmswarm_point, &point_buffer)); 627 PetscCall(DMSwarmDataExPackInitialize(de, sizeof_dmswarm_point)); 628 for (r = 0; r < size; r++) { 629 for (p = 0; p < n2collect[r]; p++) { 630 PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db, collectlist[r][p], point_buffer)); 631 /* insert point buffer into the data exchanger */ 632 PetscCall(DMSwarmDataExPackData(de, r, 1, point_buffer)); 633 } 634 } 635 PetscCall(DMSwarmDataExPackFinalize(de)); 636 /* Scatter */ 637 PetscCall(DMSwarmDataExBegin(de)); 638 PetscCall(DMSwarmDataExEnd(de)); 639 /* Collect data in DMSwarm container */ 640 PetscCall(DMSwarmDataExGetRecvData(de, &n_points_recv, (void **)&recv_points)); 641 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL)); 642 PetscCall(DMSwarmDataBucketSetSizes(swarm->db, npoints + n_points_recv, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 643 for (p = 0; p < n_points_recv; p++) { 644 void *data_p = (void *)((char *)recv_points + p * sizeof_dmswarm_point); 645 646 PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db, npoints + p, data_p)); 647 } 648 /* Release memory */ 649 for (r = 0; r < size; r++) { 650 if (collectlist[r]) PetscFree(collectlist[r]); 651 } 652 PetscCall(PetscFree(collectlist)); 653 PetscCall(PetscFree(n2collect)); 654 PetscCall(PetscFree(ctxlist)); 655 PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db, &point_buffer)); 656 PetscCall(DMSwarmDataExView(de)); 657 PetscCall(DMSwarmDataExDestroy(de)); 658 PetscFunctionReturn(0); 659 } 660