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