1df21e3a8SDave May 208056efcSDave May #include <petscdmswarm.h> 3df21e3a8SDave May #include <petsc/private/dmswarmimpl.h> /*I "petscdmswarm.h" I*/ 4df21e3a8SDave May #include "data_bucket.h" 5df21e3a8SDave May #include "data_ex.h" 6df21e3a8SDave May 7df21e3a8SDave May 8480eef7bSDave May /* 9480eef7bSDave May User loads desired location (MPI rank) into field DMSwarm_rank 10480eef7bSDave May */ 11df21e3a8SDave May #undef __FUNCT__ 12df21e3a8SDave May #define __FUNCT__ "DMSwarmMigrate_Push_Basic" 13df21e3a8SDave May PetscErrorCode DMSwarmMigrate_Push_Basic(DM dm,PetscBool remove_sent_points) 14df21e3a8SDave May { 15df21e3a8SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 16df21e3a8SDave May PetscErrorCode ierr; 17df21e3a8SDave May DataEx de; 18df21e3a8SDave May PetscInt p,npoints,*rankval,n_points_recv; 19df21e3a8SDave May PetscMPIInt rank,nrank; 20df21e3a8SDave May void *point_buffer,*recv_points; 21df21e3a8SDave May size_t sizeof_dmswarm_point; 22df21e3a8SDave May 23df21e3a8SDave May ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 24df21e3a8SDave May 25df21e3a8SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 2608056efcSDave May ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 27df21e3a8SDave May 28df21e3a8SDave May de = DataExCreate(PetscObjectComm((PetscObject)dm),0);CHKERRQ(ierr); 29df21e3a8SDave May 30df21e3a8SDave May ierr = DataExTopologyInitialize(de);CHKERRQ(ierr); 31df21e3a8SDave May for (p=0; p<npoints; p++) { 32df21e3a8SDave May nrank = rankval[p]; 33df21e3a8SDave May if (nrank != rank) { 34df21e3a8SDave May ierr = DataExTopologyAddNeighbour(de,nrank);CHKERRQ(ierr); 35df21e3a8SDave May } 36df21e3a8SDave May } 37df21e3a8SDave May ierr = DataExTopologyFinalize(de);CHKERRQ(ierr); 38df21e3a8SDave May 39df21e3a8SDave May ierr = DataExInitializeSendCount(de);CHKERRQ(ierr); 40df21e3a8SDave May for (p=0; p<npoints; p++) { 41df21e3a8SDave May nrank = rankval[p]; 42df21e3a8SDave May if (nrank != rank) { 43df21e3a8SDave May ierr = DataExAddToSendCount(de,nrank,1);CHKERRQ(ierr); 44df21e3a8SDave May } 45df21e3a8SDave May } 46df21e3a8SDave May ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr); 47df21e3a8SDave May 48df21e3a8SDave May ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr); 49df21e3a8SDave May ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr); 50df21e3a8SDave May for (p=0; p<npoints; p++) { 51df21e3a8SDave May nrank = rankval[p]; 52df21e3a8SDave May if (nrank != rank) { 53df21e3a8SDave May /* copy point into buffer */ 54df21e3a8SDave May ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr); 55df21e3a8SDave May /* insert point buffer into DataExchanger */ 56df21e3a8SDave May ierr = DataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr); 57df21e3a8SDave May } 58df21e3a8SDave May } 59df21e3a8SDave May ierr = DataExPackFinalize(de);CHKERRQ(ierr); 60df21e3a8SDave May 61df21e3a8SDave May 62df21e3a8SDave May if (remove_sent_points) { 63df21e3a8SDave May /* remove points which left processor */ 64df21e3a8SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 65df21e3a8SDave May for (p=0; p<npoints; p++) { 66df21e3a8SDave May nrank = rankval[p]; 67df21e3a8SDave May if (nrank != rank) { 68df21e3a8SDave May /* kill point */ 69df21e3a8SDave May ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr); 70df21e3a8SDave May DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */ 71df21e3a8SDave May p--; /* check replacement point */ 72df21e3a8SDave May } 73df21e3a8SDave May } 74df21e3a8SDave May } 7508056efcSDave May ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 76df21e3a8SDave May 77df21e3a8SDave May ierr = DataExBegin(de);CHKERRQ(ierr); 78df21e3a8SDave May ierr = DataExEnd(de);CHKERRQ(ierr); 79df21e3a8SDave May 80df21e3a8SDave May ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr); 81df21e3a8SDave May 82df21e3a8SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 83df21e3a8SDave May ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr); 84df21e3a8SDave May for (p=0; p<n_points_recv; p++) { 85df21e3a8SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point ); 86df21e3a8SDave May 87df21e3a8SDave May ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr); 88df21e3a8SDave May } 89df21e3a8SDave May ierr = DataExView(de);CHKERRQ(ierr); 90df21e3a8SDave May ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr); 91df21e3a8SDave May ierr = DataExDestroy(de);CHKERRQ(ierr); 92df21e3a8SDave May 93df21e3a8SDave May PetscFunctionReturn(0); 94df21e3a8SDave May } 952712d1f2SDave May 9640c453e9SDave May #undef __FUNCT__ 97889dbfe5SDave May #define __FUNCT__ "DMSwarmMigrate_DMNeighborScatter" 98889dbfe5SDave May PetscErrorCode DMSwarmMigrate_DMNeighborScatter(DM dm,DM dmcell,PetscBool remove_sent_points,PetscInt *npoints_prior_migration) 9940c453e9SDave May { 10040c453e9SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 10140c453e9SDave May PetscErrorCode ierr; 10240c453e9SDave May DataEx de; 10340c453e9SDave May PetscInt r,p,npoints,*rankval,n_points_recv; 10440c453e9SDave May PetscMPIInt rank,_rank; 10540c453e9SDave May const PetscMPIInt *neighbourranks; 10640c453e9SDave May void *point_buffer,*recv_points; 10740c453e9SDave May size_t sizeof_dmswarm_point; 10840c453e9SDave May PetscInt nneighbors; 109*7c6d1d28SDave May PetscMPIInt mynneigh,*myneigh; 11040c453e9SDave May 11140c453e9SDave May ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 11240c453e9SDave May 11340c453e9SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 11408056efcSDave May ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 11540c453e9SDave May 11640c453e9SDave May de = DataExCreate(PetscObjectComm((PetscObject)dm),0);CHKERRQ(ierr); 11740c453e9SDave May 11840c453e9SDave May ierr = DMGetNeighbors(dmcell,&nneighbors,&neighbourranks);CHKERRQ(ierr); 11940c453e9SDave May ierr = DataExTopologyInitialize(de);CHKERRQ(ierr); 12040c453e9SDave May for (r=0; r<nneighbors; r++) { 12140c453e9SDave May _rank = neighbourranks[r]; 12240c453e9SDave May if ((_rank != rank) && (_rank > 0)) { 12340c453e9SDave May ierr = DataExTopologyAddNeighbour(de,_rank);CHKERRQ(ierr); 12440c453e9SDave May } 12540c453e9SDave May } 12640c453e9SDave May ierr = DataExTopologyFinalize(de);CHKERRQ(ierr); 127*7c6d1d28SDave May ierr = DataExTopologyGetNeighbours(de,&mynneigh,&myneigh);CHKERRQ(ierr); 12840c453e9SDave May 12940c453e9SDave May ierr = DataExInitializeSendCount(de);CHKERRQ(ierr); 13040c453e9SDave May for (p=0; p<npoints; p++) { 13140c453e9SDave May if (rankval[p] == -1) { 132*7c6d1d28SDave May for (r=0; r<mynneigh; r++) { 133*7c6d1d28SDave May _rank = myneigh[r]; 13440c453e9SDave May ierr = DataExAddToSendCount(de,_rank,1);CHKERRQ(ierr); 13540c453e9SDave May } 13640c453e9SDave May } 13740c453e9SDave May } 13840c453e9SDave May ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr); 13940c453e9SDave May 14040c453e9SDave May ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr); 14140c453e9SDave May ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr); 14240c453e9SDave May for (p=0; p<npoints; p++) { 14340c453e9SDave May if (rankval[p] == -1) { 144*7c6d1d28SDave May for (r=0; r<mynneigh; r++) { 145*7c6d1d28SDave May _rank = myneigh[r]; 14640c453e9SDave May /* copy point into buffer */ 14740c453e9SDave May ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr); 14840c453e9SDave May /* insert point buffer into DataExchanger */ 14940c453e9SDave May ierr = DataExPackData(de,_rank,1,point_buffer);CHKERRQ(ierr); 15040c453e9SDave May } 15140c453e9SDave May } 15240c453e9SDave May } 15340c453e9SDave May ierr = DataExPackFinalize(de);CHKERRQ(ierr); 154*7c6d1d28SDave May ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 15540c453e9SDave May 15640c453e9SDave May if (remove_sent_points) { 157*7c6d1d28SDave May DataField PField; 158*7c6d1d28SDave May 159*7c6d1d28SDave May ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr); 160*7c6d1d28SDave May ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); 161*7c6d1d28SDave May 16240c453e9SDave May /* remove points which left processor */ 16340c453e9SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 16440c453e9SDave May for (p=0; p<npoints; p++) { 16540c453e9SDave May if (rankval[p] == -1) { 16640c453e9SDave May /* kill point */ 16740c453e9SDave May ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr); 168*7c6d1d28SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */ 169*7c6d1d28SDave May ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */ 17040c453e9SDave May p--; /* check replacement point */ 17140c453e9SDave May } 17240c453e9SDave May } 17340c453e9SDave May } 17440c453e9SDave May ierr = DataBucketGetSizes(swarm->db,npoints_prior_migration,NULL,NULL);CHKERRQ(ierr); 17540c453e9SDave May 17640c453e9SDave May ierr = DataExBegin(de);CHKERRQ(ierr); 17740c453e9SDave May ierr = DataExEnd(de);CHKERRQ(ierr); 17840c453e9SDave May 17940c453e9SDave May ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr); 18040c453e9SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 18140c453e9SDave May ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr); 18240c453e9SDave May for (p=0; p<n_points_recv; p++) { 18340c453e9SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point ); 18440c453e9SDave May 18540c453e9SDave May ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr); 18640c453e9SDave May } 187*7c6d1d28SDave May 188*7c6d1d28SDave May //ierr = DataExView(de);CHKERRQ(ierr); 18940c453e9SDave May ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr); 19040c453e9SDave May ierr = DataExDestroy(de);CHKERRQ(ierr); 19140c453e9SDave May 19240c453e9SDave May PetscFunctionReturn(0); 19340c453e9SDave May } 194480eef7bSDave May 195480eef7bSDave May #undef __FUNCT__ 19608056efcSDave May #define __FUNCT__ "DMSwarmMigrate_CellDMScatter" 19708056efcSDave May PetscErrorCode DMSwarmMigrate_CellDMScatter(DM dm,PetscBool remove_sent_points) 198480eef7bSDave May { 199480eef7bSDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 200480eef7bSDave May PetscErrorCode ierr; 201*7c6d1d28SDave May PetscInt p,npoints,npointsg=0,npoints2,npoints2g,len,*rankval,npoints_prior_migration; 202480eef7bSDave May const PetscInt *LA_iscell; 203480eef7bSDave May DM dmcell; 204480eef7bSDave May IS iscell; 205480eef7bSDave May Vec pos; 20640c453e9SDave May PetscBool error_check = swarm->migrate_error_on_missing_point; 207*7c6d1d28SDave May PetscMPIInt commsize,rank; 208480eef7bSDave May 209480eef7bSDave May ierr = DMSwarmGetCellDM(dm,&dmcell);CHKERRQ(ierr); 210480eef7bSDave May if (!dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided"); 211480eef7bSDave May 212*7c6d1d28SDave May ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr); 213*7c6d1d28SDave May ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 214*7c6d1d28SDave May if (commsize == 1) PetscFunctionReturn(0); 215*7c6d1d28SDave May 21608056efcSDave May ierr = DMSwarmCreateGlobalVectorFromField(dm,DMSwarmPICField_coor,&pos);CHKERRQ(ierr); 217480eef7bSDave May ierr = DMLocatePoints(dmcell,pos,&iscell);CHKERRQ(ierr); 21808056efcSDave May ierr = DMSwarmDestroyGlobalVectorFromField(dm,DMSwarmPICField_coor,&pos);CHKERRQ(ierr); 219480eef7bSDave May 22040c453e9SDave May if (error_check) { 22140c453e9SDave May ierr = DMSwarmGetSize(dm,&npointsg);CHKERRQ(ierr); 22240c453e9SDave May } 223480eef7bSDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 22408056efcSDave May ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 225480eef7bSDave May ierr = ISGetIndices(iscell,&LA_iscell);CHKERRQ(ierr); 226480eef7bSDave May for (p=0; p<npoints; p++) { 227*7c6d1d28SDave May rankval[p] = LA_iscell[p]; 228480eef7bSDave May } 229480eef7bSDave May ierr = ISRestoreIndices(iscell,&LA_iscell);CHKERRQ(ierr); 23008056efcSDave May ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 231480eef7bSDave May ierr = ISDestroy(&iscell);CHKERRQ(ierr); 232480eef7bSDave May 233889dbfe5SDave May ierr = DMSwarmMigrate_DMNeighborScatter(dm,dmcell,remove_sent_points,&npoints_prior_migration);CHKERRQ(ierr); 234480eef7bSDave May 23540c453e9SDave May /* locate points newly recevied */ 23640c453e9SDave May ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); 237*7c6d1d28SDave May #if 0 23840c453e9SDave May len = npoints2 - npoints_prior_migration; 23940c453e9SDave May if (len > 0) { 24040c453e9SDave May PetscScalar *LA_coor; 241*7c6d1d28SDave May PetscInt bs; 24240c453e9SDave May 243*7c6d1d28SDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr); 244*7c6d1d28SDave May ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*len,(const PetscScalar*)&LA_coor[bs*npoints_prior_migration],&pos);CHKERRQ(ierr); 24540c453e9SDave May 24640c453e9SDave May ierr = DMLocatePoints(dmcell,pos,&iscell);CHKERRQ(ierr); 24740c453e9SDave May 24840c453e9SDave May ierr = VecDestroy(&pos);CHKERRQ(ierr); 249*7c6d1d28SDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr); 25040c453e9SDave May 25140c453e9SDave May ierr = ISGetIndices(iscell,&LA_iscell);CHKERRQ(ierr); 25208056efcSDave May ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 25340c453e9SDave May for (p=0; p<len; p++) { 25440c453e9SDave May if (LA_iscell[p] == -1) { 25540c453e9SDave May rankval[npoints_prior_migration+p] = -1; 25640c453e9SDave May } 25740c453e9SDave May } 25840c453e9SDave May ierr = ISRestoreIndices(iscell,&LA_iscell);CHKERRQ(ierr); 25940c453e9SDave May ierr = ISDestroy(&iscell);CHKERRQ(ierr); 260*7c6d1d28SDave May ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 26140c453e9SDave May 26240c453e9SDave May /* remove points which left processor */ 263*7c6d1d28SDave May ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr); 264*7c6d1d28SDave May ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); 265*7c6d1d28SDave May 26640c453e9SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 26740c453e9SDave May for (p=npoints_prior_migration; p<npoints2; p++) { 26840c453e9SDave May if (rankval[p] == -1) { 26940c453e9SDave May /* kill point */ 27040c453e9SDave May ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr); 271*7c6d1d28SDave May ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */ 272*7c6d1d28SDave May ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */ 27340c453e9SDave May p--; /* check replacement point */ 27440c453e9SDave May } 27540c453e9SDave May } 276*7c6d1d28SDave May 277*7c6d1d28SDave May } 278*7c6d1d28SDave May #endif 279*7c6d1d28SDave May 280*7c6d1d28SDave May { 281*7c6d1d28SDave May PetscScalar *LA_coor; 282*7c6d1d28SDave May PetscInt bs; 283*7c6d1d28SDave May DataField PField; 284*7c6d1d28SDave May 285*7c6d1d28SDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr); 286*7c6d1d28SDave May ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*npoints2,(const PetscScalar*)LA_coor,&pos);CHKERRQ(ierr); 287*7c6d1d28SDave May ierr = DMLocatePoints(dmcell,pos,&iscell);CHKERRQ(ierr); 288*7c6d1d28SDave May 289*7c6d1d28SDave May ierr = VecDestroy(&pos);CHKERRQ(ierr); 290*7c6d1d28SDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr); 291*7c6d1d28SDave May 292*7c6d1d28SDave May ierr = ISGetIndices(iscell,&LA_iscell);CHKERRQ(ierr); 293*7c6d1d28SDave May ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 294*7c6d1d28SDave May for (p=0; p<npoints2; p++) { 295*7c6d1d28SDave May rankval[p] = LA_iscell[p]; 296*7c6d1d28SDave May } 297*7c6d1d28SDave May ierr = ISRestoreIndices(iscell,&LA_iscell);CHKERRQ(ierr); 298*7c6d1d28SDave May ierr = ISDestroy(&iscell);CHKERRQ(ierr); 29908056efcSDave May ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 30040c453e9SDave May 301*7c6d1d28SDave May /* remove points which left processor */ 302*7c6d1d28SDave May ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr); 303*7c6d1d28SDave May ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); 304*7c6d1d28SDave May 305*7c6d1d28SDave May ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); 306*7c6d1d28SDave May for (p=0; p<npoints2; p++) { 307*7c6d1d28SDave May if (rankval[p] == -1) { 308*7c6d1d28SDave May /* kill point */ 309*7c6d1d28SDave May ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr); 310*7c6d1d28SDave May ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */ 311*7c6d1d28SDave May ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */ 312*7c6d1d28SDave May p--; /* check replacement point */ 313*7c6d1d28SDave May } 314*7c6d1d28SDave May } 315*7c6d1d28SDave May 31640c453e9SDave May } 31740c453e9SDave May 31840c453e9SDave May /* check for error on removed points */ 31940c453e9SDave May if (error_check) { 32040c453e9SDave May ierr = DMSwarmGetSize(dm,&npoints2g);CHKERRQ(ierr); 32140c453e9SDave May 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); 32240c453e9SDave May } 323480eef7bSDave May PetscFunctionReturn(0); 324480eef7bSDave May } 325480eef7bSDave May 32608056efcSDave May #undef __FUNCT__ 32708056efcSDave May #define __FUNCT__ "DMSwarmMigrate_CellDMExact" 32808056efcSDave May PetscErrorCode DMSwarmMigrate_CellDMExact(DM dm,PetscBool remove_sent_points) 32908056efcSDave May { 33008056efcSDave May // DM_Swarm *swarm = (DM_Swarm*)dm->data; 33108056efcSDave May // PetscErrorCode ierr; 33208056efcSDave May 33308056efcSDave May PetscFunctionReturn(0); 33408056efcSDave May } 33508056efcSDave May 336480eef7bSDave May /* 337480eef7bSDave May Redundant as this assumes points can only be sent to a single rank 338480eef7bSDave May */ 3392712d1f2SDave May #undef __FUNCT__ 3402712d1f2SDave May #define __FUNCT__ "DMSwarmMigrate_GlobalToLocal_Basic" 3412712d1f2SDave May PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt *globalsize) 3422712d1f2SDave May { 3432712d1f2SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 3442712d1f2SDave May PetscErrorCode ierr; 3452712d1f2SDave May DataEx de; 3462712d1f2SDave May PetscInt p,npoints,*rankval,n_points_recv; 3472712d1f2SDave May PetscMPIInt rank,nrank,negrank; 3482712d1f2SDave May void *point_buffer,*recv_points; 3492712d1f2SDave May size_t sizeof_dmswarm_point; 3502712d1f2SDave May 3512712d1f2SDave May ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 3522712d1f2SDave May 3532712d1f2SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 3542712d1f2SDave May *globalsize = npoints; 35508056efcSDave May ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 3562712d1f2SDave May 3572712d1f2SDave May de = DataExCreate(PetscObjectComm((PetscObject)dm),0);CHKERRQ(ierr); 3582712d1f2SDave May 3592712d1f2SDave May ierr = DataExTopologyInitialize(de);CHKERRQ(ierr); 3602712d1f2SDave May for (p=0; p<npoints; p++) { 3612712d1f2SDave May negrank = rankval[p]; 3622712d1f2SDave May if (negrank < 0) { 3632712d1f2SDave May nrank = -negrank - 1; 3642712d1f2SDave May ierr = DataExTopologyAddNeighbour(de,nrank);CHKERRQ(ierr); 3652712d1f2SDave May } 3662712d1f2SDave May } 3672712d1f2SDave May ierr = DataExTopologyFinalize(de);CHKERRQ(ierr); 3682712d1f2SDave May 3692712d1f2SDave May ierr = DataExInitializeSendCount(de);CHKERRQ(ierr); 3702712d1f2SDave May for (p=0; p<npoints; p++) { 3712712d1f2SDave May negrank = rankval[p]; 3722712d1f2SDave May if (negrank < 0) { 3732712d1f2SDave May nrank = -negrank - 1; 3742712d1f2SDave May ierr = DataExAddToSendCount(de,nrank,1);CHKERRQ(ierr); 3752712d1f2SDave May } 3762712d1f2SDave May } 3772712d1f2SDave May ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr); 3782712d1f2SDave May 3792712d1f2SDave May ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr); 3802712d1f2SDave May ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr); 3812712d1f2SDave May for (p=0; p<npoints; p++) { 3822712d1f2SDave May negrank = rankval[p]; 3832712d1f2SDave May if (negrank < 0) { 3842712d1f2SDave May nrank = -negrank - 1; 3852712d1f2SDave May rankval[p] = nrank; 3862712d1f2SDave May /* copy point into buffer */ 3872712d1f2SDave May ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr); 3882712d1f2SDave May /* insert point buffer into DataExchanger */ 3892712d1f2SDave May ierr = DataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr); 3902712d1f2SDave May rankval[p] = negrank; 3912712d1f2SDave May } 3922712d1f2SDave May } 3932712d1f2SDave May ierr = DataExPackFinalize(de);CHKERRQ(ierr); 3942712d1f2SDave May 39508056efcSDave May ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 3962712d1f2SDave May 3972712d1f2SDave May ierr = DataExBegin(de);CHKERRQ(ierr); 3982712d1f2SDave May ierr = DataExEnd(de);CHKERRQ(ierr); 3992712d1f2SDave May 4002712d1f2SDave May ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr); 4012712d1f2SDave May 4022712d1f2SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 4032712d1f2SDave May ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr); 4042712d1f2SDave May for (p=0; p<n_points_recv; p++) { 4052712d1f2SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point ); 4062712d1f2SDave May 4072712d1f2SDave May ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr); 4082712d1f2SDave May } 4092712d1f2SDave May ierr = DataExView(de);CHKERRQ(ierr); 4102712d1f2SDave May ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr); 4112712d1f2SDave May ierr = DataExDestroy(de);CHKERRQ(ierr); 4122712d1f2SDave May 4132712d1f2SDave May PetscFunctionReturn(0); 4142712d1f2SDave May } 415b16650c8SDave May 416b16650c8SDave May typedef struct { 417b16650c8SDave May PetscMPIInt owner_rank; 418b16650c8SDave May PetscReal min[3],max[3]; 419b16650c8SDave May } CollectBBox; 420b16650c8SDave May 421b16650c8SDave May #undef __FUNCT__ 422fe39f135SDave May #define __FUNCT__ "DMSwarmCollect_DMDABoundingBox" 423fe39f135SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM dm,PetscInt *globalsize) 424b16650c8SDave May { 425b16650c8SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 426b16650c8SDave May PetscErrorCode ierr; 427b16650c8SDave May DataEx de; 428b16650c8SDave May PetscInt p,pk,npoints,*rankval,n_points_recv,n_bbox_recv,dim,neighbour_cells; 429b16650c8SDave May PetscMPIInt rank,nrank; 430b16650c8SDave May void *point_buffer,*recv_points; 431b16650c8SDave May size_t sizeof_dmswarm_point,sizeof_bbox_ctx; 432b16650c8SDave May PetscBool isdmda; 433b16650c8SDave May CollectBBox *bbox,*recv_bbox; 434b16650c8SDave May const PetscMPIInt *dmneighborranks; 435b16650c8SDave May DM dmcell; 436b16650c8SDave May 437b16650c8SDave May ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 438b16650c8SDave May 439fe39f135SDave May ierr = DMSwarmGetCellDM(dm,&dmcell);CHKERRQ(ierr); 440fe39f135SDave May if (!dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided"); 441b16650c8SDave May 442b16650c8SDave May isdmda = PETSC_FALSE; 443b16650c8SDave May PetscObjectTypeCompare((PetscObject)dmcell,DMDA,&isdmda); 444b16650c8SDave May if (!isdmda) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only DMDA support for CollectBoundingBox"); 445b16650c8SDave May 4468dbd68bcSDave May ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 447fe39f135SDave May 448b16650c8SDave May sizeof_bbox_ctx = sizeof(CollectBBox); 449b16650c8SDave May PetscMalloc1(1,&bbox); 450b16650c8SDave May bbox->owner_rank = rank; 451b16650c8SDave May 452b16650c8SDave May /* compute the bounding box based on the overlapping / stenctil size */ 453b16650c8SDave May //ierr = DMDAGetLocalBoundingBox(dmcell,bbox->min,bbox->max);CHKERRQ(ierr); 454b16650c8SDave May { 455b16650c8SDave May Vec lcoor; 456b16650c8SDave May 457b16650c8SDave May ierr = DMGetCoordinatesLocal(dmcell,&lcoor);CHKERRQ(ierr); 458b16650c8SDave May 459fe39f135SDave May if (dim >= 1) { 460b16650c8SDave May ierr = VecStrideMin(lcoor,0,NULL,&bbox->min[0]);CHKERRQ(ierr); 461b16650c8SDave May ierr = VecStrideMax(lcoor,0,NULL,&bbox->max[0]);CHKERRQ(ierr); 462fe39f135SDave May } 463fe39f135SDave May 464fe39f135SDave May if (dim >= 2) { 465fe39f135SDave May ierr = VecStrideMin(lcoor,1,NULL,&bbox->min[1]);CHKERRQ(ierr); 466b16650c8SDave May ierr = VecStrideMax(lcoor,1,NULL,&bbox->max[1]);CHKERRQ(ierr); 467b16650c8SDave May } 468b16650c8SDave May 469fe39f135SDave May if (dim == 3) { 470fe39f135SDave May ierr = VecStrideMin(lcoor,2,NULL,&bbox->min[2]);CHKERRQ(ierr); 471fe39f135SDave May ierr = VecStrideMax(lcoor,2,NULL,&bbox->max[2]);CHKERRQ(ierr); 472fe39f135SDave May } 473fe39f135SDave May } 474fe39f135SDave May 475b16650c8SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 476b16650c8SDave May *globalsize = npoints; 47708056efcSDave May ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 478b16650c8SDave May 479b16650c8SDave May de = DataExCreate(PetscObjectComm((PetscObject)dm),0);CHKERRQ(ierr); 480b16650c8SDave May 481b16650c8SDave May /* use DMDA neighbours */ 482b16650c8SDave May ierr = DMDAGetNeighbors(dmcell,&dmneighborranks);CHKERRQ(ierr); 4838dbd68bcSDave May if (dim == 1) { 4848dbd68bcSDave May neighbour_cells = 3; 4858dbd68bcSDave May } else if (dim == 2) { 4868dbd68bcSDave May neighbour_cells = 9; 4878dbd68bcSDave May } else { 4888dbd68bcSDave May neighbour_cells = 27; 4898dbd68bcSDave May } 490b16650c8SDave May 491b16650c8SDave May ierr = DataExTopologyInitialize(de);CHKERRQ(ierr); 492b16650c8SDave May for (p=0; p<neighbour_cells; p++) { 493b16650c8SDave May if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) { 494b16650c8SDave May ierr = DataExTopologyAddNeighbour(de,dmneighborranks[p]);CHKERRQ(ierr); 495b16650c8SDave May } 496b16650c8SDave May } 497b16650c8SDave May ierr = DataExTopologyFinalize(de);CHKERRQ(ierr); 498b16650c8SDave May 499b16650c8SDave May ierr = DataExInitializeSendCount(de);CHKERRQ(ierr); 500b16650c8SDave May for (p=0; p<neighbour_cells; p++) { 501b16650c8SDave May if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) { 502b16650c8SDave May ierr = DataExAddToSendCount(de,dmneighborranks[p],1);CHKERRQ(ierr); 503b16650c8SDave May } 504b16650c8SDave May } 505b16650c8SDave May ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr); 506b16650c8SDave May 507b16650c8SDave May 508b16650c8SDave May /* send bounding boxes */ 509b16650c8SDave May ierr = DataExPackInitialize(de,sizeof_bbox_ctx);CHKERRQ(ierr); 510b16650c8SDave May for (p=0; p<neighbour_cells; p++) { 511b16650c8SDave May nrank = dmneighborranks[p]; 512b16650c8SDave May if ( (nrank >= 0) && (nrank != rank) ) { 513b16650c8SDave May /* insert bbox buffer into DataExchanger */ 514b16650c8SDave May ierr = DataExPackData(de,nrank,1,bbox);CHKERRQ(ierr); 515b16650c8SDave May } 516b16650c8SDave May } 517b16650c8SDave May ierr = DataExPackFinalize(de);CHKERRQ(ierr); 518b16650c8SDave May 519b16650c8SDave May /* recv bounding boxes */ 520b16650c8SDave May ierr = DataExBegin(de);CHKERRQ(ierr); 521b16650c8SDave May ierr = DataExEnd(de);CHKERRQ(ierr); 522b16650c8SDave May 523b16650c8SDave May ierr = DataExGetRecvData(de,&n_bbox_recv,(void**)&recv_bbox);CHKERRQ(ierr); 524b16650c8SDave May 525b16650c8SDave May 526b16650c8SDave May for (p=0; p<n_bbox_recv; p++) { 527b16650c8SDave May printf("[rank %d]: box from %d : range[%+1.4e,%+1.4e]x[%+1.4e,%+1.4e]\n",rank,recv_bbox[p].owner_rank, 528b16650c8SDave May recv_bbox[p].min[0],recv_bbox[p].max[0],recv_bbox[p].min[1],recv_bbox[p].max[1]); 529b16650c8SDave May } 530b16650c8SDave May 531b16650c8SDave May /* of course this is stupid as this "generic" function should have a better way to know what the coordinates are called */ 532b16650c8SDave May ierr = DataExInitializeSendCount(de);CHKERRQ(ierr); 533b16650c8SDave May for (pk=0; pk<n_bbox_recv; pk++) { 534b16650c8SDave May PetscReal *array_x,*array_y; 535b16650c8SDave May 536b16650c8SDave May ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 537b16650c8SDave May ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 538b16650c8SDave May 539b16650c8SDave May for (p=0; p<npoints; p++) { 540b16650c8SDave May if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) { 541b16650c8SDave May if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) { 542b16650c8SDave May 543b16650c8SDave May ierr = DataExAddToSendCount(de,recv_bbox[pk].owner_rank,1);CHKERRQ(ierr); 544b16650c8SDave May 545b16650c8SDave May } 546b16650c8SDave May } 547b16650c8SDave May } 548b16650c8SDave May ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 549b16650c8SDave May ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 550b16650c8SDave May } 551b16650c8SDave May ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr); 552b16650c8SDave May 553b16650c8SDave May 554b16650c8SDave May 555b16650c8SDave May ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr); 556b16650c8SDave May ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr); 557b16650c8SDave May for (pk=0; pk<n_bbox_recv; pk++) { 558b16650c8SDave May PetscReal *array_x,*array_y; 559b16650c8SDave May 560b16650c8SDave May ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 561b16650c8SDave May ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 562b16650c8SDave May 563b16650c8SDave May for (p=0; p<npoints; p++) { 564b16650c8SDave May if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) { 565b16650c8SDave May if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) { 566b16650c8SDave May // copy point into buffer // 567b16650c8SDave May ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr); 568b16650c8SDave May // insert point buffer into DataExchanger // 569b16650c8SDave May ierr = DataExPackData(de,recv_bbox[pk].owner_rank,1,point_buffer);CHKERRQ(ierr); 570b16650c8SDave May } 571b16650c8SDave May } 572b16650c8SDave May } 573b16650c8SDave May 574b16650c8SDave May ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 575b16650c8SDave May ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 576b16650c8SDave May } 577b16650c8SDave May 578b16650c8SDave May ierr = DataExPackFinalize(de);CHKERRQ(ierr); 579b16650c8SDave May 58008056efcSDave May ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 581b16650c8SDave May 582b16650c8SDave May ierr = DataExBegin(de);CHKERRQ(ierr); 583b16650c8SDave May ierr = DataExEnd(de);CHKERRQ(ierr); 584b16650c8SDave May 585b16650c8SDave May ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr); 586b16650c8SDave May 587b16650c8SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 588b16650c8SDave May ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr); 589b16650c8SDave May for (p=0; p<n_points_recv; p++) { 590b16650c8SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point ); 591b16650c8SDave May 592b16650c8SDave May ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr); 593b16650c8SDave May } 594b16650c8SDave May 595b16650c8SDave May ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr); 596b16650c8SDave May PetscFree(bbox); 597b16650c8SDave May ierr = DataExView(de);CHKERRQ(ierr); 598b16650c8SDave May ierr = DataExDestroy(de);CHKERRQ(ierr); 599b16650c8SDave May 600b16650c8SDave May PetscFunctionReturn(0); 601b16650c8SDave May } 602a9fd7477SDave May 603a9fd7477SDave May 604a9fd7477SDave May /* General collection when no order, or neighbour information is provided */ 605a9fd7477SDave May /* 606a9fd7477SDave May User provides context and collect() method 607a9fd7477SDave May Broadcast user context 608a9fd7477SDave May 609a9fd7477SDave May for each context / rank { 610a9fd7477SDave May collect(swarm,context,n,list) 611a9fd7477SDave May } 612a9fd7477SDave May 613a9fd7477SDave May */ 614a9fd7477SDave May #undef __FUNCT__ 615a9fd7477SDave May #define __FUNCT__ "DMSwarmCollect_General" 616a9fd7477SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM dm,PetscErrorCode (*collect)(DM,void*,PetscInt*,PetscInt**),size_t ctx_size,void *ctx,PetscInt *globalsize) 617a9fd7477SDave May { 618a9fd7477SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 619a9fd7477SDave May PetscErrorCode ierr; 620a9fd7477SDave May DataEx de; 621a9fd7477SDave May PetscInt p,r,npoints,n_points_recv; 622a9fd7477SDave May PetscMPIInt commsize,rank; 623a9fd7477SDave May void *point_buffer,*recv_points; 624a9fd7477SDave May void *ctxlist; 625a9fd7477SDave May PetscInt *n2collect,**collectlist; 626a9fd7477SDave May size_t sizeof_dmswarm_point; 627a9fd7477SDave May 628a9fd7477SDave May ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr); 629a9fd7477SDave May ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 630a9fd7477SDave May 631a9fd7477SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 632a9fd7477SDave May *globalsize = npoints; 633a9fd7477SDave May 634a9fd7477SDave May /* Broadcast user context */ 635a9fd7477SDave May PetscMalloc(ctx_size*commsize,&ctxlist); 636a9fd7477SDave May ierr = MPI_Allgather(ctx,ctx_size,MPI_CHAR,ctxlist,ctx_size,MPI_CHAR,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 637a9fd7477SDave May 638a9fd7477SDave May PetscMalloc1(commsize,&n2collect); 639a9fd7477SDave May PetscMalloc1(commsize,&collectlist); 640a9fd7477SDave May 641a9fd7477SDave May for (r=0; r<commsize; r++) { 642a9fd7477SDave May PetscInt _n2collect; 643a9fd7477SDave May PetscInt *_collectlist; 644a9fd7477SDave May void *_ctx_r; 645a9fd7477SDave May 646a9fd7477SDave May _n2collect = 0; 647a9fd7477SDave May _collectlist = NULL; 648a9fd7477SDave May 649a9fd7477SDave May if (r != rank) { /* don't collect data from yourself */ 650a9fd7477SDave May _ctx_r = (void*)( (char*)ctxlist + r * ctx_size ); 651a9fd7477SDave May ierr = collect(dm,_ctx_r,&_n2collect,&_collectlist);CHKERRQ(ierr); 652a9fd7477SDave May } 653a9fd7477SDave May 654a9fd7477SDave May n2collect[r] = _n2collect; 655a9fd7477SDave May collectlist[r] = _collectlist; 656a9fd7477SDave May } 657a9fd7477SDave May 658a9fd7477SDave May de = DataExCreate(PetscObjectComm((PetscObject)dm),0);CHKERRQ(ierr); 659a9fd7477SDave May 660a9fd7477SDave May /* Define topology */ 661a9fd7477SDave May ierr = DataExTopologyInitialize(de);CHKERRQ(ierr); 662a9fd7477SDave May for (r=0; r<commsize; r++) { 663a9fd7477SDave May if (n2collect[r] > 0) { 664a9fd7477SDave May ierr = DataExTopologyAddNeighbour(de,(PetscMPIInt)r);CHKERRQ(ierr); 665a9fd7477SDave May } 666a9fd7477SDave May } 667a9fd7477SDave May ierr = DataExTopologyFinalize(de);CHKERRQ(ierr); 668a9fd7477SDave May 669a9fd7477SDave May /* Define send counts */ 670a9fd7477SDave May ierr = DataExInitializeSendCount(de);CHKERRQ(ierr); 671a9fd7477SDave May for (r=0; r<commsize; r++) { 672a9fd7477SDave May if (n2collect[r] > 0) { 673a9fd7477SDave May ierr = DataExAddToSendCount(de,r,n2collect[r]);CHKERRQ(ierr); 674a9fd7477SDave May } 675a9fd7477SDave May } 676a9fd7477SDave May ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr); 677a9fd7477SDave May 678a9fd7477SDave May /* Pack data */ 679a9fd7477SDave May ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr); 680a9fd7477SDave May 681a9fd7477SDave May ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr); 682a9fd7477SDave May for (r=0; r<commsize; r++) { 683a9fd7477SDave May 684a9fd7477SDave May for (p=0; p<n2collect[r]; p++) { 685a9fd7477SDave May ierr = DataBucketFillPackedArray(swarm->db,collectlist[r][p],point_buffer);CHKERRQ(ierr); 686a9fd7477SDave May /* insert point buffer into the data exchanger */ 687a9fd7477SDave May ierr = DataExPackData(de,r,1,point_buffer);CHKERRQ(ierr); 688a9fd7477SDave May } 689a9fd7477SDave May } 690a9fd7477SDave May 691a9fd7477SDave May ierr = DataExPackFinalize(de);CHKERRQ(ierr); 692a9fd7477SDave May 693a9fd7477SDave May /* Scatter */ 694a9fd7477SDave May ierr = DataExBegin(de);CHKERRQ(ierr); 695a9fd7477SDave May ierr = DataExEnd(de);CHKERRQ(ierr); 696a9fd7477SDave May 697a9fd7477SDave May /* Collect data in DMSwarm container */ 698a9fd7477SDave May ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr); 699a9fd7477SDave May 700a9fd7477SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 701a9fd7477SDave May ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr); 702a9fd7477SDave May for (p=0; p<n_points_recv; p++) { 703a9fd7477SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point ); 704a9fd7477SDave May 705a9fd7477SDave May ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr); 706a9fd7477SDave May } 707a9fd7477SDave May 708a9fd7477SDave May /* Release memory */ 709a9fd7477SDave May for (r=0; r<commsize; r++) { 710a9fd7477SDave May if (collectlist[r]) PetscFree(collectlist[r]); 711a9fd7477SDave May } 712a9fd7477SDave May PetscFree(collectlist); 713a9fd7477SDave May PetscFree(n2collect); 714a9fd7477SDave May PetscFree(ctxlist); 715a9fd7477SDave May ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr); 716a9fd7477SDave May ierr = DataExView(de);CHKERRQ(ierr); 717a9fd7477SDave May ierr = DataExDestroy(de);CHKERRQ(ierr); 718a9fd7477SDave May 719a9fd7477SDave May PetscFunctionReturn(0); 720a9fd7477SDave May } 721a9fd7477SDave May 722