1dfc14de9SMatthew G. Knepley #include <petscsf.h> 208056efcSDave May #include <petscdmswarm.h> 35917a6f0SStefano Zampini #include <petscdmda.h> 4df21e3a8SDave May #include <petsc/private/dmswarmimpl.h> /*I "petscdmswarm.h" I*/ 5279f676cSBarry Smith #include "../src/dm/impls/swarm/data_bucket.h" 6279f676cSBarry Smith #include "../src/dm/impls/swarm/data_ex.h" 7df21e3a8SDave May 8480eef7bSDave May /* 9480eef7bSDave May User loads desired location (MPI rank) into field DMSwarm_rank 10480eef7bSDave May */ 11df21e3a8SDave May PetscErrorCode DMSwarmMigrate_Push_Basic(DM dm,PetscBool remove_sent_points) 12df21e3a8SDave May { 13df21e3a8SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 1477048351SPatrick Sanan DMSwarmDataEx de; 15df21e3a8SDave May PetscInt p,npoints,*rankval,n_points_recv; 16df21e3a8SDave May PetscMPIInt rank,nrank; 17df21e3a8SDave May void *point_buffer,*recv_points; 18df21e3a8SDave May size_t sizeof_dmswarm_point; 19*356ed814SMatthew G. Knepley PetscBool debug = PETSC_FALSE; 20df21e3a8SDave May 21521f74f9SMatthew G. Knepley PetscFunctionBegin; 229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank)); 23df21e3a8SDave May 249566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL)); 259566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval)); 269566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm),0, &de)); 279566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExTopologyInitialize(de)); 28521f74f9SMatthew G. Knepley for (p = 0; p < npoints; ++p) { 29df21e3a8SDave May nrank = rankval[p]; 30df21e3a8SDave May if (nrank != rank) { 319566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExTopologyAddNeighbour(de,nrank)); 32df21e3a8SDave May } 33df21e3a8SDave May } 349566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExTopologyFinalize(de)); 359566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExInitializeSendCount(de)); 36df21e3a8SDave May for (p=0; p<npoints; p++) { 37df21e3a8SDave May nrank = rankval[p]; 38df21e3a8SDave May if (nrank != rank) { 399566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExAddToSendCount(de,nrank,1)); 40df21e3a8SDave May } 41df21e3a8SDave May } 429566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExFinalizeSendCount(de)); 439566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer)); 449566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackInitialize(de,sizeof_dmswarm_point)); 45df21e3a8SDave May for (p=0; p<npoints; p++) { 46df21e3a8SDave May nrank = rankval[p]; 47df21e3a8SDave May if (nrank != rank) { 48df21e3a8SDave May /* copy point into buffer */ 499566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db,p,point_buffer)); 5077048351SPatrick Sanan /* insert point buffer into DMSwarmDataExchanger */ 519566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackData(de,nrank,1,point_buffer)); 52df21e3a8SDave May } 53df21e3a8SDave May } 549566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackFinalize(de)); 559566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval)); 56df21e3a8SDave May 57df21e3a8SDave May if (remove_sent_points) { 5877048351SPatrick Sanan DMSwarmDataField gfield; 5922a417f9SDave May 609566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,DMSwarmField_rank,&gfield)); 619566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldGetAccess(gfield)); 629566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldGetEntries(gfield,(void**)&rankval)); 6322a417f9SDave May 64df21e3a8SDave May /* remove points which left processor */ 659566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL)); 66df21e3a8SDave May for (p=0; p<npoints; p++) { 67df21e3a8SDave May nrank = rankval[p]; 68df21e3a8SDave May if (nrank != rank) { 69df21e3a8SDave May /* kill point */ 709566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldRestoreAccess(gfield)); 7122a417f9SDave May 729566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db,p)); 7322a417f9SDave May 749566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL)); /* you need to update npoints as the list size decreases! */ 759566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldGetAccess(gfield)); 769566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldGetEntries(gfield,(void**)&rankval)); 77df21e3a8SDave May p--; /* check replacement point */ 78df21e3a8SDave May } 79df21e3a8SDave May } 809566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldRestoreEntries(gfield,(void**)&rankval)); 819566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldRestoreAccess(gfield)); 82df21e3a8SDave May } 839566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExBegin(de)); 849566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExEnd(de)); 859566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExGetRecvData(de,&n_points_recv,(void**)&recv_points)); 869566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL)); 879566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(swarm->db,npoints + n_points_recv,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 88df21e3a8SDave May for (p=0; p<n_points_recv; p++) { 89df21e3a8SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point); 90df21e3a8SDave May 919566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db,npoints+p,data_p)); 92df21e3a8SDave May } 93*356ed814SMatthew G. Knepley if (debug) PetscCall(DMSwarmDataExView(de)); 949566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db,&point_buffer)); 959566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExDestroy(de)); 96df21e3a8SDave May PetscFunctionReturn(0); 97df21e3a8SDave May } 982712d1f2SDave May 99889dbfe5SDave May PetscErrorCode DMSwarmMigrate_DMNeighborScatter(DM dm,DM dmcell,PetscBool remove_sent_points,PetscInt *npoints_prior_migration) 10040c453e9SDave May { 10140c453e9SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 10277048351SPatrick Sanan DMSwarmDataEx 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; 1097c6d1d28SDave May PetscMPIInt mynneigh,*myneigh; 11040c453e9SDave May 111521f74f9SMatthew G. Knepley PetscFunctionBegin; 1129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank)); 1139566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL)); 1149566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval)); 1159566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm),0,&de)); 1169566063dSJacob Faibussowitsch PetscCall(DMGetNeighbors(dmcell,&nneighbors,&neighbourranks)); 1179566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExTopologyInitialize(de)); 11840c453e9SDave May for (r=0; r<nneighbors; r++) { 11940c453e9SDave May _rank = neighbourranks[r]; 12040c453e9SDave May if ((_rank != rank) && (_rank > 0)) { 1219566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExTopologyAddNeighbour(de,_rank)); 12240c453e9SDave May } 12340c453e9SDave May } 1249566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExTopologyFinalize(de)); 1259566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExTopologyGetNeighbours(de,&mynneigh,&myneigh)); 1269566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExInitializeSendCount(de)); 12740c453e9SDave May for (p=0; p<npoints; p++) { 128f954cb40SDave May if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) { 1297c6d1d28SDave May for (r=0; r<mynneigh; r++) { 1307c6d1d28SDave May _rank = myneigh[r]; 1319566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExAddToSendCount(de,_rank,1)); 13240c453e9SDave May } 13340c453e9SDave May } 13440c453e9SDave May } 1359566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExFinalizeSendCount(de)); 1369566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer)); 1379566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackInitialize(de,sizeof_dmswarm_point)); 13840c453e9SDave May for (p=0; p<npoints; p++) { 139f954cb40SDave May if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) { 1407c6d1d28SDave May for (r=0; r<mynneigh; r++) { 1417c6d1d28SDave May _rank = myneigh[r]; 14240c453e9SDave May /* copy point into buffer */ 1439566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db,p,point_buffer)); 14477048351SPatrick Sanan /* insert point buffer into DMSwarmDataExchanger */ 1459566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackData(de,_rank,1,point_buffer)); 14640c453e9SDave May } 14740c453e9SDave May } 14840c453e9SDave May } 1499566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackFinalize(de)); 1509566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval)); 15140c453e9SDave May if (remove_sent_points) { 15277048351SPatrick Sanan DMSwarmDataField PField; 1537c6d1d28SDave May 1549566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,DMSwarmField_rank,&PField)); 1559566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldGetEntries(PField,(void**)&rankval)); 15640c453e9SDave May /* remove points which left processor */ 1579566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL)); 15840c453e9SDave May for (p=0; p<npoints; p++) { 159f954cb40SDave May if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) { 16040c453e9SDave May /* kill point */ 1619566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db,p)); 1629566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL)); /* you need to update npoints as the list size decreases! */ 1639566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldGetEntries(PField,(void**)&rankval)); /* update date point increase realloc performed */ 16440c453e9SDave May p--; /* check replacement point */ 16540c453e9SDave May } 16640c453e9SDave May } 16740c453e9SDave May } 1689566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,npoints_prior_migration,NULL,NULL)); 1699566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExBegin(de)); 1709566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExEnd(de)); 1719566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExGetRecvData(de,&n_points_recv,(void**)&recv_points)); 1729566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL)); 1739566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(swarm->db,npoints + n_points_recv,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 17440c453e9SDave May for (p=0; p<n_points_recv; p++) { 17540c453e9SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point); 17640c453e9SDave May 1779566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db,npoints+p,data_p)); 17840c453e9SDave May } 1799566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db,&point_buffer)); 1809566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExDestroy(de)); 18140c453e9SDave May PetscFunctionReturn(0); 18240c453e9SDave May } 183480eef7bSDave May 18408056efcSDave May PetscErrorCode DMSwarmMigrate_CellDMScatter(DM dm,PetscBool remove_sent_points) 185480eef7bSDave May { 186480eef7bSDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 1876fbf25f8SDave May PetscInt p,npoints,npointsg=0,npoints2,npoints2g,*rankval,npoints_prior_migration; 188bbe8250bSMatthew G. Knepley PetscSF sfcell = NULL; 189dfc14de9SMatthew G. Knepley const PetscSFNode *LA_sfcell; 190480eef7bSDave May DM dmcell; 191480eef7bSDave May Vec pos; 19240c453e9SDave May PetscBool error_check = swarm->migrate_error_on_missing_point; 193e4fbd051SBarry Smith PetscMPIInt size,rank; 194480eef7bSDave May 195521f74f9SMatthew G. Knepley PetscFunctionBegin; 1969566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm,&dmcell)); 19728b400f6SJacob Faibussowitsch PetscCheck(dmcell,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided"); 198480eef7bSDave May 1999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size)); 2009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank)); 2017c6d1d28SDave May 20243a82f2bSDave May #if 1 20343a82f2bSDave May { 20443a82f2bSDave May PetscInt *p_cellid; 20543a82f2bSDave May PetscInt npoints_curr,range = 0; 20643a82f2bSDave May PetscSFNode *sf_cells; 20743a82f2bSDave May 2089566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints_curr,NULL,NULL)); 2099566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints_curr, &sf_cells)); 21043a82f2bSDave May 2119566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval)); 2129566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&p_cellid)); 21343a82f2bSDave May for (p=0; p<npoints_curr; p++) { 21443a82f2bSDave May 21543a82f2bSDave May sf_cells[p].rank = 0; 21643a82f2bSDave May sf_cells[p].index = p_cellid[p]; 21743a82f2bSDave May if (p_cellid[p] > range) { 21843a82f2bSDave May range = p_cellid[p]; 21943a82f2bSDave May } 22043a82f2bSDave May } 2219566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&p_cellid)); 2229566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval)); 22343a82f2bSDave May 2249566063dSJacob Faibussowitsch /* PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm),&sfcell)); */ 2259566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(PETSC_COMM_SELF,&sfcell)); 2269566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sfcell, range, npoints_curr, NULL, PETSC_OWN_POINTER, sf_cells, PETSC_OWN_POINTER)); 22743a82f2bSDave May } 22843a82f2bSDave May #endif 22943a82f2bSDave May 2309566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateLocalVectorFromField(dm, DMSwarmPICField_coor, &pos)); 2319566063dSJacob Faibussowitsch PetscCall(DMLocatePoints(dmcell, pos, DM_POINTLOCATION_NONE, &sfcell)); 2329566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyLocalVectorFromField(dm, DMSwarmPICField_coor, &pos)); 233480eef7bSDave May 23440c453e9SDave May if (error_check) { 2359566063dSJacob Faibussowitsch PetscCall(DMSwarmGetSize(dm,&npointsg)); 23640c453e9SDave May } 2379566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL)); 2389566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval)); 2399566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell)); 240480eef7bSDave May for (p=0; p<npoints; p++) { 241dfc14de9SMatthew G. Knepley rankval[p] = LA_sfcell[p].index; 242480eef7bSDave May } 2439566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval)); 2449566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfcell)); 245480eef7bSDave May 246e4fbd051SBarry Smith if (size > 1) { 2479566063dSJacob Faibussowitsch PetscCall(DMSwarmMigrate_DMNeighborScatter(dm,dmcell,remove_sent_points,&npoints_prior_migration)); 2486fbf25f8SDave May } else { 24977048351SPatrick Sanan DMSwarmDataField PField; 2500ed23c7fSDave May PetscInt npoints_curr; 2510ed23c7fSDave May 2520ed23c7fSDave May /* remove points which the domain */ 2539566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,DMSwarmField_rank,&PField)); 2549566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldGetEntries(PField,(void**)&rankval)); 2550ed23c7fSDave May 2569566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints_curr,NULL,NULL)); 2570ed23c7fSDave May for (p=0; p<npoints_curr; p++) { 2580ed23c7fSDave May if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) { 2590ed23c7fSDave May /* kill point */ 2609566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db,p)); 2619566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints_curr,NULL,NULL)); /* you need to update npoints as the list size decreases! */ 2629566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldGetEntries(PField,(void**)&rankval)); /* update date point increase realloc performed */ 2630ed23c7fSDave May p--; /* check replacement point */ 2640ed23c7fSDave May } 2650ed23c7fSDave May } 2669566063dSJacob Faibussowitsch PetscCall(DMSwarmGetSize(dm,&npoints_prior_migration)); 2670ed23c7fSDave May 2686fbf25f8SDave May } 269480eef7bSDave May 2702d4ee042Sprj- /* locate points newly received */ 2719566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL)); 272009b43efSDave May 2737c6d1d28SDave May #if 0 2742d4ee042Sprj- { /* safe alternative - however this performs two point locations on: (i) the initial points set and; (ii) the (initial + received) point set */ 2757c6d1d28SDave May PetscScalar *LA_coor; 2767c6d1d28SDave May PetscInt bs; 27777048351SPatrick Sanan DMSwarmDataField PField; 2787c6d1d28SDave May 2799566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor)); 2809566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*npoints2,(const PetscScalar*)LA_coor,&pos)); 2819566063dSJacob Faibussowitsch PetscCall(DMLocatePoints(dmcell,pos,DM_POINTLOCATION_NONE,&sfcell)); 2827c6d1d28SDave May 2839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pos)); 2849566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor)); 2857c6d1d28SDave May 2869566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell)); 2879566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval)); 2887c6d1d28SDave May for (p=0; p<npoints2; p++) { 289dfc14de9SMatthew G. Knepley rankval[p] = LA_sfcell[p].index; 2907c6d1d28SDave May } 2919566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfcell)); 2929566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval)); 29340c453e9SDave May 2947c6d1d28SDave May /* remove points which left processor */ 2959566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,DMSwarmField_rank,&PField)); 2969566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldGetEntries(PField,(void**)&rankval)); 2977c6d1d28SDave May 2989566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL)); 2997c6d1d28SDave May for (p=0; p<npoints2; p++) { 300f954cb40SDave May if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) { 3017c6d1d28SDave May /* kill point */ 3029566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db,p)); 3039566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL)); /* you need to update npoints as the list size decreases! */ 3049566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldGetEntries(PField,(void**)&rankval)); /* update date point increase realloc performed */ 3057c6d1d28SDave May p--; /* check replacement point */ 3067c6d1d28SDave May } 3077c6d1d28SDave May } 30840c453e9SDave May } 309009b43efSDave May #endif 310009b43efSDave May 3115627991aSBarry Smith { /* perform two point locations: (i) on the initial points set prior to communication; and (ii) on the new (received) points */ 312009b43efSDave May PetscScalar *LA_coor; 313009b43efSDave May PetscInt npoints_from_neighbours,bs; 31477048351SPatrick Sanan DMSwarmDataField PField; 315009b43efSDave May 316009b43efSDave May npoints_from_neighbours = npoints2 - npoints_prior_migration; 317009b43efSDave May 3189566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor)); 3199566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*npoints_from_neighbours,(const PetscScalar*)&LA_coor[bs*npoints_prior_migration],&pos)); 320009b43efSDave May 3219566063dSJacob Faibussowitsch PetscCall(DMLocatePoints(dmcell,pos,DM_POINTLOCATION_NONE,&sfcell)); 322009b43efSDave May 3239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pos)); 3249566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor)); 325009b43efSDave May 3269566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell)); 3279566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval)); 328009b43efSDave May for (p=0; p<npoints_from_neighbours; p++) { 329009b43efSDave May rankval[npoints_prior_migration + p] = LA_sfcell[p].index; 330009b43efSDave May } 3319566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval)); 3329566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfcell)); 333009b43efSDave May 334009b43efSDave May /* remove points which left processor */ 3359566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,DMSwarmField_rank,&PField)); 3369566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldGetEntries(PField,(void**)&rankval)); 337009b43efSDave May 3389566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL)); 339009b43efSDave May for (p=npoints_prior_migration; p<npoints2; p++) { 340009b43efSDave May if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) { 341009b43efSDave May /* kill point */ 3429566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db,p)); 3439566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL)); /* you need to update npoints as the list size decreases! */ 3449566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldGetEntries(PField,(void**)&rankval)); /* update date point increase realloc performed */ 345009b43efSDave May p--; /* check replacement point */ 346009b43efSDave May } 347009b43efSDave May } 348009b43efSDave May } 349009b43efSDave May 3503151f1d5SDave May { 3513151f1d5SDave May PetscInt *p_cellid; 3523151f1d5SDave May 3539566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL)); 3549566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval)); 3559566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&p_cellid)); 3563151f1d5SDave May for (p=0; p<npoints2; p++) { 3573151f1d5SDave May p_cellid[p] = rankval[p]; 3583151f1d5SDave May } 3599566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&p_cellid)); 3609566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval)); 3613151f1d5SDave May } 3623151f1d5SDave May 36340c453e9SDave May /* check for error on removed points */ 36440c453e9SDave May if (error_check) { 3659566063dSJacob Faibussowitsch PetscCall(DMSwarmGetSize(dm,&npoints2g)); 3662c71b3e2SJacob Faibussowitsch PetscCheckFalse(npointsg != npoints2g,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Points from the DMSwarm must remain constant during migration (initial %D - final %D)",npointsg,npoints2g); 36740c453e9SDave May } 368480eef7bSDave May PetscFunctionReturn(0); 369480eef7bSDave May } 370480eef7bSDave May 37108056efcSDave May PetscErrorCode DMSwarmMigrate_CellDMExact(DM dm,PetscBool remove_sent_points) 37208056efcSDave May { 373521f74f9SMatthew G. Knepley PetscFunctionBegin; 37408056efcSDave May PetscFunctionReturn(0); 37508056efcSDave May } 37608056efcSDave May 377480eef7bSDave May /* 378480eef7bSDave May Redundant as this assumes points can only be sent to a single rank 379480eef7bSDave May */ 3802712d1f2SDave May PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt *globalsize) 3812712d1f2SDave May { 3822712d1f2SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 38377048351SPatrick Sanan DMSwarmDataEx de; 3842712d1f2SDave May PetscInt p,npoints,*rankval,n_points_recv; 3852712d1f2SDave May PetscMPIInt rank,nrank,negrank; 3862712d1f2SDave May void *point_buffer,*recv_points; 3872712d1f2SDave May size_t sizeof_dmswarm_point; 3882712d1f2SDave May 389521f74f9SMatthew G. Knepley PetscFunctionBegin; 3909566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank)); 3919566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL)); 3922712d1f2SDave May *globalsize = npoints; 3939566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval)); 3949566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm),0,&de)); 3959566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExTopologyInitialize(de)); 3962712d1f2SDave May for (p=0; p<npoints; p++) { 3972712d1f2SDave May negrank = rankval[p]; 3982712d1f2SDave May if (negrank < 0) { 3992712d1f2SDave May nrank = -negrank - 1; 4009566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExTopologyAddNeighbour(de,nrank)); 4012712d1f2SDave May } 4022712d1f2SDave May } 4039566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExTopologyFinalize(de)); 4049566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExInitializeSendCount(de)); 4052712d1f2SDave May for (p=0; p<npoints; p++) { 4062712d1f2SDave May negrank = rankval[p]; 4072712d1f2SDave May if (negrank < 0) { 4082712d1f2SDave May nrank = -negrank - 1; 4099566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExAddToSendCount(de,nrank,1)); 4102712d1f2SDave May } 4112712d1f2SDave May } 4129566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExFinalizeSendCount(de)); 4139566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer)); 4149566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackInitialize(de,sizeof_dmswarm_point)); 4152712d1f2SDave May for (p=0; p<npoints; p++) { 4162712d1f2SDave May negrank = rankval[p]; 4172712d1f2SDave May if (negrank < 0) { 4182712d1f2SDave May nrank = -negrank - 1; 4192712d1f2SDave May rankval[p] = nrank; 4202712d1f2SDave May /* copy point into buffer */ 4219566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db,p,point_buffer)); 42277048351SPatrick Sanan /* insert point buffer into DMSwarmDataExchanger */ 4239566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackData(de,nrank,1,point_buffer)); 4242712d1f2SDave May rankval[p] = negrank; 4252712d1f2SDave May } 4262712d1f2SDave May } 4279566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackFinalize(de)); 4289566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval)); 4299566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExBegin(de)); 4309566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExEnd(de)); 4319566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExGetRecvData(de,&n_points_recv,(void**)&recv_points)); 4329566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL)); 4339566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(swarm->db,npoints + n_points_recv,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 4342712d1f2SDave May for (p=0; p<n_points_recv; p++) { 4352712d1f2SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point); 4362712d1f2SDave May 4379566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db,npoints+p,data_p)); 4382712d1f2SDave May } 4399566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExView(de)); 4409566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db,&point_buffer)); 4419566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExDestroy(de)); 4422712d1f2SDave May PetscFunctionReturn(0); 4432712d1f2SDave May } 444b16650c8SDave May 445b16650c8SDave May typedef struct { 446b16650c8SDave May PetscMPIInt owner_rank; 447b16650c8SDave May PetscReal min[3],max[3]; 448b16650c8SDave May } CollectBBox; 449b16650c8SDave May 450fe39f135SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM dm,PetscInt *globalsize) 451b16650c8SDave May { 452b16650c8SDave May DM_Swarm * swarm = (DM_Swarm*)dm->data; 45377048351SPatrick Sanan DMSwarmDataEx de; 454b16650c8SDave May PetscInt p,pk,npoints,*rankval,n_points_recv,n_bbox_recv,dim,neighbour_cells; 455b16650c8SDave May PetscMPIInt rank,nrank; 456b16650c8SDave May void *point_buffer,*recv_points; 457b16650c8SDave May size_t sizeof_dmswarm_point,sizeof_bbox_ctx; 458b16650c8SDave May PetscBool isdmda; 459b16650c8SDave May CollectBBox *bbox,*recv_bbox; 460b16650c8SDave May const PetscMPIInt *dmneighborranks; 461b16650c8SDave May DM dmcell; 462b16650c8SDave May 463521f74f9SMatthew G. Knepley PetscFunctionBegin; 4649566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank)); 465b16650c8SDave May 4669566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm,&dmcell)); 46728b400f6SJacob Faibussowitsch PetscCheck(dmcell,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided"); 468b16650c8SDave May isdmda = PETSC_FALSE; 469b16650c8SDave May PetscObjectTypeCompare((PetscObject)dmcell,DMDA,&isdmda); 47028b400f6SJacob Faibussowitsch PetscCheck(isdmda,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only DMDA support for CollectBoundingBox"); 471b16650c8SDave May 4729566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm,&dim)); 473b16650c8SDave May sizeof_bbox_ctx = sizeof(CollectBBox); 474b16650c8SDave May PetscMalloc1(1,&bbox); 475b16650c8SDave May bbox->owner_rank = rank; 476b16650c8SDave May 477b16650c8SDave May /* compute the bounding box based on the overlapping / stenctil size */ 478b16650c8SDave May { 479b16650c8SDave May Vec lcoor; 480b16650c8SDave May 4819566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmcell,&lcoor)); 482fe39f135SDave May if (dim >= 1) { 4839566063dSJacob Faibussowitsch PetscCall(VecStrideMin(lcoor,0,NULL,&bbox->min[0])); 4849566063dSJacob Faibussowitsch PetscCall(VecStrideMax(lcoor,0,NULL,&bbox->max[0])); 485fe39f135SDave May } 486fe39f135SDave May if (dim >= 2) { 4879566063dSJacob Faibussowitsch PetscCall(VecStrideMin(lcoor,1,NULL,&bbox->min[1])); 4889566063dSJacob Faibussowitsch PetscCall(VecStrideMax(lcoor,1,NULL,&bbox->max[1])); 489b16650c8SDave May } 490fe39f135SDave May if (dim == 3) { 4919566063dSJacob Faibussowitsch PetscCall(VecStrideMin(lcoor,2,NULL,&bbox->min[2])); 4929566063dSJacob Faibussowitsch PetscCall(VecStrideMax(lcoor,2,NULL,&bbox->max[2])); 493fe39f135SDave May } 494fe39f135SDave May } 4959566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL)); 496b16650c8SDave May *globalsize = npoints; 4979566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval)); 4989566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm),0,&de)); 499b16650c8SDave May /* use DMDA neighbours */ 5009566063dSJacob Faibussowitsch PetscCall(DMDAGetNeighbors(dmcell,&dmneighborranks)); 5018dbd68bcSDave May if (dim == 1) { 5028dbd68bcSDave May neighbour_cells = 3; 5038dbd68bcSDave May } else if (dim == 2) { 5048dbd68bcSDave May neighbour_cells = 9; 5058dbd68bcSDave May } else { 5068dbd68bcSDave May neighbour_cells = 27; 5078dbd68bcSDave May } 5089566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExTopologyInitialize(de)); 509b16650c8SDave May for (p=0; p<neighbour_cells; p++) { 510b16650c8SDave May if ((dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank)) { 5119566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExTopologyAddNeighbour(de,dmneighborranks[p])); 512b16650c8SDave May } 513b16650c8SDave May } 5149566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExTopologyFinalize(de)); 5159566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExInitializeSendCount(de)); 516b16650c8SDave May for (p=0; p<neighbour_cells; p++) { 517b16650c8SDave May if ((dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank)) { 5189566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExAddToSendCount(de,dmneighborranks[p],1)); 519b16650c8SDave May } 520b16650c8SDave May } 5219566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExFinalizeSendCount(de)); 522b16650c8SDave May /* send bounding boxes */ 5239566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackInitialize(de,sizeof_bbox_ctx)); 524b16650c8SDave May for (p=0; p<neighbour_cells; p++) { 525b16650c8SDave May nrank = dmneighborranks[p]; 526b16650c8SDave May if ((nrank >= 0) && (nrank != rank)) { 52777048351SPatrick Sanan /* insert bbox buffer into DMSwarmDataExchanger */ 5289566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackData(de,nrank,1,bbox)); 529b16650c8SDave May } 530b16650c8SDave May } 5319566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackFinalize(de)); 532b16650c8SDave May /* recv bounding boxes */ 5339566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExBegin(de)); 5349566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExEnd(de)); 5359566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExGetRecvData(de,&n_bbox_recv,(void**)&recv_bbox)); 536298827fbSBarry Smith /* Wrong, should not be using PETSC_COMM_WORLD */ 537b16650c8SDave May for (p=0; p<n_bbox_recv; p++) { 5389566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[rank %d]: box from %d : range[%+1.4e,%+1.4e]x[%+1.4e,%+1.4e]\n",rank,recv_bbox[p].owner_rank, 539b122ec5aSJacob Faibussowitsch (double)recv_bbox[p].min[0],(double)recv_bbox[p].max[0],(double)recv_bbox[p].min[1],(double)recv_bbox[p].max[1])); 540b16650c8SDave May } 5419566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,stdout)); 542b16650c8SDave May /* of course this is stupid as this "generic" function should have a better way to know what the coordinates are called */ 5439566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExInitializeSendCount(de)); 544b16650c8SDave May for (pk=0; pk<n_bbox_recv; pk++) { 545b16650c8SDave May PetscReal *array_x,*array_y; 546b16650c8SDave May 5479566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x)); 5489566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y)); 549b16650c8SDave May for (p=0; p<npoints; p++) { 550b16650c8SDave May if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0])) { 551b16650c8SDave May if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1])) { 5529566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExAddToSendCount(de,recv_bbox[pk].owner_rank,1)); 553b16650c8SDave May } 554b16650c8SDave May } 555b16650c8SDave May } 5569566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y)); 5579566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x)); 558b16650c8SDave May } 5599566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExFinalizeSendCount(de)); 5609566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer)); 5619566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackInitialize(de,sizeof_dmswarm_point)); 562b16650c8SDave May for (pk=0; pk<n_bbox_recv; pk++) { 563b16650c8SDave May PetscReal *array_x,*array_y; 564b16650c8SDave May 5659566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x)); 5669566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y)); 567b16650c8SDave May for (p=0; p<npoints; p++) { 568b16650c8SDave May if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0])) { 569b16650c8SDave May if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1])) { 570521f74f9SMatthew G. Knepley /* copy point into buffer */ 5719566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db,p,point_buffer)); 57277048351SPatrick Sanan /* insert point buffer into DMSwarmDataExchanger */ 5739566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackData(de,recv_bbox[pk].owner_rank,1,point_buffer)); 574b16650c8SDave May } 575b16650c8SDave May } 576b16650c8SDave May } 5779566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y)); 5789566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x)); 579b16650c8SDave May } 5809566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackFinalize(de)); 5819566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval)); 5829566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExBegin(de)); 5839566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExEnd(de)); 5849566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExGetRecvData(de,&n_points_recv,(void**)&recv_points)); 5859566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL)); 5869566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(swarm->db,npoints + n_points_recv,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 587b16650c8SDave May for (p=0; p<n_points_recv; p++) { 588b16650c8SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point); 589b16650c8SDave May 5909566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db,npoints+p,data_p)); 591b16650c8SDave May } 5929566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db,&point_buffer)); 593b16650c8SDave May PetscFree(bbox); 5949566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExView(de)); 5959566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExDestroy(de)); 596b16650c8SDave May PetscFunctionReturn(0); 597b16650c8SDave May } 598a9fd7477SDave May 599a9fd7477SDave May /* General collection when no order, or neighbour information is provided */ 600a9fd7477SDave May /* 601a9fd7477SDave May User provides context and collect() method 602a9fd7477SDave May Broadcast user context 603a9fd7477SDave May 604a9fd7477SDave May for each context / rank { 605a9fd7477SDave May collect(swarm,context,n,list) 606a9fd7477SDave May } 607a9fd7477SDave May */ 608a9fd7477SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM dm,PetscErrorCode (*collect)(DM,void*,PetscInt*,PetscInt**),size_t ctx_size,void *ctx,PetscInt *globalsize) 609a9fd7477SDave May { 610a9fd7477SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 61177048351SPatrick Sanan DMSwarmDataEx de; 612a9fd7477SDave May PetscInt p,r,npoints,n_points_recv; 613e4fbd051SBarry Smith PetscMPIInt size,rank; 614a9fd7477SDave May void *point_buffer,*recv_points; 615a9fd7477SDave May void *ctxlist; 616a9fd7477SDave May PetscInt *n2collect,**collectlist; 617a9fd7477SDave May size_t sizeof_dmswarm_point; 618a9fd7477SDave May 619521f74f9SMatthew G. Knepley PetscFunctionBegin; 6209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size)); 6219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank)); 6229566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL)); 623a9fd7477SDave May *globalsize = npoints; 624a9fd7477SDave May /* Broadcast user context */ 625e4fbd051SBarry Smith PetscMalloc(ctx_size*size,&ctxlist); 6269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(ctx,ctx_size,MPI_CHAR,ctxlist,ctx_size,MPI_CHAR,PetscObjectComm((PetscObject)dm))); 6279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size,&n2collect)); 6289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size,&collectlist)); 629e4fbd051SBarry Smith for (r=0; r<size; r++) { 630a9fd7477SDave May PetscInt _n2collect; 631a9fd7477SDave May PetscInt *_collectlist; 632a9fd7477SDave May void *_ctx_r; 633a9fd7477SDave May 634a9fd7477SDave May _n2collect = 0; 635a9fd7477SDave May _collectlist = NULL; 636a9fd7477SDave May if (r != rank) { /* don't collect data from yourself */ 637a9fd7477SDave May _ctx_r = (void*)( (char*)ctxlist + r * ctx_size); 6389566063dSJacob Faibussowitsch PetscCall(collect(dm,_ctx_r,&_n2collect,&_collectlist)); 639a9fd7477SDave May } 640a9fd7477SDave May n2collect[r] = _n2collect; 641a9fd7477SDave May collectlist[r] = _collectlist; 642a9fd7477SDave May } 6439566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm),0,&de)); 644a9fd7477SDave May /* Define topology */ 6459566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExTopologyInitialize(de)); 646e4fbd051SBarry Smith for (r=0; r<size; r++) { 647a9fd7477SDave May if (n2collect[r] > 0) { 6489566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExTopologyAddNeighbour(de,(PetscMPIInt)r)); 649a9fd7477SDave May } 650a9fd7477SDave May } 6519566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExTopologyFinalize(de)); 652a9fd7477SDave May /* Define send counts */ 6539566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExInitializeSendCount(de)); 654e4fbd051SBarry Smith for (r=0; r<size; r++) { 655a9fd7477SDave May if (n2collect[r] > 0) { 6569566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExAddToSendCount(de,r,n2collect[r])); 657a9fd7477SDave May } 658a9fd7477SDave May } 6599566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExFinalizeSendCount(de)); 660a9fd7477SDave May /* Pack data */ 6619566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer)); 6629566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackInitialize(de,sizeof_dmswarm_point)); 663e4fbd051SBarry Smith for (r=0; r<size; r++) { 664a9fd7477SDave May for (p=0; p<n2collect[r]; p++) { 6659566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db,collectlist[r][p],point_buffer)); 666a9fd7477SDave May /* insert point buffer into the data exchanger */ 6679566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackData(de,r,1,point_buffer)); 668a9fd7477SDave May } 669a9fd7477SDave May } 6709566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackFinalize(de)); 671a9fd7477SDave May /* Scatter */ 6729566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExBegin(de)); 6739566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExEnd(de)); 674a9fd7477SDave May /* Collect data in DMSwarm container */ 6759566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExGetRecvData(de,&n_points_recv,(void**)&recv_points)); 6769566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL)); 6779566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(swarm->db,npoints + n_points_recv,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 678a9fd7477SDave May for (p=0; p<n_points_recv; p++) { 679a9fd7477SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point); 680a9fd7477SDave May 6819566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db,npoints+p,data_p)); 682a9fd7477SDave May } 683a9fd7477SDave May /* Release memory */ 684e4fbd051SBarry Smith for (r=0; r<size; r++) { 685a9fd7477SDave May if (collectlist[r]) PetscFree(collectlist[r]); 686a9fd7477SDave May } 6879566063dSJacob Faibussowitsch PetscCall(PetscFree(collectlist)); 6889566063dSJacob Faibussowitsch PetscCall(PetscFree(n2collect)); 6899566063dSJacob Faibussowitsch PetscCall(PetscFree(ctxlist)); 6909566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db,&point_buffer)); 6919566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExView(de)); 6929566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExDestroy(de)); 693a9fd7477SDave May PetscFunctionReturn(0); 694a9fd7477SDave May } 695