1df21e3a8SDave May 2dfc14de9SMatthew G. Knepley #include <petscsf.h> 308056efcSDave May #include <petscdmswarm.h> 4df21e3a8SDave May #include <petsc/private/dmswarmimpl.h> /*I "petscdmswarm.h" I*/ 5df21e3a8SDave May #include "data_bucket.h" 6df21e3a8SDave May #include "data_ex.h" 7df21e3a8SDave May 8df21e3a8SDave May 9480eef7bSDave May /* 10480eef7bSDave May User loads desired location (MPI rank) into field DMSwarm_rank 11480eef7bSDave May */ 12df21e3a8SDave May #undef __FUNCT__ 13df21e3a8SDave May #define __FUNCT__ "DMSwarmMigrate_Push_Basic" 14df21e3a8SDave May PetscErrorCode DMSwarmMigrate_Push_Basic(DM dm,PetscBool remove_sent_points) 15df21e3a8SDave May { 16df21e3a8SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 17df21e3a8SDave May PetscErrorCode ierr; 18df21e3a8SDave May DataEx de; 19df21e3a8SDave May PetscInt p,npoints,*rankval,n_points_recv; 20df21e3a8SDave May PetscMPIInt rank,nrank; 21df21e3a8SDave May void *point_buffer,*recv_points; 22df21e3a8SDave May size_t sizeof_dmswarm_point; 23df21e3a8SDave May 24521f74f9SMatthew G. Knepley PetscFunctionBegin; 25df21e3a8SDave May ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 26df21e3a8SDave May 27df21e3a8SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 2808056efcSDave May ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 29521f74f9SMatthew G. Knepley ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0, &de);CHKERRQ(ierr); 30df21e3a8SDave May ierr = DataExTopologyInitialize(de);CHKERRQ(ierr); 31521f74f9SMatthew G. Knepley 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 ierr = DataExInitializeSendCount(de);CHKERRQ(ierr); 39df21e3a8SDave May for (p=0; p<npoints; p++) { 40df21e3a8SDave May nrank = rankval[p]; 41df21e3a8SDave May if (nrank != rank) { 42df21e3a8SDave May ierr = DataExAddToSendCount(de,nrank,1);CHKERRQ(ierr); 43df21e3a8SDave May } 44df21e3a8SDave May } 45df21e3a8SDave May ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr); 46df21e3a8SDave May ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr); 47df21e3a8SDave May ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr); 48df21e3a8SDave May for (p=0; p<npoints; p++) { 49df21e3a8SDave May nrank = rankval[p]; 50df21e3a8SDave May if (nrank != rank) { 51df21e3a8SDave May /* copy point into buffer */ 52df21e3a8SDave May ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr); 53df21e3a8SDave May /* insert point buffer into DataExchanger */ 54df21e3a8SDave May ierr = DataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr); 55df21e3a8SDave May } 56df21e3a8SDave May } 57df21e3a8SDave May ierr = DataExPackFinalize(de);CHKERRQ(ierr); 58*22a417f9SDave May ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 59df21e3a8SDave May 60df21e3a8SDave May if (remove_sent_points) { 61*22a417f9SDave May DataField gfield; 62*22a417f9SDave May 63*22a417f9SDave May ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&gfield);CHKERRQ(ierr); 64*22a417f9SDave May ierr = DataFieldGetAccess(gfield);CHKERRQ(ierr); 65*22a417f9SDave May ierr = DataFieldGetEntries(gfield,(void**)&rankval);CHKERRQ(ierr); 66*22a417f9SDave May 67df21e3a8SDave May /* remove points which left processor */ 68df21e3a8SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 69df21e3a8SDave May for (p=0; p<npoints; p++) { 70df21e3a8SDave May nrank = rankval[p]; 71df21e3a8SDave May if (nrank != rank) { 72df21e3a8SDave May /* kill point */ 73*22a417f9SDave May ierr = DataFieldRestoreAccess(gfield);CHKERRQ(ierr); 74*22a417f9SDave May 75df21e3a8SDave May ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr); 76*22a417f9SDave May 77*22a417f9SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */ 78*22a417f9SDave May ierr = DataFieldGetAccess(gfield);CHKERRQ(ierr); 79*22a417f9SDave May ierr = DataFieldGetEntries(gfield,(void**)&rankval);CHKERRQ(ierr); 80df21e3a8SDave May p--; /* check replacement point */ 81df21e3a8SDave May } 82df21e3a8SDave May } 83*22a417f9SDave May ierr = DataFieldRestoreEntries(gfield,(void**)&rankval);CHKERRQ(ierr); 84*22a417f9SDave May ierr = DataFieldRestoreAccess(gfield);CHKERRQ(ierr); 85df21e3a8SDave May } 86df21e3a8SDave May ierr = DataExBegin(de);CHKERRQ(ierr); 87df21e3a8SDave May ierr = DataExEnd(de);CHKERRQ(ierr); 88df21e3a8SDave May ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr); 89df21e3a8SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 90df21e3a8SDave May ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr); 91df21e3a8SDave May for (p=0; p<n_points_recv; p++) { 92df21e3a8SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point ); 93df21e3a8SDave May 94df21e3a8SDave May ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr); 95df21e3a8SDave May } 96df21e3a8SDave May ierr = DataExView(de);CHKERRQ(ierr); 97df21e3a8SDave May ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr); 98df21e3a8SDave May ierr = DataExDestroy(de);CHKERRQ(ierr); 99df21e3a8SDave May PetscFunctionReturn(0); 100df21e3a8SDave May } 1012712d1f2SDave May 10240c453e9SDave May #undef __FUNCT__ 103889dbfe5SDave May #define __FUNCT__ "DMSwarmMigrate_DMNeighborScatter" 104889dbfe5SDave May PetscErrorCode DMSwarmMigrate_DMNeighborScatter(DM dm,DM dmcell,PetscBool remove_sent_points,PetscInt *npoints_prior_migration) 10540c453e9SDave May { 10640c453e9SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 10740c453e9SDave May PetscErrorCode ierr; 10840c453e9SDave May DataEx de; 10940c453e9SDave May PetscInt r,p,npoints,*rankval,n_points_recv; 11040c453e9SDave May PetscMPIInt rank,_rank; 11140c453e9SDave May const PetscMPIInt *neighbourranks; 11240c453e9SDave May void *point_buffer,*recv_points; 11340c453e9SDave May size_t sizeof_dmswarm_point; 11440c453e9SDave May PetscInt nneighbors; 1157c6d1d28SDave May PetscMPIInt mynneigh,*myneigh; 11640c453e9SDave May 117521f74f9SMatthew G. Knepley PetscFunctionBegin; 11840c453e9SDave May ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 11940c453e9SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 12008056efcSDave May ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 121521f74f9SMatthew G. Knepley ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr); 12240c453e9SDave May ierr = DMGetNeighbors(dmcell,&nneighbors,&neighbourranks);CHKERRQ(ierr); 12340c453e9SDave May ierr = DataExTopologyInitialize(de);CHKERRQ(ierr); 12440c453e9SDave May for (r=0; r<nneighbors; r++) { 12540c453e9SDave May _rank = neighbourranks[r]; 12640c453e9SDave May if ((_rank != rank) && (_rank > 0)) { 12740c453e9SDave May ierr = DataExTopologyAddNeighbour(de,_rank);CHKERRQ(ierr); 12840c453e9SDave May } 12940c453e9SDave May } 13040c453e9SDave May ierr = DataExTopologyFinalize(de);CHKERRQ(ierr); 1317c6d1d28SDave May ierr = DataExTopologyGetNeighbours(de,&mynneigh,&myneigh);CHKERRQ(ierr); 13240c453e9SDave May ierr = DataExInitializeSendCount(de);CHKERRQ(ierr); 13340c453e9SDave May for (p=0; p<npoints; p++) { 134f954cb40SDave May if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) { 1357c6d1d28SDave May for (r=0; r<mynneigh; r++) { 1367c6d1d28SDave May _rank = myneigh[r]; 13740c453e9SDave May ierr = DataExAddToSendCount(de,_rank,1);CHKERRQ(ierr); 13840c453e9SDave May } 13940c453e9SDave May } 14040c453e9SDave May } 14140c453e9SDave May ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr); 14240c453e9SDave May ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr); 14340c453e9SDave May ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr); 14440c453e9SDave May for (p=0; p<npoints; p++) { 145f954cb40SDave May if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) { 1467c6d1d28SDave May for (r=0; r<mynneigh; r++) { 1477c6d1d28SDave May _rank = myneigh[r]; 14840c453e9SDave May /* copy point into buffer */ 14940c453e9SDave May ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr); 15040c453e9SDave May /* insert point buffer into DataExchanger */ 15140c453e9SDave May ierr = DataExPackData(de,_rank,1,point_buffer);CHKERRQ(ierr); 15240c453e9SDave May } 15340c453e9SDave May } 15440c453e9SDave May } 15540c453e9SDave May ierr = DataExPackFinalize(de);CHKERRQ(ierr); 1567c6d1d28SDave May ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 15740c453e9SDave May if (remove_sent_points) { 1587c6d1d28SDave May DataField PField; 1597c6d1d28SDave May 1607c6d1d28SDave May ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr); 1617c6d1d28SDave May ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); 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++) { 165f954cb40SDave May if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) { 16640c453e9SDave May /* kill point */ 16740c453e9SDave May ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr); 1687c6d1d28SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */ 1697c6d1d28SDave 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 ierr = DataExBegin(de);CHKERRQ(ierr); 17640c453e9SDave May ierr = DataExEnd(de);CHKERRQ(ierr); 17740c453e9SDave May ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr); 17840c453e9SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 17940c453e9SDave May ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr); 18040c453e9SDave May for (p=0; p<n_points_recv; p++) { 18140c453e9SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point ); 18240c453e9SDave May 18340c453e9SDave May ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr); 18440c453e9SDave May } 18540c453e9SDave May ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr); 18640c453e9SDave May ierr = DataExDestroy(de);CHKERRQ(ierr); 18740c453e9SDave May PetscFunctionReturn(0); 18840c453e9SDave May } 189480eef7bSDave May 190480eef7bSDave May #undef __FUNCT__ 19108056efcSDave May #define __FUNCT__ "DMSwarmMigrate_CellDMScatter" 19208056efcSDave May PetscErrorCode DMSwarmMigrate_CellDMScatter(DM dm,PetscBool remove_sent_points) 193480eef7bSDave May { 194480eef7bSDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 195480eef7bSDave May PetscErrorCode ierr; 1966fbf25f8SDave May PetscInt p,npoints,npointsg=0,npoints2,npoints2g,*rankval,npoints_prior_migration; 197bbe8250bSMatthew G. Knepley PetscSF sfcell = NULL; 198dfc14de9SMatthew G. Knepley const PetscSFNode *LA_sfcell; 199480eef7bSDave May DM dmcell; 200480eef7bSDave May Vec pos; 20140c453e9SDave May PetscBool error_check = swarm->migrate_error_on_missing_point; 2027c6d1d28SDave May PetscMPIInt commsize,rank; 203480eef7bSDave May 204521f74f9SMatthew G. Knepley PetscFunctionBegin; 205480eef7bSDave May ierr = DMSwarmGetCellDM(dm,&dmcell);CHKERRQ(ierr); 206480eef7bSDave May if (!dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided"); 207480eef7bSDave May 2087c6d1d28SDave May ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr); 2097c6d1d28SDave May ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 2107c6d1d28SDave May 211fb1bcc12SMatthew G. Knepley ierr = DMSwarmCreateLocalVectorFromField(dm, DMSwarmPICField_coor, &pos);CHKERRQ(ierr); 212dfc14de9SMatthew G. Knepley ierr = DMLocatePoints(dmcell, pos, DM_POINTLOCATION_NONE, &sfcell);CHKERRQ(ierr); 213fb1bcc12SMatthew G. Knepley ierr = DMSwarmDestroyLocalVectorFromField(dm, DMSwarmPICField_coor, &pos);CHKERRQ(ierr); 214480eef7bSDave May 21540c453e9SDave May if (error_check) { 21640c453e9SDave May ierr = DMSwarmGetSize(dm,&npointsg);CHKERRQ(ierr); 21740c453e9SDave May } 218480eef7bSDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 21908056efcSDave May ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 220dfc14de9SMatthew G. Knepley ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr); 221480eef7bSDave May for (p=0; p<npoints; p++) { 222dfc14de9SMatthew G. Knepley rankval[p] = LA_sfcell[p].index; 223480eef7bSDave May } 22408056efcSDave May ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 225dfc14de9SMatthew G. Knepley ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr); 226480eef7bSDave May 2276fbf25f8SDave May if (commsize > 1) { 228889dbfe5SDave May ierr = DMSwarmMigrate_DMNeighborScatter(dm,dmcell,remove_sent_points,&npoints_prior_migration);CHKERRQ(ierr); 2296fbf25f8SDave May } else { 2306fbf25f8SDave May ierr = DMSwarmGetSize(dm,&npoints_prior_migration);CHKERRQ(ierr); 2316fbf25f8SDave May } 232480eef7bSDave May 23340c453e9SDave May /* locate points newly recevied */ 23440c453e9SDave May ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); 2357c6d1d28SDave May #if 0 23640c453e9SDave May len = npoints2 - npoints_prior_migration; 23740c453e9SDave May if (len > 0) { 23840c453e9SDave May PetscScalar *LA_coor; 2397c6d1d28SDave May PetscInt bs; 24040c453e9SDave May 2417c6d1d28SDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr); 2427c6d1d28SDave May ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*len,(const PetscScalar*)&LA_coor[bs*npoints_prior_migration],&pos);CHKERRQ(ierr); 24340c453e9SDave May 24440c453e9SDave May ierr = DMLocatePoints(dmcell,pos,&iscell);CHKERRQ(ierr); 24540c453e9SDave May 24640c453e9SDave May ierr = VecDestroy(&pos);CHKERRQ(ierr); 2477c6d1d28SDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr); 24840c453e9SDave May 24940c453e9SDave May ierr = ISGetIndices(iscell,&LA_iscell);CHKERRQ(ierr); 25008056efcSDave May ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 25140c453e9SDave May for (p=0; p<len; p++) { 252f954cb40SDave May rankval[npoints_prior_migration+p] = LA_iscell[p]; 25340c453e9SDave May } 25440c453e9SDave May ierr = ISRestoreIndices(iscell,&LA_iscell);CHKERRQ(ierr); 25540c453e9SDave May ierr = ISDestroy(&iscell);CHKERRQ(ierr); 2567c6d1d28SDave May ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 25740c453e9SDave May 25840c453e9SDave May /* remove points which left processor */ 2597c6d1d28SDave May ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr); 2607c6d1d28SDave May ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); 2617c6d1d28SDave May 26240c453e9SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 26340c453e9SDave May for (p=npoints_prior_migration; p<npoints2; p++) { 264f954cb40SDave May if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) { 26540c453e9SDave May /* kill point */ 26640c453e9SDave May ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr); 2677c6d1d28SDave May ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */ 2687c6d1d28SDave May ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */ 26940c453e9SDave May p--; /* check replacement point */ 27040c453e9SDave May } 27140c453e9SDave May } 2727c6d1d28SDave May 2737c6d1d28SDave May } 2747c6d1d28SDave May #endif 2757c6d1d28SDave May 2767c6d1d28SDave May { 2777c6d1d28SDave May PetscScalar *LA_coor; 2787c6d1d28SDave May PetscInt bs; 2797c6d1d28SDave May DataField PField; 2807c6d1d28SDave May 2817c6d1d28SDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr); 2827c6d1d28SDave May ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*npoints2,(const PetscScalar*)LA_coor,&pos);CHKERRQ(ierr); 283dfc14de9SMatthew G. Knepley ierr = DMLocatePoints(dmcell,pos,DM_POINTLOCATION_NONE,&sfcell);CHKERRQ(ierr); 2847c6d1d28SDave May 2857c6d1d28SDave May ierr = VecDestroy(&pos);CHKERRQ(ierr); 2867c6d1d28SDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr); 2877c6d1d28SDave May 288dfc14de9SMatthew G. Knepley ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr); 2897c6d1d28SDave May ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 2907c6d1d28SDave May for (p=0; p<npoints2; p++) { 291dfc14de9SMatthew G. Knepley rankval[p] = LA_sfcell[p].index; 2927c6d1d28SDave May } 293dfc14de9SMatthew G. Knepley ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr); 29408056efcSDave May ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 29540c453e9SDave May 2967c6d1d28SDave May /* remove points which left processor */ 2977c6d1d28SDave May ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr); 2987c6d1d28SDave May ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); 2997c6d1d28SDave May 3007c6d1d28SDave May ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); 3017c6d1d28SDave May for (p=0; p<npoints2; p++) { 302f954cb40SDave May if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) { 3037c6d1d28SDave May /* kill point */ 3047c6d1d28SDave May ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr); 3057c6d1d28SDave May ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */ 3067c6d1d28SDave May ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */ 3077c6d1d28SDave May p--; /* check replacement point */ 3087c6d1d28SDave May } 3097c6d1d28SDave May } 31040c453e9SDave May } 31140c453e9SDave May /* check for error on removed points */ 31240c453e9SDave May if (error_check) { 31340c453e9SDave May ierr = DMSwarmGetSize(dm,&npoints2g);CHKERRQ(ierr); 31440c453e9SDave 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); 31540c453e9SDave May } 316480eef7bSDave May PetscFunctionReturn(0); 317480eef7bSDave May } 318480eef7bSDave May 31908056efcSDave May #undef __FUNCT__ 32008056efcSDave May #define __FUNCT__ "DMSwarmMigrate_CellDMExact" 32108056efcSDave May PetscErrorCode DMSwarmMigrate_CellDMExact(DM dm,PetscBool remove_sent_points) 32208056efcSDave May { 323521f74f9SMatthew G. Knepley PetscFunctionBegin; 32408056efcSDave May PetscFunctionReturn(0); 32508056efcSDave May } 32608056efcSDave May 327480eef7bSDave May /* 328480eef7bSDave May Redundant as this assumes points can only be sent to a single rank 329480eef7bSDave May */ 3302712d1f2SDave May #undef __FUNCT__ 3312712d1f2SDave May #define __FUNCT__ "DMSwarmMigrate_GlobalToLocal_Basic" 3322712d1f2SDave May PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt *globalsize) 3332712d1f2SDave May { 3342712d1f2SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 3352712d1f2SDave May PetscErrorCode ierr; 3362712d1f2SDave May DataEx de; 3372712d1f2SDave May PetscInt p,npoints,*rankval,n_points_recv; 3382712d1f2SDave May PetscMPIInt rank,nrank,negrank; 3392712d1f2SDave May void *point_buffer,*recv_points; 3402712d1f2SDave May size_t sizeof_dmswarm_point; 3412712d1f2SDave May 342521f74f9SMatthew G. Knepley PetscFunctionBegin; 3432712d1f2SDave May ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 3442712d1f2SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 3452712d1f2SDave May *globalsize = npoints; 34608056efcSDave May ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 347521f74f9SMatthew G. Knepley ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr); 3482712d1f2SDave May ierr = DataExTopologyInitialize(de);CHKERRQ(ierr); 3492712d1f2SDave May for (p=0; p<npoints; p++) { 3502712d1f2SDave May negrank = rankval[p]; 3512712d1f2SDave May if (negrank < 0) { 3522712d1f2SDave May nrank = -negrank - 1; 3532712d1f2SDave May ierr = DataExTopologyAddNeighbour(de,nrank);CHKERRQ(ierr); 3542712d1f2SDave May } 3552712d1f2SDave May } 3562712d1f2SDave May ierr = DataExTopologyFinalize(de);CHKERRQ(ierr); 3572712d1f2SDave May ierr = DataExInitializeSendCount(de);CHKERRQ(ierr); 3582712d1f2SDave May for (p=0; p<npoints; p++) { 3592712d1f2SDave May negrank = rankval[p]; 3602712d1f2SDave May if (negrank < 0) { 3612712d1f2SDave May nrank = -negrank - 1; 3622712d1f2SDave May ierr = DataExAddToSendCount(de,nrank,1);CHKERRQ(ierr); 3632712d1f2SDave May } 3642712d1f2SDave May } 3652712d1f2SDave May ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr); 3662712d1f2SDave May ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr); 3672712d1f2SDave May ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr); 3682712d1f2SDave May for (p=0; p<npoints; p++) { 3692712d1f2SDave May negrank = rankval[p]; 3702712d1f2SDave May if (negrank < 0) { 3712712d1f2SDave May nrank = -negrank - 1; 3722712d1f2SDave May rankval[p] = nrank; 3732712d1f2SDave May /* copy point into buffer */ 3742712d1f2SDave May ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr); 3752712d1f2SDave May /* insert point buffer into DataExchanger */ 3762712d1f2SDave May ierr = DataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr); 3772712d1f2SDave May rankval[p] = negrank; 3782712d1f2SDave May } 3792712d1f2SDave May } 3802712d1f2SDave May ierr = DataExPackFinalize(de);CHKERRQ(ierr); 38108056efcSDave May ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 3822712d1f2SDave May ierr = DataExBegin(de);CHKERRQ(ierr); 3832712d1f2SDave May ierr = DataExEnd(de);CHKERRQ(ierr); 3842712d1f2SDave May ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr); 3852712d1f2SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 3862712d1f2SDave May ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr); 3872712d1f2SDave May for (p=0; p<n_points_recv; p++) { 3882712d1f2SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point ); 3892712d1f2SDave May 3902712d1f2SDave May ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr); 3912712d1f2SDave May } 3922712d1f2SDave May ierr = DataExView(de);CHKERRQ(ierr); 3932712d1f2SDave May ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr); 3942712d1f2SDave May ierr = DataExDestroy(de);CHKERRQ(ierr); 3952712d1f2SDave May PetscFunctionReturn(0); 3962712d1f2SDave May } 397b16650c8SDave May 398b16650c8SDave May typedef struct { 399b16650c8SDave May PetscMPIInt owner_rank; 400b16650c8SDave May PetscReal min[3],max[3]; 401b16650c8SDave May } CollectBBox; 402b16650c8SDave May 403b16650c8SDave May #undef __FUNCT__ 404fe39f135SDave May #define __FUNCT__ "DMSwarmCollect_DMDABoundingBox" 405fe39f135SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM dm,PetscInt *globalsize) 406b16650c8SDave May { 407b16650c8SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 408b16650c8SDave May PetscErrorCode ierr; 409b16650c8SDave May DataEx de; 410b16650c8SDave May PetscInt p,pk,npoints,*rankval,n_points_recv,n_bbox_recv,dim,neighbour_cells; 411b16650c8SDave May PetscMPIInt rank,nrank; 412b16650c8SDave May void *point_buffer,*recv_points; 413b16650c8SDave May size_t sizeof_dmswarm_point,sizeof_bbox_ctx; 414b16650c8SDave May PetscBool isdmda; 415b16650c8SDave May CollectBBox *bbox,*recv_bbox; 416b16650c8SDave May const PetscMPIInt *dmneighborranks; 417b16650c8SDave May DM dmcell; 418b16650c8SDave May 419521f74f9SMatthew G. Knepley PetscFunctionBegin; 420b16650c8SDave May ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 421b16650c8SDave May 422fe39f135SDave May ierr = DMSwarmGetCellDM(dm,&dmcell);CHKERRQ(ierr); 423fe39f135SDave May if (!dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided"); 424b16650c8SDave May isdmda = PETSC_FALSE; 425b16650c8SDave May PetscObjectTypeCompare((PetscObject)dmcell,DMDA,&isdmda); 426b16650c8SDave May if (!isdmda) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only DMDA support for CollectBoundingBox"); 427b16650c8SDave May 4288dbd68bcSDave May ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 429b16650c8SDave May sizeof_bbox_ctx = sizeof(CollectBBox); 430b16650c8SDave May PetscMalloc1(1,&bbox); 431b16650c8SDave May bbox->owner_rank = rank; 432b16650c8SDave May 433b16650c8SDave May /* compute the bounding box based on the overlapping / stenctil size */ 434b16650c8SDave May { 435b16650c8SDave May Vec lcoor; 436b16650c8SDave May 437b16650c8SDave May ierr = DMGetCoordinatesLocal(dmcell,&lcoor);CHKERRQ(ierr); 438fe39f135SDave May if (dim >= 1) { 439b16650c8SDave May ierr = VecStrideMin(lcoor,0,NULL,&bbox->min[0]);CHKERRQ(ierr); 440b16650c8SDave May ierr = VecStrideMax(lcoor,0,NULL,&bbox->max[0]);CHKERRQ(ierr); 441fe39f135SDave May } 442fe39f135SDave May if (dim >= 2) { 443fe39f135SDave May ierr = VecStrideMin(lcoor,1,NULL,&bbox->min[1]);CHKERRQ(ierr); 444b16650c8SDave May ierr = VecStrideMax(lcoor,1,NULL,&bbox->max[1]);CHKERRQ(ierr); 445b16650c8SDave May } 446fe39f135SDave May if (dim == 3) { 447fe39f135SDave May ierr = VecStrideMin(lcoor,2,NULL,&bbox->min[2]);CHKERRQ(ierr); 448fe39f135SDave May ierr = VecStrideMax(lcoor,2,NULL,&bbox->max[2]);CHKERRQ(ierr); 449fe39f135SDave May } 450fe39f135SDave May } 451b16650c8SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 452b16650c8SDave May *globalsize = npoints; 45308056efcSDave May ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 454521f74f9SMatthew G. Knepley ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr); 455b16650c8SDave May /* use DMDA neighbours */ 456b16650c8SDave May ierr = DMDAGetNeighbors(dmcell,&dmneighborranks);CHKERRQ(ierr); 4578dbd68bcSDave May if (dim == 1) { 4588dbd68bcSDave May neighbour_cells = 3; 4598dbd68bcSDave May } else if (dim == 2) { 4608dbd68bcSDave May neighbour_cells = 9; 4618dbd68bcSDave May } else { 4628dbd68bcSDave May neighbour_cells = 27; 4638dbd68bcSDave May } 464b16650c8SDave May ierr = DataExTopologyInitialize(de);CHKERRQ(ierr); 465b16650c8SDave May for (p=0; p<neighbour_cells; p++) { 466b16650c8SDave May if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) { 467b16650c8SDave May ierr = DataExTopologyAddNeighbour(de,dmneighborranks[p]);CHKERRQ(ierr); 468b16650c8SDave May } 469b16650c8SDave May } 470b16650c8SDave May ierr = DataExTopologyFinalize(de);CHKERRQ(ierr); 471b16650c8SDave May ierr = DataExInitializeSendCount(de);CHKERRQ(ierr); 472b16650c8SDave May for (p=0; p<neighbour_cells; p++) { 473b16650c8SDave May if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) { 474b16650c8SDave May ierr = DataExAddToSendCount(de,dmneighborranks[p],1);CHKERRQ(ierr); 475b16650c8SDave May } 476b16650c8SDave May } 477b16650c8SDave May ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr); 478b16650c8SDave May /* send bounding boxes */ 479b16650c8SDave May ierr = DataExPackInitialize(de,sizeof_bbox_ctx);CHKERRQ(ierr); 480b16650c8SDave May for (p=0; p<neighbour_cells; p++) { 481b16650c8SDave May nrank = dmneighborranks[p]; 482b16650c8SDave May if ( (nrank >= 0) && (nrank != rank) ) { 483b16650c8SDave May /* insert bbox buffer into DataExchanger */ 484b16650c8SDave May ierr = DataExPackData(de,nrank,1,bbox);CHKERRQ(ierr); 485b16650c8SDave May } 486b16650c8SDave May } 487b16650c8SDave May ierr = DataExPackFinalize(de);CHKERRQ(ierr); 488b16650c8SDave May /* recv bounding boxes */ 489b16650c8SDave May ierr = DataExBegin(de);CHKERRQ(ierr); 490b16650c8SDave May ierr = DataExEnd(de);CHKERRQ(ierr); 491b16650c8SDave May ierr = DataExGetRecvData(de,&n_bbox_recv,(void**)&recv_bbox);CHKERRQ(ierr); 492b16650c8SDave May for (p=0; p<n_bbox_recv; p++) { 493b4b02483SDave May PetscPrintf(PETSC_COMM_SELF,"[rank %d]: box from %d : range[%+1.4e,%+1.4e]x[%+1.4e,%+1.4e]\n",rank,recv_bbox[p].owner_rank, 494b4b02483SDave May (double)recv_bbox[p].min[0],(double)recv_bbox[p].max[0],(double)recv_bbox[p].min[1],(double)recv_bbox[p].max[1]); 495b16650c8SDave May } 496b16650c8SDave May /* of course this is stupid as this "generic" function should have a better way to know what the coordinates are called */ 497b16650c8SDave May ierr = DataExInitializeSendCount(de);CHKERRQ(ierr); 498b16650c8SDave May for (pk=0; pk<n_bbox_recv; pk++) { 499b16650c8SDave May PetscReal *array_x,*array_y; 500b16650c8SDave May 501b16650c8SDave May ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 502b16650c8SDave May ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 503b16650c8SDave May for (p=0; p<npoints; p++) { 504b16650c8SDave May if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) { 505b16650c8SDave May if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) { 506b16650c8SDave May ierr = DataExAddToSendCount(de,recv_bbox[pk].owner_rank,1);CHKERRQ(ierr); 507b16650c8SDave May } 508b16650c8SDave May } 509b16650c8SDave May } 510b16650c8SDave May ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 511b16650c8SDave May ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 512b16650c8SDave May } 513b16650c8SDave May ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr); 514b16650c8SDave May ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr); 515b16650c8SDave May ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr); 516b16650c8SDave May for (pk=0; pk<n_bbox_recv; pk++) { 517b16650c8SDave May PetscReal *array_x,*array_y; 518b16650c8SDave May 519b16650c8SDave May ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 520b16650c8SDave May ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 521b16650c8SDave May for (p=0; p<npoints; p++) { 522b16650c8SDave May if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) { 523b16650c8SDave May if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) { 524521f74f9SMatthew G. Knepley /* copy point into buffer */ 525b16650c8SDave May ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr); 526521f74f9SMatthew G. Knepley /* insert point buffer into DataExchanger */ 527b16650c8SDave May ierr = DataExPackData(de,recv_bbox[pk].owner_rank,1,point_buffer);CHKERRQ(ierr); 528b16650c8SDave May } 529b16650c8SDave May } 530b16650c8SDave May } 531b16650c8SDave May ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 532b16650c8SDave May ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 533b16650c8SDave May } 534b16650c8SDave May ierr = DataExPackFinalize(de);CHKERRQ(ierr); 53508056efcSDave May ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 536b16650c8SDave May ierr = DataExBegin(de);CHKERRQ(ierr); 537b16650c8SDave May ierr = DataExEnd(de);CHKERRQ(ierr); 538b16650c8SDave May ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr); 539b16650c8SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 540b16650c8SDave May ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr); 541b16650c8SDave May for (p=0; p<n_points_recv; p++) { 542b16650c8SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point ); 543b16650c8SDave May 544b16650c8SDave May ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr); 545b16650c8SDave May } 546b16650c8SDave May ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr); 547b16650c8SDave May PetscFree(bbox); 548b16650c8SDave May ierr = DataExView(de);CHKERRQ(ierr); 549b16650c8SDave May ierr = DataExDestroy(de);CHKERRQ(ierr); 550b16650c8SDave May PetscFunctionReturn(0); 551b16650c8SDave May } 552a9fd7477SDave May 553a9fd7477SDave May 554a9fd7477SDave May /* General collection when no order, or neighbour information is provided */ 555a9fd7477SDave May /* 556a9fd7477SDave May User provides context and collect() method 557a9fd7477SDave May Broadcast user context 558a9fd7477SDave May 559a9fd7477SDave May for each context / rank { 560a9fd7477SDave May collect(swarm,context,n,list) 561a9fd7477SDave May } 562a9fd7477SDave May */ 563a9fd7477SDave May #undef __FUNCT__ 564a9fd7477SDave May #define __FUNCT__ "DMSwarmCollect_General" 565a9fd7477SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM dm,PetscErrorCode (*collect)(DM,void*,PetscInt*,PetscInt**),size_t ctx_size,void *ctx,PetscInt *globalsize) 566a9fd7477SDave May { 567a9fd7477SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 568a9fd7477SDave May PetscErrorCode ierr; 569a9fd7477SDave May DataEx de; 570a9fd7477SDave May PetscInt p,r,npoints,n_points_recv; 571a9fd7477SDave May PetscMPIInt commsize,rank; 572a9fd7477SDave May void *point_buffer,*recv_points; 573a9fd7477SDave May void *ctxlist; 574a9fd7477SDave May PetscInt *n2collect,**collectlist; 575a9fd7477SDave May size_t sizeof_dmswarm_point; 576a9fd7477SDave May 577521f74f9SMatthew G. Knepley PetscFunctionBegin; 578a9fd7477SDave May ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr); 579a9fd7477SDave May ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 580a9fd7477SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 581a9fd7477SDave May *globalsize = npoints; 582a9fd7477SDave May /* Broadcast user context */ 583a9fd7477SDave May PetscMalloc(ctx_size*commsize,&ctxlist); 584a9fd7477SDave May ierr = MPI_Allgather(ctx,ctx_size,MPI_CHAR,ctxlist,ctx_size,MPI_CHAR,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 585521f74f9SMatthew G. Knepley ierr = PetscMalloc1(commsize,&n2collect);CHKERRQ(ierr); 586521f74f9SMatthew G. Knepley ierr = PetscMalloc1(commsize,&collectlist);CHKERRQ(ierr); 587a9fd7477SDave May for (r=0; r<commsize; r++) { 588a9fd7477SDave May PetscInt _n2collect; 589a9fd7477SDave May PetscInt *_collectlist; 590a9fd7477SDave May void *_ctx_r; 591a9fd7477SDave May 592a9fd7477SDave May _n2collect = 0; 593a9fd7477SDave May _collectlist = NULL; 594a9fd7477SDave May if (r != rank) { /* don't collect data from yourself */ 595a9fd7477SDave May _ctx_r = (void*)( (char*)ctxlist + r * ctx_size ); 596a9fd7477SDave May ierr = collect(dm,_ctx_r,&_n2collect,&_collectlist);CHKERRQ(ierr); 597a9fd7477SDave May } 598a9fd7477SDave May n2collect[r] = _n2collect; 599a9fd7477SDave May collectlist[r] = _collectlist; 600a9fd7477SDave May } 601521f74f9SMatthew G. Knepley ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr); 602a9fd7477SDave May /* Define topology */ 603a9fd7477SDave May ierr = DataExTopologyInitialize(de);CHKERRQ(ierr); 604a9fd7477SDave May for (r=0; r<commsize; r++) { 605a9fd7477SDave May if (n2collect[r] > 0) { 606a9fd7477SDave May ierr = DataExTopologyAddNeighbour(de,(PetscMPIInt)r);CHKERRQ(ierr); 607a9fd7477SDave May } 608a9fd7477SDave May } 609a9fd7477SDave May ierr = DataExTopologyFinalize(de);CHKERRQ(ierr); 610a9fd7477SDave May /* Define send counts */ 611a9fd7477SDave May ierr = DataExInitializeSendCount(de);CHKERRQ(ierr); 612a9fd7477SDave May for (r=0; r<commsize; r++) { 613a9fd7477SDave May if (n2collect[r] > 0) { 614a9fd7477SDave May ierr = DataExAddToSendCount(de,r,n2collect[r]);CHKERRQ(ierr); 615a9fd7477SDave May } 616a9fd7477SDave May } 617a9fd7477SDave May ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr); 618a9fd7477SDave May /* Pack data */ 619a9fd7477SDave May ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr); 620a9fd7477SDave May ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr); 621a9fd7477SDave May for (r=0; r<commsize; r++) { 622a9fd7477SDave May for (p=0; p<n2collect[r]; p++) { 623a9fd7477SDave May ierr = DataBucketFillPackedArray(swarm->db,collectlist[r][p],point_buffer);CHKERRQ(ierr); 624a9fd7477SDave May /* insert point buffer into the data exchanger */ 625a9fd7477SDave May ierr = DataExPackData(de,r,1,point_buffer);CHKERRQ(ierr); 626a9fd7477SDave May } 627a9fd7477SDave May } 628a9fd7477SDave May ierr = DataExPackFinalize(de);CHKERRQ(ierr); 629a9fd7477SDave May /* Scatter */ 630a9fd7477SDave May ierr = DataExBegin(de);CHKERRQ(ierr); 631a9fd7477SDave May ierr = DataExEnd(de);CHKERRQ(ierr); 632a9fd7477SDave May /* Collect data in DMSwarm container */ 633a9fd7477SDave May ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr); 634a9fd7477SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 635a9fd7477SDave May ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr); 636a9fd7477SDave May for (p=0; p<n_points_recv; p++) { 637a9fd7477SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point ); 638a9fd7477SDave May 639a9fd7477SDave May ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr); 640a9fd7477SDave May } 641a9fd7477SDave May /* Release memory */ 642a9fd7477SDave May for (r=0; r<commsize; r++) { 643a9fd7477SDave May if (collectlist[r]) PetscFree(collectlist[r]); 644a9fd7477SDave May } 645521f74f9SMatthew G. Knepley ierr = PetscFree(collectlist);CHKERRQ(ierr); 646521f74f9SMatthew G. Knepley ierr = PetscFree(n2collect);CHKERRQ(ierr); 647521f74f9SMatthew G. Knepley ierr = PetscFree(ctxlist);CHKERRQ(ierr); 648a9fd7477SDave May ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr); 649a9fd7477SDave May ierr = DataExView(de);CHKERRQ(ierr); 650a9fd7477SDave May ierr = DataExDestroy(de);CHKERRQ(ierr); 651a9fd7477SDave May PetscFunctionReturn(0); 652a9fd7477SDave May } 653a9fd7477SDave May 654