1df21e3a8SDave May 2df21e3a8SDave May #include <petsc/private/dmswarmimpl.h> /*I "petscdmswarm.h" I*/ 3df21e3a8SDave May #include "data_bucket.h" 4df21e3a8SDave May #include "data_ex.h" 5df21e3a8SDave May 6df21e3a8SDave May 7df21e3a8SDave May #undef __FUNCT__ 8df21e3a8SDave May #define __FUNCT__ "DMSwarmMigrate_Push_Basic" 9df21e3a8SDave May PetscErrorCode DMSwarmMigrate_Push_Basic(DM dm,PetscBool remove_sent_points) 10df21e3a8SDave May { 11df21e3a8SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 12df21e3a8SDave May PetscErrorCode ierr; 13df21e3a8SDave May DataEx de; 14df21e3a8SDave May PetscInt p,npoints,*rankval,n_points_recv; 15df21e3a8SDave May PetscMPIInt rank,nrank; 16df21e3a8SDave May void *point_buffer,*recv_points; 17df21e3a8SDave May size_t sizeof_dmswarm_point; 18df21e3a8SDave May 19df21e3a8SDave May ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 20df21e3a8SDave May 21df21e3a8SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 22df21e3a8SDave May ierr = DMSwarmGetField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 23df21e3a8SDave May 24df21e3a8SDave May de = DataExCreate(PetscObjectComm((PetscObject)dm),0);CHKERRQ(ierr); 25df21e3a8SDave May 26df21e3a8SDave May ierr = DataExTopologyInitialize(de);CHKERRQ(ierr); 27df21e3a8SDave May for (p=0; p<npoints; p++) { 28df21e3a8SDave May nrank = rankval[p]; 29df21e3a8SDave May if (nrank != rank) { 30df21e3a8SDave May ierr = DataExTopologyAddNeighbour(de,nrank);CHKERRQ(ierr); 31df21e3a8SDave May } 32df21e3a8SDave May } 33df21e3a8SDave May ierr = DataExTopologyFinalize(de);CHKERRQ(ierr); 34df21e3a8SDave May 35df21e3a8SDave May ierr = DataExInitializeSendCount(de);CHKERRQ(ierr); 36df21e3a8SDave May for (p=0; p<npoints; p++) { 37df21e3a8SDave May nrank = rankval[p]; 38df21e3a8SDave May if (nrank != rank) { 39df21e3a8SDave May ierr = DataExAddToSendCount(de,nrank,1);CHKERRQ(ierr); 40df21e3a8SDave May } 41df21e3a8SDave May } 42df21e3a8SDave May ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr); 43df21e3a8SDave May 44df21e3a8SDave May ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr); 45df21e3a8SDave May ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr); 46df21e3a8SDave May for (p=0; p<npoints; p++) { 47df21e3a8SDave May nrank = rankval[p]; 48df21e3a8SDave May if (nrank != rank) { 49df21e3a8SDave May /* copy point into buffer */ 50df21e3a8SDave May ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr); 51df21e3a8SDave May /* insert point buffer into DataExchanger */ 52df21e3a8SDave May ierr = DataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr); 53df21e3a8SDave May } 54df21e3a8SDave May } 55df21e3a8SDave May ierr = DataExPackFinalize(de);CHKERRQ(ierr); 56df21e3a8SDave May 57df21e3a8SDave May 58df21e3a8SDave May if (remove_sent_points) { 59df21e3a8SDave May /* remove points which left processor */ 60df21e3a8SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 61df21e3a8SDave May for (p=0; p<npoints; p++) { 62df21e3a8SDave May nrank = rankval[p]; 63df21e3a8SDave May if (nrank != rank) { 64df21e3a8SDave May /* kill point */ 65df21e3a8SDave May ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr); 66df21e3a8SDave May DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */ 67df21e3a8SDave May p--; /* check replacement point */ 68df21e3a8SDave May } 69df21e3a8SDave May } 70df21e3a8SDave May } 71df21e3a8SDave May ierr = DMSwarmRestoreField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 72df21e3a8SDave May 73df21e3a8SDave May ierr = DataExBegin(de);CHKERRQ(ierr); 74df21e3a8SDave May ierr = DataExEnd(de);CHKERRQ(ierr); 75df21e3a8SDave May 76df21e3a8SDave May ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr); 77df21e3a8SDave May 78df21e3a8SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 79df21e3a8SDave May ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr); 80df21e3a8SDave May for (p=0; p<n_points_recv; p++) { 81df21e3a8SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point ); 82df21e3a8SDave May 83df21e3a8SDave May ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr); 84df21e3a8SDave May } 85df21e3a8SDave May ierr = DataExView(de);CHKERRQ(ierr); 86df21e3a8SDave May ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr); 87df21e3a8SDave May ierr = DataExDestroy(de);CHKERRQ(ierr); 88df21e3a8SDave May 89df21e3a8SDave May PetscFunctionReturn(0); 90df21e3a8SDave May } 912712d1f2SDave May 922712d1f2SDave May #undef __FUNCT__ 932712d1f2SDave May #define __FUNCT__ "DMSwarmMigrate_GlobalToLocal_Basic" 942712d1f2SDave May PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt *globalsize) 952712d1f2SDave May { 962712d1f2SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 972712d1f2SDave May PetscErrorCode ierr; 982712d1f2SDave May DataEx de; 992712d1f2SDave May PetscInt p,npoints,*rankval,n_points_recv; 1002712d1f2SDave May PetscMPIInt rank,nrank,negrank; 1012712d1f2SDave May void *point_buffer,*recv_points; 1022712d1f2SDave May size_t sizeof_dmswarm_point; 1032712d1f2SDave May 1042712d1f2SDave May ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 1052712d1f2SDave May 1062712d1f2SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 1072712d1f2SDave May *globalsize = npoints; 1082712d1f2SDave May ierr = DMSwarmGetField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 1092712d1f2SDave May 1102712d1f2SDave May de = DataExCreate(PetscObjectComm((PetscObject)dm),0);CHKERRQ(ierr); 1112712d1f2SDave May 1122712d1f2SDave May ierr = DataExTopologyInitialize(de);CHKERRQ(ierr); 1132712d1f2SDave May for (p=0; p<npoints; p++) { 1142712d1f2SDave May negrank = rankval[p]; 1152712d1f2SDave May if (negrank < 0) { 1162712d1f2SDave May nrank = -negrank - 1; 1172712d1f2SDave May ierr = DataExTopologyAddNeighbour(de,nrank);CHKERRQ(ierr); 1182712d1f2SDave May } 1192712d1f2SDave May } 1202712d1f2SDave May ierr = DataExTopologyFinalize(de);CHKERRQ(ierr); 1212712d1f2SDave May 1222712d1f2SDave May ierr = DataExInitializeSendCount(de);CHKERRQ(ierr); 1232712d1f2SDave May for (p=0; p<npoints; p++) { 1242712d1f2SDave May negrank = rankval[p]; 1252712d1f2SDave May if (negrank < 0) { 1262712d1f2SDave May nrank = -negrank - 1; 1272712d1f2SDave May ierr = DataExAddToSendCount(de,nrank,1);CHKERRQ(ierr); 1282712d1f2SDave May } 1292712d1f2SDave May } 1302712d1f2SDave May ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr); 1312712d1f2SDave May 1322712d1f2SDave May ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr); 1332712d1f2SDave May ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr); 1342712d1f2SDave May for (p=0; p<npoints; p++) { 1352712d1f2SDave May negrank = rankval[p]; 1362712d1f2SDave May if (negrank < 0) { 1372712d1f2SDave May nrank = -negrank - 1; 1382712d1f2SDave May rankval[p] = nrank; 1392712d1f2SDave May /* copy point into buffer */ 1402712d1f2SDave May ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr); 1412712d1f2SDave May /* insert point buffer into DataExchanger */ 1422712d1f2SDave May ierr = DataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr); 1432712d1f2SDave May rankval[p] = negrank; 1442712d1f2SDave May } 1452712d1f2SDave May } 1462712d1f2SDave May ierr = DataExPackFinalize(de);CHKERRQ(ierr); 1472712d1f2SDave May 1482712d1f2SDave May ierr = DMSwarmRestoreField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 1492712d1f2SDave May 1502712d1f2SDave May ierr = DataExBegin(de);CHKERRQ(ierr); 1512712d1f2SDave May ierr = DataExEnd(de);CHKERRQ(ierr); 1522712d1f2SDave May 1532712d1f2SDave May ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr); 1542712d1f2SDave May 1552712d1f2SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 1562712d1f2SDave May ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr); 1572712d1f2SDave May for (p=0; p<n_points_recv; p++) { 1582712d1f2SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point ); 1592712d1f2SDave May 1602712d1f2SDave May ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr); 1612712d1f2SDave May } 1622712d1f2SDave May ierr = DataExView(de);CHKERRQ(ierr); 1632712d1f2SDave May ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr); 1642712d1f2SDave May ierr = DataExDestroy(de);CHKERRQ(ierr); 1652712d1f2SDave May 1662712d1f2SDave May PetscFunctionReturn(0); 1672712d1f2SDave May } 168b16650c8SDave May 169b16650c8SDave May typedef struct { 170b16650c8SDave May PetscMPIInt owner_rank; 171b16650c8SDave May PetscReal min[3],max[3]; 172b16650c8SDave May } CollectBBox; 173b16650c8SDave May 174b16650c8SDave May #undef __FUNCT__ 175*fe39f135SDave May #define __FUNCT__ "DMSwarmCollect_DMDABoundingBox" 176*fe39f135SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM dm,PetscInt *globalsize) 177b16650c8SDave May { 178b16650c8SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 179b16650c8SDave May PetscErrorCode ierr; 180b16650c8SDave May DataEx de; 181b16650c8SDave May PetscInt p,pk,npoints,*rankval,n_points_recv,n_bbox_recv,dim,neighbour_cells; 182b16650c8SDave May PetscMPIInt rank,nrank; 183b16650c8SDave May void *point_buffer,*recv_points; 184b16650c8SDave May size_t sizeof_dmswarm_point,sizeof_bbox_ctx; 185b16650c8SDave May PetscBool isdmda; 186b16650c8SDave May CollectBBox *bbox,*recv_bbox; 187b16650c8SDave May const PetscMPIInt *dmneighborranks; 188b16650c8SDave May DM dmcell; 189b16650c8SDave May 190b16650c8SDave May ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 191b16650c8SDave May 192*fe39f135SDave May ierr = DMSwarmGetCellDM(dm,&dmcell);CHKERRQ(ierr); 193*fe39f135SDave May if (!dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided"); 194b16650c8SDave May 195b16650c8SDave May isdmda = PETSC_FALSE; 196b16650c8SDave May PetscObjectTypeCompare((PetscObject)dmcell,DMDA,&isdmda); 197b16650c8SDave May if (!isdmda) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only DMDA support for CollectBoundingBox"); 198b16650c8SDave May 199*fe39f135SDave May ierr = DMDAGetInfo(dm,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 200*fe39f135SDave May if (dim == 1) { 201*fe39f135SDave May neighbour_cells = 3; 202*fe39f135SDave May } else if (dim == 2) { 203*fe39f135SDave May neighbour_cells = 9; 204*fe39f135SDave May } else { 205*fe39f135SDave May neighbour_cells = 27; 206*fe39f135SDave May } 207*fe39f135SDave May 208b16650c8SDave May sizeof_bbox_ctx = sizeof(CollectBBox); 209b16650c8SDave May PetscMalloc1(1,&bbox); 210b16650c8SDave May bbox->owner_rank = rank; 211b16650c8SDave May 212b16650c8SDave May /* compute the bounding box based on the overlapping / stenctil size */ 213b16650c8SDave May //ierr = DMDAGetLocalBoundingBox(dmcell,bbox->min,bbox->max);CHKERRQ(ierr); 214b16650c8SDave May { 215b16650c8SDave May Vec lcoor; 216b16650c8SDave May 217b16650c8SDave May ierr = DMGetCoordinatesLocal(dmcell,&lcoor);CHKERRQ(ierr); 218b16650c8SDave May 219*fe39f135SDave May if (dim >= 1) { 220b16650c8SDave May ierr = VecStrideMin(lcoor,0,NULL,&bbox->min[0]);CHKERRQ(ierr); 221b16650c8SDave May ierr = VecStrideMax(lcoor,0,NULL,&bbox->max[0]);CHKERRQ(ierr); 222*fe39f135SDave May } 223*fe39f135SDave May 224*fe39f135SDave May if (dim >= 2) { 225*fe39f135SDave May ierr = VecStrideMin(lcoor,1,NULL,&bbox->min[1]);CHKERRQ(ierr); 226b16650c8SDave May ierr = VecStrideMax(lcoor,1,NULL,&bbox->max[1]);CHKERRQ(ierr); 227b16650c8SDave May } 228b16650c8SDave May 229*fe39f135SDave May if (dim == 3) { 230*fe39f135SDave May ierr = VecStrideMin(lcoor,2,NULL,&bbox->min[2]);CHKERRQ(ierr); 231*fe39f135SDave May ierr = VecStrideMax(lcoor,2,NULL,&bbox->max[2]);CHKERRQ(ierr); 232*fe39f135SDave May } 233*fe39f135SDave May } 234*fe39f135SDave May 235b16650c8SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 236b16650c8SDave May *globalsize = npoints; 237b16650c8SDave May ierr = DMSwarmGetField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 238b16650c8SDave May 239b16650c8SDave May de = DataExCreate(PetscObjectComm((PetscObject)dm),0);CHKERRQ(ierr); 240b16650c8SDave May 241b16650c8SDave May /* use DMDA neighbours */ 242b16650c8SDave May ierr = DMDAGetNeighbors(dmcell,&dmneighborranks);CHKERRQ(ierr); 243b16650c8SDave May 244b16650c8SDave May ierr = DataExTopologyInitialize(de);CHKERRQ(ierr); 245b16650c8SDave May for (p=0; p<neighbour_cells; p++) { 246b16650c8SDave May if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) { 247b16650c8SDave May ierr = DataExTopologyAddNeighbour(de,dmneighborranks[p]);CHKERRQ(ierr); 248b16650c8SDave May } 249b16650c8SDave May } 250b16650c8SDave May ierr = DataExTopologyFinalize(de);CHKERRQ(ierr); 251b16650c8SDave May 252b16650c8SDave May ierr = DataExInitializeSendCount(de);CHKERRQ(ierr); 253b16650c8SDave May for (p=0; p<neighbour_cells; p++) { 254b16650c8SDave May if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) { 255b16650c8SDave May ierr = DataExAddToSendCount(de,dmneighborranks[p],1);CHKERRQ(ierr); 256b16650c8SDave May } 257b16650c8SDave May } 258b16650c8SDave May ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr); 259b16650c8SDave May 260b16650c8SDave May 261b16650c8SDave May /* send bounding boxes */ 262b16650c8SDave May ierr = DataExPackInitialize(de,sizeof_bbox_ctx);CHKERRQ(ierr); 263b16650c8SDave May for (p=0; p<neighbour_cells; p++) { 264b16650c8SDave May nrank = dmneighborranks[p]; 265b16650c8SDave May if ( (nrank >= 0) && (nrank != rank) ) { 266b16650c8SDave May /* insert bbox buffer into DataExchanger */ 267b16650c8SDave May ierr = DataExPackData(de,nrank,1,bbox);CHKERRQ(ierr); 268b16650c8SDave May } 269b16650c8SDave May } 270b16650c8SDave May ierr = DataExPackFinalize(de);CHKERRQ(ierr); 271b16650c8SDave May 272b16650c8SDave May /* recv bounding boxes */ 273b16650c8SDave May ierr = DataExBegin(de);CHKERRQ(ierr); 274b16650c8SDave May ierr = DataExEnd(de);CHKERRQ(ierr); 275b16650c8SDave May 276b16650c8SDave May ierr = DataExGetRecvData(de,&n_bbox_recv,(void**)&recv_bbox);CHKERRQ(ierr); 277b16650c8SDave May 278b16650c8SDave May 279b16650c8SDave May for (p=0; p<n_bbox_recv; p++) { 280b16650c8SDave May printf("[rank %d]: box from %d : range[%+1.4e,%+1.4e]x[%+1.4e,%+1.4e]\n",rank,recv_bbox[p].owner_rank, 281b16650c8SDave May recv_bbox[p].min[0],recv_bbox[p].max[0],recv_bbox[p].min[1],recv_bbox[p].max[1]); 282b16650c8SDave May } 283b16650c8SDave May 284b16650c8SDave May /* of course this is stupid as this "generic" function should have a better way to know what the coordinates are called */ 285b16650c8SDave May ierr = DataExInitializeSendCount(de);CHKERRQ(ierr); 286b16650c8SDave May for (pk=0; pk<n_bbox_recv; pk++) { 287b16650c8SDave May PetscReal *array_x,*array_y; 288b16650c8SDave May 289b16650c8SDave May ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 290b16650c8SDave May ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 291b16650c8SDave May 292b16650c8SDave May for (p=0; p<npoints; p++) { 293b16650c8SDave May if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) { 294b16650c8SDave May if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) { 295b16650c8SDave May 296b16650c8SDave May ierr = DataExAddToSendCount(de,recv_bbox[pk].owner_rank,1);CHKERRQ(ierr); 297b16650c8SDave May 298b16650c8SDave May } 299b16650c8SDave May } 300b16650c8SDave May } 301b16650c8SDave May ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 302b16650c8SDave May ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 303b16650c8SDave May } 304b16650c8SDave May ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr); 305b16650c8SDave May 306b16650c8SDave May 307b16650c8SDave May 308b16650c8SDave May ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr); 309b16650c8SDave May ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr); 310b16650c8SDave May for (pk=0; pk<n_bbox_recv; pk++) { 311b16650c8SDave May PetscReal *array_x,*array_y; 312b16650c8SDave May 313b16650c8SDave May ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 314b16650c8SDave May ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 315b16650c8SDave May 316b16650c8SDave May for (p=0; p<npoints; p++) { 317b16650c8SDave May if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) { 318b16650c8SDave May if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) { 319b16650c8SDave May // copy point into buffer // 320b16650c8SDave May ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr); 321b16650c8SDave May // insert point buffer into DataExchanger // 322b16650c8SDave May ierr = DataExPackData(de,recv_bbox[pk].owner_rank,1,point_buffer);CHKERRQ(ierr); 323b16650c8SDave May } 324b16650c8SDave May } 325b16650c8SDave May } 326b16650c8SDave May 327b16650c8SDave May ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 328b16650c8SDave May ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 329b16650c8SDave May } 330b16650c8SDave May 331b16650c8SDave May ierr = DataExPackFinalize(de);CHKERRQ(ierr); 332b16650c8SDave May 333b16650c8SDave May ierr = DMSwarmRestoreField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 334b16650c8SDave May 335b16650c8SDave May ierr = DataExBegin(de);CHKERRQ(ierr); 336b16650c8SDave May ierr = DataExEnd(de);CHKERRQ(ierr); 337b16650c8SDave May 338b16650c8SDave May ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr); 339b16650c8SDave May 340b16650c8SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 341b16650c8SDave May ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr); 342b16650c8SDave May for (p=0; p<n_points_recv; p++) { 343b16650c8SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point ); 344b16650c8SDave May 345b16650c8SDave May ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr); 346b16650c8SDave May } 347b16650c8SDave May 348b16650c8SDave May ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr); 349b16650c8SDave May PetscFree(bbox); 350b16650c8SDave May ierr = DataExView(de);CHKERRQ(ierr); 351b16650c8SDave May ierr = DataExDestroy(de);CHKERRQ(ierr); 352b16650c8SDave May 353b16650c8SDave May PetscFunctionReturn(0); 354b16650c8SDave May } 355a9fd7477SDave May 356a9fd7477SDave May 357a9fd7477SDave May /* General collection when no order, or neighbour information is provided */ 358a9fd7477SDave May /* 359a9fd7477SDave May User provides context and collect() method 360a9fd7477SDave May Broadcast user context 361a9fd7477SDave May 362a9fd7477SDave May for each context / rank { 363a9fd7477SDave May collect(swarm,context,n,list) 364a9fd7477SDave May } 365a9fd7477SDave May 366a9fd7477SDave May */ 367a9fd7477SDave May #undef __FUNCT__ 368a9fd7477SDave May #define __FUNCT__ "DMSwarmCollect_General" 369a9fd7477SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM dm,PetscErrorCode (*collect)(DM,void*,PetscInt*,PetscInt**),size_t ctx_size,void *ctx,PetscInt *globalsize) 370a9fd7477SDave May { 371a9fd7477SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 372a9fd7477SDave May PetscErrorCode ierr; 373a9fd7477SDave May DataEx de; 374a9fd7477SDave May PetscInt p,r,npoints,n_points_recv; 375a9fd7477SDave May PetscMPIInt commsize,rank; 376a9fd7477SDave May void *point_buffer,*recv_points; 377a9fd7477SDave May void *ctxlist; 378a9fd7477SDave May PetscInt *n2collect,**collectlist; 379a9fd7477SDave May size_t sizeof_dmswarm_point; 380a9fd7477SDave May 381a9fd7477SDave May ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr); 382a9fd7477SDave May ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 383a9fd7477SDave May 384a9fd7477SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 385a9fd7477SDave May *globalsize = npoints; 386a9fd7477SDave May 387a9fd7477SDave May /* Broadcast user context */ 388a9fd7477SDave May PetscMalloc(ctx_size*commsize,&ctxlist); 389a9fd7477SDave May ierr = MPI_Allgather(ctx,ctx_size,MPI_CHAR,ctxlist,ctx_size,MPI_CHAR,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 390a9fd7477SDave May 391a9fd7477SDave May PetscMalloc1(commsize,&n2collect); 392a9fd7477SDave May PetscMalloc1(commsize,&collectlist); 393a9fd7477SDave May 394a9fd7477SDave May for (r=0; r<commsize; r++) { 395a9fd7477SDave May PetscInt _n2collect; 396a9fd7477SDave May PetscInt *_collectlist; 397a9fd7477SDave May void *_ctx_r; 398a9fd7477SDave May 399a9fd7477SDave May _n2collect = 0; 400a9fd7477SDave May _collectlist = NULL; 401a9fd7477SDave May 402a9fd7477SDave May if (r != rank) { /* don't collect data from yourself */ 403a9fd7477SDave May _ctx_r = (void*)( (char*)ctxlist + r * ctx_size ); 404a9fd7477SDave May ierr = collect(dm,_ctx_r,&_n2collect,&_collectlist);CHKERRQ(ierr); 405a9fd7477SDave May } 406a9fd7477SDave May 407a9fd7477SDave May n2collect[r] = _n2collect; 408a9fd7477SDave May collectlist[r] = _collectlist; 409a9fd7477SDave May } 410a9fd7477SDave May 411a9fd7477SDave May de = DataExCreate(PetscObjectComm((PetscObject)dm),0);CHKERRQ(ierr); 412a9fd7477SDave May 413a9fd7477SDave May /* Define topology */ 414a9fd7477SDave May ierr = DataExTopologyInitialize(de);CHKERRQ(ierr); 415a9fd7477SDave May for (r=0; r<commsize; r++) { 416a9fd7477SDave May if (n2collect[r] > 0) { 417a9fd7477SDave May ierr = DataExTopologyAddNeighbour(de,(PetscMPIInt)r);CHKERRQ(ierr); 418a9fd7477SDave May } 419a9fd7477SDave May } 420a9fd7477SDave May ierr = DataExTopologyFinalize(de);CHKERRQ(ierr); 421a9fd7477SDave May 422a9fd7477SDave May /* Define send counts */ 423a9fd7477SDave May ierr = DataExInitializeSendCount(de);CHKERRQ(ierr); 424a9fd7477SDave May for (r=0; r<commsize; r++) { 425a9fd7477SDave May if (n2collect[r] > 0) { 426a9fd7477SDave May ierr = DataExAddToSendCount(de,r,n2collect[r]);CHKERRQ(ierr); 427a9fd7477SDave May } 428a9fd7477SDave May } 429a9fd7477SDave May ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr); 430a9fd7477SDave May 431a9fd7477SDave May /* Pack data */ 432a9fd7477SDave May ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr); 433a9fd7477SDave May 434a9fd7477SDave May ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr); 435a9fd7477SDave May for (r=0; r<commsize; r++) { 436a9fd7477SDave May 437a9fd7477SDave May for (p=0; p<n2collect[r]; p++) { 438a9fd7477SDave May ierr = DataBucketFillPackedArray(swarm->db,collectlist[r][p],point_buffer);CHKERRQ(ierr); 439a9fd7477SDave May /* insert point buffer into the data exchanger */ 440a9fd7477SDave May ierr = DataExPackData(de,r,1,point_buffer);CHKERRQ(ierr); 441a9fd7477SDave May } 442a9fd7477SDave May } 443a9fd7477SDave May 444a9fd7477SDave May ierr = DataExPackFinalize(de);CHKERRQ(ierr); 445a9fd7477SDave May 446a9fd7477SDave May /* Scatter */ 447a9fd7477SDave May ierr = DataExBegin(de);CHKERRQ(ierr); 448a9fd7477SDave May ierr = DataExEnd(de);CHKERRQ(ierr); 449a9fd7477SDave May 450a9fd7477SDave May /* Collect data in DMSwarm container */ 451a9fd7477SDave May ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr); 452a9fd7477SDave May 453a9fd7477SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 454a9fd7477SDave May ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr); 455a9fd7477SDave May for (p=0; p<n_points_recv; p++) { 456a9fd7477SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point ); 457a9fd7477SDave May 458a9fd7477SDave May ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr); 459a9fd7477SDave May } 460a9fd7477SDave May 461a9fd7477SDave May /* Release memory */ 462a9fd7477SDave May for (r=0; r<commsize; r++) { 463a9fd7477SDave May if (collectlist[r]) PetscFree(collectlist[r]); 464a9fd7477SDave May } 465a9fd7477SDave May PetscFree(collectlist); 466a9fd7477SDave May PetscFree(n2collect); 467a9fd7477SDave May PetscFree(ctxlist); 468a9fd7477SDave May ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr); 469a9fd7477SDave May ierr = DataExView(de);CHKERRQ(ierr); 470a9fd7477SDave May ierr = DataExDestroy(de);CHKERRQ(ierr); 471a9fd7477SDave May 472a9fd7477SDave May PetscFunctionReturn(0); 473a9fd7477SDave May } 474a9fd7477SDave May 475