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; 14df21e3a8SDave May PetscErrorCode ierr; 1577048351SPatrick Sanan DMSwarmDataEx de; 16df21e3a8SDave May PetscInt p,npoints,*rankval,n_points_recv; 17df21e3a8SDave May PetscMPIInt rank,nrank; 18df21e3a8SDave May void *point_buffer,*recv_points; 19df21e3a8SDave May size_t sizeof_dmswarm_point; 20df21e3a8SDave May 21521f74f9SMatthew G. Knepley PetscFunctionBegin; 22df21e3a8SDave May ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 23df21e3a8SDave May 2477048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 2508056efcSDave May ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 2677048351SPatrick Sanan ierr = DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm),0, &de);CHKERRQ(ierr); 2777048351SPatrick Sanan ierr = DMSwarmDataExTopologyInitialize(de);CHKERRQ(ierr); 28521f74f9SMatthew G. Knepley for (p = 0; p < npoints; ++p) { 29df21e3a8SDave May nrank = rankval[p]; 30df21e3a8SDave May if (nrank != rank) { 3177048351SPatrick Sanan ierr = DMSwarmDataExTopologyAddNeighbour(de,nrank);CHKERRQ(ierr); 32df21e3a8SDave May } 33df21e3a8SDave May } 3477048351SPatrick Sanan ierr = DMSwarmDataExTopologyFinalize(de);CHKERRQ(ierr); 3577048351SPatrick Sanan ierr = DMSwarmDataExInitializeSendCount(de);CHKERRQ(ierr); 36df21e3a8SDave May for (p=0; p<npoints; p++) { 37df21e3a8SDave May nrank = rankval[p]; 38df21e3a8SDave May if (nrank != rank) { 3977048351SPatrick Sanan ierr = DMSwarmDataExAddToSendCount(de,nrank,1);CHKERRQ(ierr); 40df21e3a8SDave May } 41df21e3a8SDave May } 4277048351SPatrick Sanan ierr = DMSwarmDataExFinalizeSendCount(de);CHKERRQ(ierr); 4377048351SPatrick Sanan ierr = DMSwarmDataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr); 4477048351SPatrick Sanan ierr = DMSwarmDataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr); 45df21e3a8SDave May for (p=0; p<npoints; p++) { 46df21e3a8SDave May nrank = rankval[p]; 47df21e3a8SDave May if (nrank != rank) { 48df21e3a8SDave May /* copy point into buffer */ 4977048351SPatrick Sanan ierr = DMSwarmDataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr); 5077048351SPatrick Sanan /* insert point buffer into DMSwarmDataExchanger */ 5177048351SPatrick Sanan ierr = DMSwarmDataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr); 52df21e3a8SDave May } 53df21e3a8SDave May } 5477048351SPatrick Sanan ierr = DMSwarmDataExPackFinalize(de);CHKERRQ(ierr); 5522a417f9SDave May ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 56df21e3a8SDave May 57df21e3a8SDave May if (remove_sent_points) { 5877048351SPatrick Sanan DMSwarmDataField gfield; 5922a417f9SDave May 6077048351SPatrick Sanan ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,DMSwarmField_rank,&gfield);CHKERRQ(ierr); 6177048351SPatrick Sanan ierr = DMSwarmDataFieldGetAccess(gfield);CHKERRQ(ierr); 6277048351SPatrick Sanan ierr = DMSwarmDataFieldGetEntries(gfield,(void**)&rankval);CHKERRQ(ierr); 6322a417f9SDave May 64df21e3a8SDave May /* remove points which left processor */ 6577048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 66df21e3a8SDave May for (p=0; p<npoints; p++) { 67df21e3a8SDave May nrank = rankval[p]; 68df21e3a8SDave May if (nrank != rank) { 69df21e3a8SDave May /* kill point */ 7077048351SPatrick Sanan ierr = DMSwarmDataFieldRestoreAccess(gfield);CHKERRQ(ierr); 7122a417f9SDave May 7277048351SPatrick Sanan ierr = DMSwarmDataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr); 7322a417f9SDave May 7477048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */ 7577048351SPatrick Sanan ierr = DMSwarmDataFieldGetAccess(gfield);CHKERRQ(ierr); 7677048351SPatrick Sanan ierr = DMSwarmDataFieldGetEntries(gfield,(void**)&rankval);CHKERRQ(ierr); 77df21e3a8SDave May p--; /* check replacement point */ 78df21e3a8SDave May } 79df21e3a8SDave May } 8077048351SPatrick Sanan ierr = DMSwarmDataFieldRestoreEntries(gfield,(void**)&rankval);CHKERRQ(ierr); 8177048351SPatrick Sanan ierr = DMSwarmDataFieldRestoreAccess(gfield);CHKERRQ(ierr); 82df21e3a8SDave May } 8377048351SPatrick Sanan ierr = DMSwarmDataExBegin(de);CHKERRQ(ierr); 8477048351SPatrick Sanan ierr = DMSwarmDataExEnd(de);CHKERRQ(ierr); 8577048351SPatrick Sanan ierr = DMSwarmDataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr); 8677048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 8777048351SPatrick Sanan ierr = DMSwarmDataBucketSetSizes(swarm->db,npoints + n_points_recv,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr); 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 9177048351SPatrick Sanan ierr = DMSwarmDataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr); 92df21e3a8SDave May } 9377048351SPatrick Sanan ierr = DMSwarmDataExView(de);CHKERRQ(ierr); 9477048351SPatrick Sanan ierr = DMSwarmDataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr); 9577048351SPatrick Sanan ierr = DMSwarmDataExDestroy(de);CHKERRQ(ierr); 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; 10240c453e9SDave May PetscErrorCode ierr; 10377048351SPatrick Sanan DMSwarmDataEx de; 10440c453e9SDave May PetscInt r,p,npoints,*rankval,n_points_recv; 10540c453e9SDave May PetscMPIInt rank,_rank; 10640c453e9SDave May const PetscMPIInt *neighbourranks; 10740c453e9SDave May void *point_buffer,*recv_points; 10840c453e9SDave May size_t sizeof_dmswarm_point; 10940c453e9SDave May PetscInt nneighbors; 1107c6d1d28SDave May PetscMPIInt mynneigh,*myneigh; 11140c453e9SDave May 112521f74f9SMatthew G. Knepley PetscFunctionBegin; 11340c453e9SDave May ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 11477048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 11508056efcSDave May ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 11677048351SPatrick Sanan ierr = DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr); 11740c453e9SDave May ierr = DMGetNeighbors(dmcell,&nneighbors,&neighbourranks);CHKERRQ(ierr); 11877048351SPatrick Sanan ierr = DMSwarmDataExTopologyInitialize(de);CHKERRQ(ierr); 11940c453e9SDave May for (r=0; r<nneighbors; r++) { 12040c453e9SDave May _rank = neighbourranks[r]; 12140c453e9SDave May if ((_rank != rank) && (_rank > 0)) { 12277048351SPatrick Sanan ierr = DMSwarmDataExTopologyAddNeighbour(de,_rank);CHKERRQ(ierr); 12340c453e9SDave May } 12440c453e9SDave May } 12577048351SPatrick Sanan ierr = DMSwarmDataExTopologyFinalize(de);CHKERRQ(ierr); 12677048351SPatrick Sanan ierr = DMSwarmDataExTopologyGetNeighbours(de,&mynneigh,&myneigh);CHKERRQ(ierr); 12777048351SPatrick Sanan ierr = DMSwarmDataExInitializeSendCount(de);CHKERRQ(ierr); 12840c453e9SDave May for (p=0; p<npoints; p++) { 129f954cb40SDave May if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) { 1307c6d1d28SDave May for (r=0; r<mynneigh; r++) { 1317c6d1d28SDave May _rank = myneigh[r]; 13277048351SPatrick Sanan ierr = DMSwarmDataExAddToSendCount(de,_rank,1);CHKERRQ(ierr); 13340c453e9SDave May } 13440c453e9SDave May } 13540c453e9SDave May } 13677048351SPatrick Sanan ierr = DMSwarmDataExFinalizeSendCount(de);CHKERRQ(ierr); 13777048351SPatrick Sanan ierr = DMSwarmDataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr); 13877048351SPatrick Sanan ierr = DMSwarmDataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr); 13940c453e9SDave May for (p=0; p<npoints; p++) { 140f954cb40SDave May if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) { 1417c6d1d28SDave May for (r=0; r<mynneigh; r++) { 1427c6d1d28SDave May _rank = myneigh[r]; 14340c453e9SDave May /* copy point into buffer */ 14477048351SPatrick Sanan ierr = DMSwarmDataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr); 14577048351SPatrick Sanan /* insert point buffer into DMSwarmDataExchanger */ 14677048351SPatrick Sanan ierr = DMSwarmDataExPackData(de,_rank,1,point_buffer);CHKERRQ(ierr); 14740c453e9SDave May } 14840c453e9SDave May } 14940c453e9SDave May } 15077048351SPatrick Sanan ierr = DMSwarmDataExPackFinalize(de);CHKERRQ(ierr); 1517c6d1d28SDave May ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 15240c453e9SDave May if (remove_sent_points) { 15377048351SPatrick Sanan DMSwarmDataField PField; 1547c6d1d28SDave May 15577048351SPatrick Sanan ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr); 15677048351SPatrick Sanan ierr = DMSwarmDataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); 15740c453e9SDave May /* remove points which left processor */ 15877048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 15940c453e9SDave May for (p=0; p<npoints; p++) { 160f954cb40SDave May if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) { 16140c453e9SDave May /* kill point */ 16277048351SPatrick Sanan ierr = DMSwarmDataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr); 16377048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */ 16477048351SPatrick Sanan ierr = DMSwarmDataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */ 16540c453e9SDave May p--; /* check replacement point */ 16640c453e9SDave May } 16740c453e9SDave May } 16840c453e9SDave May } 16977048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,npoints_prior_migration,NULL,NULL);CHKERRQ(ierr); 17077048351SPatrick Sanan ierr = DMSwarmDataExBegin(de);CHKERRQ(ierr); 17177048351SPatrick Sanan ierr = DMSwarmDataExEnd(de);CHKERRQ(ierr); 17277048351SPatrick Sanan ierr = DMSwarmDataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr); 17377048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 17477048351SPatrick Sanan ierr = DMSwarmDataBucketSetSizes(swarm->db,npoints + n_points_recv,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr); 17540c453e9SDave May for (p=0; p<n_points_recv; p++) { 17640c453e9SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point ); 17740c453e9SDave May 17877048351SPatrick Sanan ierr = DMSwarmDataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr); 17940c453e9SDave May } 18077048351SPatrick Sanan ierr = DMSwarmDataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr); 18177048351SPatrick Sanan ierr = DMSwarmDataExDestroy(de);CHKERRQ(ierr); 18240c453e9SDave May PetscFunctionReturn(0); 18340c453e9SDave May } 184480eef7bSDave May 18508056efcSDave May PetscErrorCode DMSwarmMigrate_CellDMScatter(DM dm,PetscBool remove_sent_points) 186480eef7bSDave May { 187480eef7bSDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 188480eef7bSDave May PetscErrorCode ierr; 1896fbf25f8SDave May PetscInt p,npoints,npointsg=0,npoints2,npoints2g,*rankval,npoints_prior_migration; 190bbe8250bSMatthew G. Knepley PetscSF sfcell = NULL; 191dfc14de9SMatthew G. Knepley const PetscSFNode *LA_sfcell; 192480eef7bSDave May DM dmcell; 193480eef7bSDave May Vec pos; 19440c453e9SDave May PetscBool error_check = swarm->migrate_error_on_missing_point; 195e4fbd051SBarry Smith PetscMPIInt size,rank; 196480eef7bSDave May 197521f74f9SMatthew G. Knepley PetscFunctionBegin; 198480eef7bSDave May ierr = DMSwarmGetCellDM(dm,&dmcell);CHKERRQ(ierr); 199480eef7bSDave May if (!dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided"); 200480eef7bSDave May 201e4fbd051SBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr); 2027c6d1d28SDave May ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 2037c6d1d28SDave May 20443a82f2bSDave May #if 1 20543a82f2bSDave May { 20643a82f2bSDave May PetscInt *p_cellid; 20743a82f2bSDave May PetscInt npoints_curr,range = 0; 20843a82f2bSDave May PetscSFNode *sf_cells; 20943a82f2bSDave May 21043a82f2bSDave May 21177048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints_curr,NULL,NULL);CHKERRQ(ierr); 21243a82f2bSDave May ierr = PetscMalloc1(npoints_curr, &sf_cells);CHKERRQ(ierr); 21343a82f2bSDave May 21443a82f2bSDave May ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 21543a82f2bSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&p_cellid);CHKERRQ(ierr); 21643a82f2bSDave May for (p=0; p<npoints_curr; p++) { 21743a82f2bSDave May 21843a82f2bSDave May sf_cells[p].rank = 0; 21943a82f2bSDave May sf_cells[p].index = p_cellid[p]; 22043a82f2bSDave May if (p_cellid[p] > range) { 22143a82f2bSDave May range = p_cellid[p]; 22243a82f2bSDave May } 22343a82f2bSDave May } 22443a82f2bSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&p_cellid);CHKERRQ(ierr); 22543a82f2bSDave May ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 22643a82f2bSDave May 227c28a3c15SDave May /*ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm),&sfcell);CHKERRQ(ierr);*/ 228c28a3c15SDave May ierr = PetscSFCreate(PETSC_COMM_SELF,&sfcell);CHKERRQ(ierr); 22943a82f2bSDave May ierr = PetscSFSetGraph(sfcell, range, npoints_curr, NULL, PETSC_OWN_POINTER, sf_cells, PETSC_OWN_POINTER);CHKERRQ(ierr); 23043a82f2bSDave May } 23143a82f2bSDave May #endif 23243a82f2bSDave May 233fb1bcc12SMatthew G. Knepley ierr = DMSwarmCreateLocalVectorFromField(dm, DMSwarmPICField_coor, &pos);CHKERRQ(ierr); 234dfc14de9SMatthew G. Knepley ierr = DMLocatePoints(dmcell, pos, DM_POINTLOCATION_NONE, &sfcell);CHKERRQ(ierr); 235fb1bcc12SMatthew G. Knepley ierr = DMSwarmDestroyLocalVectorFromField(dm, DMSwarmPICField_coor, &pos);CHKERRQ(ierr); 236480eef7bSDave May 23740c453e9SDave May if (error_check) { 23840c453e9SDave May ierr = DMSwarmGetSize(dm,&npointsg);CHKERRQ(ierr); 23940c453e9SDave May } 24077048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 24108056efcSDave May ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 242dfc14de9SMatthew G. Knepley ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr); 243480eef7bSDave May for (p=0; p<npoints; p++) { 244dfc14de9SMatthew G. Knepley rankval[p] = LA_sfcell[p].index; 245480eef7bSDave May } 24608056efcSDave May ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 247dfc14de9SMatthew G. Knepley ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr); 248480eef7bSDave May 249e4fbd051SBarry Smith if (size > 1) { 250889dbfe5SDave May ierr = DMSwarmMigrate_DMNeighborScatter(dm,dmcell,remove_sent_points,&npoints_prior_migration);CHKERRQ(ierr); 2516fbf25f8SDave May } else { 25277048351SPatrick Sanan DMSwarmDataField PField; 2530ed23c7fSDave May PetscInt npoints_curr; 2540ed23c7fSDave May 2550ed23c7fSDave May /* remove points which the domain */ 25677048351SPatrick Sanan ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr); 25777048351SPatrick Sanan ierr = DMSwarmDataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); 2580ed23c7fSDave May 25977048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints_curr,NULL,NULL);CHKERRQ(ierr); 2600ed23c7fSDave May for (p=0; p<npoints_curr; p++) { 2610ed23c7fSDave May if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) { 2620ed23c7fSDave May /* kill point */ 26377048351SPatrick Sanan ierr = DMSwarmDataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr); 26477048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints_curr,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */ 26577048351SPatrick Sanan ierr = DMSwarmDataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */ 2660ed23c7fSDave May p--; /* check replacement point */ 2670ed23c7fSDave May } 2680ed23c7fSDave May } 2696fbf25f8SDave May ierr = DMSwarmGetSize(dm,&npoints_prior_migration);CHKERRQ(ierr); 2700ed23c7fSDave May 2716fbf25f8SDave May } 272480eef7bSDave May 273*2d4ee042Sprj- /* locate points newly received */ 27477048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); 275009b43efSDave May 2767c6d1d28SDave May #if 0 277*2d4ee042Sprj- { /* safe alternative - however this performs two point locations on: (i) the initial points set and; (ii) the (initial + received) point set */ 2787c6d1d28SDave May PetscScalar *LA_coor; 2797c6d1d28SDave May PetscInt bs; 28077048351SPatrick Sanan DMSwarmDataField PField; 2817c6d1d28SDave May 2827c6d1d28SDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr); 2837c6d1d28SDave May ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*npoints2,(const PetscScalar*)LA_coor,&pos);CHKERRQ(ierr); 284dfc14de9SMatthew G. Knepley ierr = DMLocatePoints(dmcell,pos,DM_POINTLOCATION_NONE,&sfcell);CHKERRQ(ierr); 2857c6d1d28SDave May 2867c6d1d28SDave May ierr = VecDestroy(&pos);CHKERRQ(ierr); 2877c6d1d28SDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr); 2887c6d1d28SDave May 289dfc14de9SMatthew G. Knepley ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr); 2907c6d1d28SDave May ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 2917c6d1d28SDave May for (p=0; p<npoints2; p++) { 292dfc14de9SMatthew G. Knepley rankval[p] = LA_sfcell[p].index; 2937c6d1d28SDave May } 294dfc14de9SMatthew G. Knepley ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr); 29508056efcSDave May ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 29640c453e9SDave May 2977c6d1d28SDave May /* remove points which left processor */ 29877048351SPatrick Sanan ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr); 29977048351SPatrick Sanan ierr = DMSwarmDataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); 3007c6d1d28SDave May 30177048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); 3027c6d1d28SDave May for (p=0; p<npoints2; p++) { 303f954cb40SDave May if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) { 3047c6d1d28SDave May /* kill point */ 30577048351SPatrick Sanan ierr = DMSwarmDataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr); 30677048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */ 30777048351SPatrick Sanan ierr = DMSwarmDataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */ 3087c6d1d28SDave May p--; /* check replacement point */ 3097c6d1d28SDave May } 3107c6d1d28SDave May } 31140c453e9SDave May } 312009b43efSDave May #endif 313009b43efSDave May 314*2d4ee042Sprj- { /* this performs two point locations: (i) on the initial points set prior to communication; and (ii) on the new (received) points */ 315009b43efSDave May PetscScalar *LA_coor; 316009b43efSDave May PetscInt npoints_from_neighbours,bs; 31777048351SPatrick Sanan DMSwarmDataField PField; 318009b43efSDave May 319009b43efSDave May npoints_from_neighbours = npoints2 - npoints_prior_migration; 320009b43efSDave May 321009b43efSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr); 322009b43efSDave May ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*npoints_from_neighbours,(const PetscScalar*)&LA_coor[bs*npoints_prior_migration],&pos);CHKERRQ(ierr); 323009b43efSDave May 324009b43efSDave May ierr = DMLocatePoints(dmcell,pos,DM_POINTLOCATION_NONE,&sfcell);CHKERRQ(ierr); 325009b43efSDave May 326009b43efSDave May ierr = VecDestroy(&pos);CHKERRQ(ierr); 327009b43efSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr); 328009b43efSDave May 329009b43efSDave May ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr); 330009b43efSDave May ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 331009b43efSDave May for (p=0; p<npoints_from_neighbours; p++) { 332009b43efSDave May rankval[npoints_prior_migration + p] = LA_sfcell[p].index; 333009b43efSDave May } 334009b43efSDave May ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 335009b43efSDave May ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr); 336009b43efSDave May 337009b43efSDave May /* remove points which left processor */ 33877048351SPatrick Sanan ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr); 33977048351SPatrick Sanan ierr = DMSwarmDataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); 340009b43efSDave May 34177048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); 342009b43efSDave May for (p=npoints_prior_migration; p<npoints2; p++) { 343009b43efSDave May if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) { 344009b43efSDave May /* kill point */ 34577048351SPatrick Sanan ierr = DMSwarmDataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr); 34677048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */ 34777048351SPatrick Sanan ierr = DMSwarmDataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */ 348009b43efSDave May p--; /* check replacement point */ 349009b43efSDave May } 350009b43efSDave May } 351009b43efSDave May } 352009b43efSDave May 3533151f1d5SDave May { 3543151f1d5SDave May PetscInt *p_cellid; 3553151f1d5SDave May 35677048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); 3573151f1d5SDave May ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 3583151f1d5SDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&p_cellid);CHKERRQ(ierr); 3593151f1d5SDave May for (p=0; p<npoints2; p++) { 3603151f1d5SDave May p_cellid[p] = rankval[p]; 3613151f1d5SDave May } 3623151f1d5SDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&p_cellid);CHKERRQ(ierr); 3633151f1d5SDave May ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 3643151f1d5SDave May } 3653151f1d5SDave May 36640c453e9SDave May /* check for error on removed points */ 36740c453e9SDave May if (error_check) { 36840c453e9SDave May ierr = DMSwarmGetSize(dm,&npoints2g);CHKERRQ(ierr); 36940c453e9SDave 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); 37040c453e9SDave May } 371480eef7bSDave May PetscFunctionReturn(0); 372480eef7bSDave May } 373480eef7bSDave May 37408056efcSDave May PetscErrorCode DMSwarmMigrate_CellDMExact(DM dm,PetscBool remove_sent_points) 37508056efcSDave May { 376521f74f9SMatthew G. Knepley PetscFunctionBegin; 37708056efcSDave May PetscFunctionReturn(0); 37808056efcSDave May } 37908056efcSDave May 380480eef7bSDave May /* 381480eef7bSDave May Redundant as this assumes points can only be sent to a single rank 382480eef7bSDave May */ 3832712d1f2SDave May PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt *globalsize) 3842712d1f2SDave May { 3852712d1f2SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 3862712d1f2SDave May PetscErrorCode ierr; 38777048351SPatrick Sanan DMSwarmDataEx de; 3882712d1f2SDave May PetscInt p,npoints,*rankval,n_points_recv; 3892712d1f2SDave May PetscMPIInt rank,nrank,negrank; 3902712d1f2SDave May void *point_buffer,*recv_points; 3912712d1f2SDave May size_t sizeof_dmswarm_point; 3922712d1f2SDave May 393521f74f9SMatthew G. Knepley PetscFunctionBegin; 3942712d1f2SDave May ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 39577048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 3962712d1f2SDave May *globalsize = npoints; 39708056efcSDave May ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 39877048351SPatrick Sanan ierr = DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr); 39977048351SPatrick Sanan ierr = DMSwarmDataExTopologyInitialize(de);CHKERRQ(ierr); 4002712d1f2SDave May for (p=0; p<npoints; p++) { 4012712d1f2SDave May negrank = rankval[p]; 4022712d1f2SDave May if (negrank < 0) { 4032712d1f2SDave May nrank = -negrank - 1; 40477048351SPatrick Sanan ierr = DMSwarmDataExTopologyAddNeighbour(de,nrank);CHKERRQ(ierr); 4052712d1f2SDave May } 4062712d1f2SDave May } 40777048351SPatrick Sanan ierr = DMSwarmDataExTopologyFinalize(de);CHKERRQ(ierr); 40877048351SPatrick Sanan ierr = DMSwarmDataExInitializeSendCount(de);CHKERRQ(ierr); 4092712d1f2SDave May for (p=0; p<npoints; p++) { 4102712d1f2SDave May negrank = rankval[p]; 4112712d1f2SDave May if (negrank < 0) { 4122712d1f2SDave May nrank = -negrank - 1; 41377048351SPatrick Sanan ierr = DMSwarmDataExAddToSendCount(de,nrank,1);CHKERRQ(ierr); 4142712d1f2SDave May } 4152712d1f2SDave May } 41677048351SPatrick Sanan ierr = DMSwarmDataExFinalizeSendCount(de);CHKERRQ(ierr); 41777048351SPatrick Sanan ierr = DMSwarmDataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr); 41877048351SPatrick Sanan ierr = DMSwarmDataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr); 4192712d1f2SDave May for (p=0; p<npoints; p++) { 4202712d1f2SDave May negrank = rankval[p]; 4212712d1f2SDave May if (negrank < 0) { 4222712d1f2SDave May nrank = -negrank - 1; 4232712d1f2SDave May rankval[p] = nrank; 4242712d1f2SDave May /* copy point into buffer */ 42577048351SPatrick Sanan ierr = DMSwarmDataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr); 42677048351SPatrick Sanan /* insert point buffer into DMSwarmDataExchanger */ 42777048351SPatrick Sanan ierr = DMSwarmDataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr); 4282712d1f2SDave May rankval[p] = negrank; 4292712d1f2SDave May } 4302712d1f2SDave May } 43177048351SPatrick Sanan ierr = DMSwarmDataExPackFinalize(de);CHKERRQ(ierr); 43208056efcSDave May ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 43377048351SPatrick Sanan ierr = DMSwarmDataExBegin(de);CHKERRQ(ierr); 43477048351SPatrick Sanan ierr = DMSwarmDataExEnd(de);CHKERRQ(ierr); 43577048351SPatrick Sanan ierr = DMSwarmDataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr); 43677048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 43777048351SPatrick Sanan ierr = DMSwarmDataBucketSetSizes(swarm->db,npoints + n_points_recv,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr); 4382712d1f2SDave May for (p=0; p<n_points_recv; p++) { 4392712d1f2SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point ); 4402712d1f2SDave May 44177048351SPatrick Sanan ierr = DMSwarmDataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr); 4422712d1f2SDave May } 44377048351SPatrick Sanan ierr = DMSwarmDataExView(de);CHKERRQ(ierr); 44477048351SPatrick Sanan ierr = DMSwarmDataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr); 44577048351SPatrick Sanan ierr = DMSwarmDataExDestroy(de);CHKERRQ(ierr); 4462712d1f2SDave May PetscFunctionReturn(0); 4472712d1f2SDave May } 448b16650c8SDave May 449b16650c8SDave May typedef struct { 450b16650c8SDave May PetscMPIInt owner_rank; 451b16650c8SDave May PetscReal min[3],max[3]; 452b16650c8SDave May } CollectBBox; 453b16650c8SDave May 454fe39f135SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM dm,PetscInt *globalsize) 455b16650c8SDave May { 456b16650c8SDave May DM_Swarm * swarm = (DM_Swarm*)dm->data; 457b16650c8SDave May PetscErrorCode ierr; 45877048351SPatrick Sanan DMSwarmDataEx de; 459b16650c8SDave May PetscInt p,pk,npoints,*rankval,n_points_recv,n_bbox_recv,dim,neighbour_cells; 460b16650c8SDave May PetscMPIInt rank,nrank; 461b16650c8SDave May void *point_buffer,*recv_points; 462b16650c8SDave May size_t sizeof_dmswarm_point,sizeof_bbox_ctx; 463b16650c8SDave May PetscBool isdmda; 464b16650c8SDave May CollectBBox *bbox,*recv_bbox; 465b16650c8SDave May const PetscMPIInt *dmneighborranks; 466b16650c8SDave May DM dmcell; 467b16650c8SDave May 468521f74f9SMatthew G. Knepley PetscFunctionBegin; 469b16650c8SDave May ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 470b16650c8SDave May 471fe39f135SDave May ierr = DMSwarmGetCellDM(dm,&dmcell);CHKERRQ(ierr); 472fe39f135SDave May if (!dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided"); 473b16650c8SDave May isdmda = PETSC_FALSE; 474b16650c8SDave May PetscObjectTypeCompare((PetscObject)dmcell,DMDA,&isdmda); 475b16650c8SDave May if (!isdmda) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only DMDA support for CollectBoundingBox"); 476b16650c8SDave May 4778dbd68bcSDave May ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 478b16650c8SDave May sizeof_bbox_ctx = sizeof(CollectBBox); 479b16650c8SDave May PetscMalloc1(1,&bbox); 480b16650c8SDave May bbox->owner_rank = rank; 481b16650c8SDave May 482b16650c8SDave May /* compute the bounding box based on the overlapping / stenctil size */ 483b16650c8SDave May { 484b16650c8SDave May Vec lcoor; 485b16650c8SDave May 486b16650c8SDave May ierr = DMGetCoordinatesLocal(dmcell,&lcoor);CHKERRQ(ierr); 487fe39f135SDave May if (dim >= 1) { 488b16650c8SDave May ierr = VecStrideMin(lcoor,0,NULL,&bbox->min[0]);CHKERRQ(ierr); 489b16650c8SDave May ierr = VecStrideMax(lcoor,0,NULL,&bbox->max[0]);CHKERRQ(ierr); 490fe39f135SDave May } 491fe39f135SDave May if (dim >= 2) { 492fe39f135SDave May ierr = VecStrideMin(lcoor,1,NULL,&bbox->min[1]);CHKERRQ(ierr); 493b16650c8SDave May ierr = VecStrideMax(lcoor,1,NULL,&bbox->max[1]);CHKERRQ(ierr); 494b16650c8SDave May } 495fe39f135SDave May if (dim == 3) { 496fe39f135SDave May ierr = VecStrideMin(lcoor,2,NULL,&bbox->min[2]);CHKERRQ(ierr); 497fe39f135SDave May ierr = VecStrideMax(lcoor,2,NULL,&bbox->max[2]);CHKERRQ(ierr); 498fe39f135SDave May } 499fe39f135SDave May } 50077048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 501b16650c8SDave May *globalsize = npoints; 50208056efcSDave May ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 50377048351SPatrick Sanan ierr = DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr); 504b16650c8SDave May /* use DMDA neighbours */ 505b16650c8SDave May ierr = DMDAGetNeighbors(dmcell,&dmneighborranks);CHKERRQ(ierr); 5068dbd68bcSDave May if (dim == 1) { 5078dbd68bcSDave May neighbour_cells = 3; 5088dbd68bcSDave May } else if (dim == 2) { 5098dbd68bcSDave May neighbour_cells = 9; 5108dbd68bcSDave May } else { 5118dbd68bcSDave May neighbour_cells = 27; 5128dbd68bcSDave May } 51377048351SPatrick Sanan ierr = DMSwarmDataExTopologyInitialize(de);CHKERRQ(ierr); 514b16650c8SDave May for (p=0; p<neighbour_cells; p++) { 515b16650c8SDave May if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) { 51677048351SPatrick Sanan ierr = DMSwarmDataExTopologyAddNeighbour(de,dmneighborranks[p]);CHKERRQ(ierr); 517b16650c8SDave May } 518b16650c8SDave May } 51977048351SPatrick Sanan ierr = DMSwarmDataExTopologyFinalize(de);CHKERRQ(ierr); 52077048351SPatrick Sanan ierr = DMSwarmDataExInitializeSendCount(de);CHKERRQ(ierr); 521b16650c8SDave May for (p=0; p<neighbour_cells; p++) { 522b16650c8SDave May if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) { 52377048351SPatrick Sanan ierr = DMSwarmDataExAddToSendCount(de,dmneighborranks[p],1);CHKERRQ(ierr); 524b16650c8SDave May } 525b16650c8SDave May } 52677048351SPatrick Sanan ierr = DMSwarmDataExFinalizeSendCount(de);CHKERRQ(ierr); 527b16650c8SDave May /* send bounding boxes */ 52877048351SPatrick Sanan ierr = DMSwarmDataExPackInitialize(de,sizeof_bbox_ctx);CHKERRQ(ierr); 529b16650c8SDave May for (p=0; p<neighbour_cells; p++) { 530b16650c8SDave May nrank = dmneighborranks[p]; 531b16650c8SDave May if ( (nrank >= 0) && (nrank != rank) ) { 53277048351SPatrick Sanan /* insert bbox buffer into DMSwarmDataExchanger */ 53377048351SPatrick Sanan ierr = DMSwarmDataExPackData(de,nrank,1,bbox);CHKERRQ(ierr); 534b16650c8SDave May } 535b16650c8SDave May } 53677048351SPatrick Sanan ierr = DMSwarmDataExPackFinalize(de);CHKERRQ(ierr); 537b16650c8SDave May /* recv bounding boxes */ 53877048351SPatrick Sanan ierr = DMSwarmDataExBegin(de);CHKERRQ(ierr); 53977048351SPatrick Sanan ierr = DMSwarmDataExEnd(de);CHKERRQ(ierr); 54077048351SPatrick Sanan ierr = DMSwarmDataExGetRecvData(de,&n_bbox_recv,(void**)&recv_bbox);CHKERRQ(ierr); 541298827fbSBarry Smith /* Wrong, should not be using PETSC_COMM_WORLD */ 542b16650c8SDave May for (p=0; p<n_bbox_recv; p++) { 543298827fbSBarry Smith ierr = 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, 544298827fbSBarry Smith (double)recv_bbox[p].min[0],(double)recv_bbox[p].max[0],(double)recv_bbox[p].min[1],(double)recv_bbox[p].max[1]);CHKERRQ(ierr); 545b16650c8SDave May } 546298827fbSBarry Smith ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,stdout);CHKERRQ(ierr); 547b16650c8SDave May /* of course this is stupid as this "generic" function should have a better way to know what the coordinates are called */ 54877048351SPatrick Sanan ierr = DMSwarmDataExInitializeSendCount(de);CHKERRQ(ierr); 549b16650c8SDave May for (pk=0; pk<n_bbox_recv; pk++) { 550b16650c8SDave May PetscReal *array_x,*array_y; 551b16650c8SDave May 552b16650c8SDave May ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 553b16650c8SDave May ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 554b16650c8SDave May for (p=0; p<npoints; p++) { 555b16650c8SDave May if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) { 556b16650c8SDave May if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) { 55777048351SPatrick Sanan ierr = DMSwarmDataExAddToSendCount(de,recv_bbox[pk].owner_rank,1);CHKERRQ(ierr); 558b16650c8SDave May } 559b16650c8SDave May } 560b16650c8SDave May } 561b16650c8SDave May ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 562b16650c8SDave May ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 563b16650c8SDave May } 56477048351SPatrick Sanan ierr = DMSwarmDataExFinalizeSendCount(de);CHKERRQ(ierr); 56577048351SPatrick Sanan ierr = DMSwarmDataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr); 56677048351SPatrick Sanan ierr = DMSwarmDataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr); 567b16650c8SDave May for (pk=0; pk<n_bbox_recv; pk++) { 568b16650c8SDave May PetscReal *array_x,*array_y; 569b16650c8SDave May 570b16650c8SDave May ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 571b16650c8SDave May ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 572b16650c8SDave May for (p=0; p<npoints; p++) { 573b16650c8SDave May if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) { 574b16650c8SDave May if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) { 575521f74f9SMatthew G. Knepley /* copy point into buffer */ 57677048351SPatrick Sanan ierr = DMSwarmDataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr); 57777048351SPatrick Sanan /* insert point buffer into DMSwarmDataExchanger */ 57877048351SPatrick Sanan ierr = DMSwarmDataExPackData(de,recv_bbox[pk].owner_rank,1,point_buffer);CHKERRQ(ierr); 579b16650c8SDave May } 580b16650c8SDave May } 581b16650c8SDave May } 582b16650c8SDave May ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 583b16650c8SDave May ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 584b16650c8SDave May } 58577048351SPatrick Sanan ierr = DMSwarmDataExPackFinalize(de);CHKERRQ(ierr); 58608056efcSDave May ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 58777048351SPatrick Sanan ierr = DMSwarmDataExBegin(de);CHKERRQ(ierr); 58877048351SPatrick Sanan ierr = DMSwarmDataExEnd(de);CHKERRQ(ierr); 58977048351SPatrick Sanan ierr = DMSwarmDataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr); 59077048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 59177048351SPatrick Sanan ierr = DMSwarmDataBucketSetSizes(swarm->db,npoints + n_points_recv,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr); 592b16650c8SDave May for (p=0; p<n_points_recv; p++) { 593b16650c8SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point ); 594b16650c8SDave May 59577048351SPatrick Sanan ierr = DMSwarmDataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr); 596b16650c8SDave May } 59777048351SPatrick Sanan ierr = DMSwarmDataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr); 598b16650c8SDave May PetscFree(bbox); 59977048351SPatrick Sanan ierr = DMSwarmDataExView(de);CHKERRQ(ierr); 60077048351SPatrick Sanan ierr = DMSwarmDataExDestroy(de);CHKERRQ(ierr); 601b16650c8SDave May PetscFunctionReturn(0); 602b16650c8SDave May } 603a9fd7477SDave May 604a9fd7477SDave May 605a9fd7477SDave May /* General collection when no order, or neighbour information is provided */ 606a9fd7477SDave May /* 607a9fd7477SDave May User provides context and collect() method 608a9fd7477SDave May Broadcast user context 609a9fd7477SDave May 610a9fd7477SDave May for each context / rank { 611a9fd7477SDave May collect(swarm,context,n,list) 612a9fd7477SDave May } 613a9fd7477SDave May */ 614a9fd7477SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM dm,PetscErrorCode (*collect)(DM,void*,PetscInt*,PetscInt**),size_t ctx_size,void *ctx,PetscInt *globalsize) 615a9fd7477SDave May { 616a9fd7477SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 617a9fd7477SDave May PetscErrorCode ierr; 61877048351SPatrick Sanan DMSwarmDataEx de; 619a9fd7477SDave May PetscInt p,r,npoints,n_points_recv; 620e4fbd051SBarry Smith PetscMPIInt size,rank; 621a9fd7477SDave May void *point_buffer,*recv_points; 622a9fd7477SDave May void *ctxlist; 623a9fd7477SDave May PetscInt *n2collect,**collectlist; 624a9fd7477SDave May size_t sizeof_dmswarm_point; 625a9fd7477SDave May 626521f74f9SMatthew G. Knepley PetscFunctionBegin; 627e4fbd051SBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr); 628a9fd7477SDave May ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 62977048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 630a9fd7477SDave May *globalsize = npoints; 631a9fd7477SDave May /* Broadcast user context */ 632e4fbd051SBarry Smith PetscMalloc(ctx_size*size,&ctxlist); 633a9fd7477SDave May ierr = MPI_Allgather(ctx,ctx_size,MPI_CHAR,ctxlist,ctx_size,MPI_CHAR,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 634e4fbd051SBarry Smith ierr = PetscMalloc1(size,&n2collect);CHKERRQ(ierr); 635e4fbd051SBarry Smith ierr = PetscMalloc1(size,&collectlist);CHKERRQ(ierr); 636e4fbd051SBarry Smith for (r=0; r<size; r++) { 637a9fd7477SDave May PetscInt _n2collect; 638a9fd7477SDave May PetscInt *_collectlist; 639a9fd7477SDave May void *_ctx_r; 640a9fd7477SDave May 641a9fd7477SDave May _n2collect = 0; 642a9fd7477SDave May _collectlist = NULL; 643a9fd7477SDave May if (r != rank) { /* don't collect data from yourself */ 644a9fd7477SDave May _ctx_r = (void*)( (char*)ctxlist + r * ctx_size ); 645a9fd7477SDave May ierr = collect(dm,_ctx_r,&_n2collect,&_collectlist);CHKERRQ(ierr); 646a9fd7477SDave May } 647a9fd7477SDave May n2collect[r] = _n2collect; 648a9fd7477SDave May collectlist[r] = _collectlist; 649a9fd7477SDave May } 65077048351SPatrick Sanan ierr = DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr); 651a9fd7477SDave May /* Define topology */ 65277048351SPatrick Sanan ierr = DMSwarmDataExTopologyInitialize(de);CHKERRQ(ierr); 653e4fbd051SBarry Smith for (r=0; r<size; r++) { 654a9fd7477SDave May if (n2collect[r] > 0) { 65577048351SPatrick Sanan ierr = DMSwarmDataExTopologyAddNeighbour(de,(PetscMPIInt)r);CHKERRQ(ierr); 656a9fd7477SDave May } 657a9fd7477SDave May } 65877048351SPatrick Sanan ierr = DMSwarmDataExTopologyFinalize(de);CHKERRQ(ierr); 659a9fd7477SDave May /* Define send counts */ 66077048351SPatrick Sanan ierr = DMSwarmDataExInitializeSendCount(de);CHKERRQ(ierr); 661e4fbd051SBarry Smith for (r=0; r<size; r++) { 662a9fd7477SDave May if (n2collect[r] > 0) { 66377048351SPatrick Sanan ierr = DMSwarmDataExAddToSendCount(de,r,n2collect[r]);CHKERRQ(ierr); 664a9fd7477SDave May } 665a9fd7477SDave May } 66677048351SPatrick Sanan ierr = DMSwarmDataExFinalizeSendCount(de);CHKERRQ(ierr); 667a9fd7477SDave May /* Pack data */ 66877048351SPatrick Sanan ierr = DMSwarmDataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr); 66977048351SPatrick Sanan ierr = DMSwarmDataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr); 670e4fbd051SBarry Smith for (r=0; r<size; r++) { 671a9fd7477SDave May for (p=0; p<n2collect[r]; p++) { 67277048351SPatrick Sanan ierr = DMSwarmDataBucketFillPackedArray(swarm->db,collectlist[r][p],point_buffer);CHKERRQ(ierr); 673a9fd7477SDave May /* insert point buffer into the data exchanger */ 67477048351SPatrick Sanan ierr = DMSwarmDataExPackData(de,r,1,point_buffer);CHKERRQ(ierr); 675a9fd7477SDave May } 676a9fd7477SDave May } 67777048351SPatrick Sanan ierr = DMSwarmDataExPackFinalize(de);CHKERRQ(ierr); 678a9fd7477SDave May /* Scatter */ 67977048351SPatrick Sanan ierr = DMSwarmDataExBegin(de);CHKERRQ(ierr); 68077048351SPatrick Sanan ierr = DMSwarmDataExEnd(de);CHKERRQ(ierr); 681a9fd7477SDave May /* Collect data in DMSwarm container */ 68277048351SPatrick Sanan ierr = DMSwarmDataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr); 68377048351SPatrick Sanan ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 68477048351SPatrick Sanan ierr = DMSwarmDataBucketSetSizes(swarm->db,npoints + n_points_recv,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr); 685a9fd7477SDave May for (p=0; p<n_points_recv; p++) { 686a9fd7477SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point ); 687a9fd7477SDave May 68877048351SPatrick Sanan ierr = DMSwarmDataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr); 689a9fd7477SDave May } 690a9fd7477SDave May /* Release memory */ 691e4fbd051SBarry Smith for (r=0; r<size; r++) { 692a9fd7477SDave May if (collectlist[r]) PetscFree(collectlist[r]); 693a9fd7477SDave May } 694521f74f9SMatthew G. Knepley ierr = PetscFree(collectlist);CHKERRQ(ierr); 695521f74f9SMatthew G. Knepley ierr = PetscFree(n2collect);CHKERRQ(ierr); 696521f74f9SMatthew G. Knepley ierr = PetscFree(ctxlist);CHKERRQ(ierr); 69777048351SPatrick Sanan ierr = DMSwarmDataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr); 69877048351SPatrick Sanan ierr = DMSwarmDataExView(de);CHKERRQ(ierr); 69977048351SPatrick Sanan ierr = DMSwarmDataExDestroy(de);CHKERRQ(ierr); 700a9fd7477SDave May PetscFunctionReturn(0); 701a9fd7477SDave May } 702a9fd7477SDave May 703