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