1*df21e3a8SDave May 2*df21e3a8SDave May #include <petsc/private/dmswarmimpl.h> /*I "petscdmswarm.h" I*/ 3*df21e3a8SDave May #include "data_bucket.h" 4*df21e3a8SDave May #include "data_ex.h" 5*df21e3a8SDave May 6*df21e3a8SDave May 7*df21e3a8SDave May #undef __FUNCT__ 8*df21e3a8SDave May #define __FUNCT__ "DMSwarmMigrate_Push_Basic" 9*df21e3a8SDave May PetscErrorCode DMSwarmMigrate_Push_Basic(DM dm,PetscBool remove_sent_points) 10*df21e3a8SDave May { 11*df21e3a8SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 12*df21e3a8SDave May PetscErrorCode ierr; 13*df21e3a8SDave May DataEx de; 14*df21e3a8SDave May PetscInt p,npoints,*rankval,n_points_recv; 15*df21e3a8SDave May PetscMPIInt rank,nrank; 16*df21e3a8SDave May void *point_buffer,*recv_points; 17*df21e3a8SDave May size_t sizeof_dmswarm_point; 18*df21e3a8SDave May 19*df21e3a8SDave May ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 20*df21e3a8SDave May 21*df21e3a8SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 22*df21e3a8SDave May ierr = DMSwarmGetField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 23*df21e3a8SDave May 24*df21e3a8SDave May de = DataExCreate(PetscObjectComm((PetscObject)dm),0);CHKERRQ(ierr); 25*df21e3a8SDave May 26*df21e3a8SDave May ierr = DataExTopologyInitialize(de);CHKERRQ(ierr); 27*df21e3a8SDave May for (p=0; p<npoints; p++) { 28*df21e3a8SDave May nrank = rankval[p]; 29*df21e3a8SDave May if (nrank != rank) { 30*df21e3a8SDave May ierr = DataExTopologyAddNeighbour(de,nrank);CHKERRQ(ierr); 31*df21e3a8SDave May } 32*df21e3a8SDave May } 33*df21e3a8SDave May ierr = DataExTopologyFinalize(de);CHKERRQ(ierr); 34*df21e3a8SDave May 35*df21e3a8SDave May ierr = DataExInitializeSendCount(de);CHKERRQ(ierr); 36*df21e3a8SDave May for (p=0; p<npoints; p++) { 37*df21e3a8SDave May nrank = rankval[p]; 38*df21e3a8SDave May if (nrank != rank) { 39*df21e3a8SDave May ierr = DataExAddToSendCount(de,nrank,1);CHKERRQ(ierr); 40*df21e3a8SDave May } 41*df21e3a8SDave May } 42*df21e3a8SDave May ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr); 43*df21e3a8SDave May 44*df21e3a8SDave May ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr); 45*df21e3a8SDave May ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr); 46*df21e3a8SDave May for (p=0; p<npoints; p++) { 47*df21e3a8SDave May nrank = rankval[p]; 48*df21e3a8SDave May if (nrank != rank) { 49*df21e3a8SDave May /* copy point into buffer */ 50*df21e3a8SDave May ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr); 51*df21e3a8SDave May /* insert point buffer into DataExchanger */ 52*df21e3a8SDave May ierr = DataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr); 53*df21e3a8SDave May } 54*df21e3a8SDave May } 55*df21e3a8SDave May ierr = DataExPackFinalize(de);CHKERRQ(ierr); 56*df21e3a8SDave May 57*df21e3a8SDave May 58*df21e3a8SDave May if (remove_sent_points) { 59*df21e3a8SDave May /* remove points which left processor */ 60*df21e3a8SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 61*df21e3a8SDave May for (p=0; p<npoints; p++) { 62*df21e3a8SDave May nrank = rankval[p]; 63*df21e3a8SDave May if (nrank != rank) { 64*df21e3a8SDave May /* kill point */ 65*df21e3a8SDave May ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr); 66*df21e3a8SDave May DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */ 67*df21e3a8SDave May p--; /* check replacement point */ 68*df21e3a8SDave May } 69*df21e3a8SDave May } 70*df21e3a8SDave May } 71*df21e3a8SDave May ierr = DMSwarmRestoreField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 72*df21e3a8SDave May 73*df21e3a8SDave May ierr = DataExBegin(de);CHKERRQ(ierr); 74*df21e3a8SDave May ierr = DataExEnd(de);CHKERRQ(ierr); 75*df21e3a8SDave May 76*df21e3a8SDave May ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr); 77*df21e3a8SDave May 78*df21e3a8SDave May ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 79*df21e3a8SDave May ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr); 80*df21e3a8SDave May for (p=0; p<n_points_recv; p++) { 81*df21e3a8SDave May void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point ); 82*df21e3a8SDave May 83*df21e3a8SDave May ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr); 84*df21e3a8SDave May } 85*df21e3a8SDave May ierr = DataExView(de);CHKERRQ(ierr); 86*df21e3a8SDave May ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr); 87*df21e3a8SDave May ierr = DataExDestroy(de);CHKERRQ(ierr); 88*df21e3a8SDave May 89*df21e3a8SDave May PetscFunctionReturn(0); 90*df21e3a8SDave May } 91