| 1 | /**
|
|---|
| 2 | * @file comm_mpi_particle.cpp
|
|---|
| 3 | * @author Julian Iseringhausen <isering@ins.uni-bonn.de>
|
|---|
| 4 | * @date Tue May 8 15:27:06 2012
|
|---|
| 5 | *
|
|---|
| 6 | * @brief Class for MPI-based particle-related communication.
|
|---|
| 7 | *
|
|---|
| 8 | */
|
|---|
| 9 |
|
|---|
| 10 | #ifdef HAVE_CONFIG_H
|
|---|
| 11 | #include <config.h>
|
|---|
| 12 | #endif
|
|---|
| 13 |
|
|---|
| 14 | #ifdef HAVE_MPI
|
|---|
| 15 |
|
|---|
| 16 | #include <mpi.h>
|
|---|
| 17 | #ifdef HAVE_MARMOT
|
|---|
| 18 | #include <enhancempicalls.h>
|
|---|
| 19 | #include <sourceinfompicalls.h>
|
|---|
| 20 | #endif
|
|---|
| 21 |
|
|---|
| 22 | #include "units/particle/comm_mpi_particle.hpp"
|
|---|
| 23 | #include "units/particle/linked_cell_list.hpp"
|
|---|
| 24 | #include "units/particle/particle.hpp"
|
|---|
| 25 |
|
|---|
| 26 | using namespace VMG;
|
|---|
| 27 |
|
|---|
| 28 | void Particle::CommMPI::CommParticles(const Grid& grid, std::list<Particle>& particles)
|
|---|
| 29 | {
|
|---|
| 30 | Factory& factory = MG::GetFactory();
|
|---|
| 31 |
|
|---|
| 32 | const vmg_int& num_particles_local = factory.GetObjectStorageVal<vmg_int>("PARTICLE_NUM_LOCAL");
|
|---|
| 33 | vmg_float* x = factory.GetObjectStorageArray<vmg_float>("PARTICLE_POS_ARRAY");
|
|---|
| 34 | vmg_float* q = factory.GetObjectStorageArray<vmg_float>("PARTICLE_CHARGE_ARRAY");
|
|---|
| 35 |
|
|---|
| 36 | int rank, size;
|
|---|
| 37 | MPI_Comm_rank(comm_global, &rank);
|
|---|
| 38 | MPI_Comm_size(comm_global, &size);
|
|---|
| 39 |
|
|---|
| 40 | #ifndef VMG_ONE_SIDED
|
|---|
| 41 | vmg_int* receiver;
|
|---|
| 42 | if (!factory.TestObject("PARTICLE_RECEIVER_ARRAY"))
|
|---|
| 43 | new ObjectStorageArray<vmg_int>("PARTICLE_RECEIVER_ARRAY", size);
|
|---|
| 44 | receiver = factory.GetObjectStorageArray<vmg_int>("PARTICLE_RECEIVER_ARRAY");
|
|---|
| 45 | #endif
|
|---|
| 46 |
|
|---|
| 47 | Index index;
|
|---|
| 48 | std::vector<int> global_extent(6*size);
|
|---|
| 49 | std::vector<int> send_sizes(size);
|
|---|
| 50 | std::vector<int> recv_sizes(size);
|
|---|
| 51 | std::vector<Index> begin_remote(size);
|
|---|
| 52 | std::vector<Index> end_remote(size);
|
|---|
| 53 | std::vector< std::vector<vmg_float> > send_buffer_x(size);
|
|---|
| 54 | std::vector< std::vector<vmg_float> > send_buffer_q(size);
|
|---|
| 55 | std::vector< std::vector<vmg_int> > send_buffer_ind(size);
|
|---|
| 56 | std::vector< std::vector<vmg_float> > recv_buffer_x(size);
|
|---|
| 57 | std::vector< std::vector<vmg_float> > recv_buffer_q(size);
|
|---|
| 58 | std::vector< std::vector<vmg_int> > recv_buffer_ind(size);
|
|---|
| 59 |
|
|---|
| 60 | std::memcpy(&global_extent[6*rank], grid.Global().LocalBegin().vec(), 3*sizeof(int));
|
|---|
| 61 | std::memcpy(&global_extent[6*rank+3], grid.Global().LocalEnd().vec(), 3*sizeof(int));
|
|---|
| 62 |
|
|---|
| 63 | MPI_Allgather(MPI_IN_PLACE, 6, MPI_INT, &global_extent.front(), 6, MPI_INT, comm_global);
|
|---|
| 64 |
|
|---|
| 65 | for (int i=0; i<size; ++i) {
|
|---|
| 66 | begin_remote[i] = static_cast<Index>(&global_extent[6*i]);
|
|---|
| 67 | end_remote[i] = static_cast<Index>(&global_extent[6*i+3]);
|
|---|
| 68 | }
|
|---|
| 69 |
|
|---|
| 70 | for (int i=0; i<num_particles_local; ++i) {
|
|---|
| 71 | index = static_cast<Index>((Vector(&x[3*i]) - grid.Extent().Begin()) / grid.Extent().MeshWidth());
|
|---|
| 72 | for (int j=0; j<size; ++j)
|
|---|
| 73 | if (index.IsInBounds(begin_remote[j], end_remote[j])) {
|
|---|
| 74 | send_buffer_x[j].push_back(x[3*i+0]);
|
|---|
| 75 | send_buffer_x[j].push_back(x[3*i+1]);
|
|---|
| 76 | send_buffer_x[j].push_back(x[3*i+2]);
|
|---|
| 77 | send_buffer_q[j].push_back(q[i]);
|
|---|
| 78 | send_buffer_ind[j].push_back(i);
|
|---|
| 79 | break;
|
|---|
| 80 | }
|
|---|
| 81 | }
|
|---|
| 82 |
|
|---|
| 83 | /*
|
|---|
| 84 | * Communicate which process gets how many particles
|
|---|
| 85 | */
|
|---|
| 86 | for (int i=0; i<size; ++i)
|
|---|
| 87 | send_sizes[i] = send_buffer_q[i].size();
|
|---|
| 88 |
|
|---|
| 89 | MPI_Alltoall(&send_sizes.front(), 1, MPI_INT, &recv_sizes.front(), 1, MPI_INT, comm_global);
|
|---|
| 90 |
|
|---|
| 91 | assert(RequestsPending() == 0);
|
|---|
| 92 |
|
|---|
| 93 | /*
|
|---|
| 94 | * Send particles
|
|---|
| 95 | */
|
|---|
| 96 | for (int i=0; i<size; ++i) {
|
|---|
| 97 |
|
|---|
| 98 | if (!send_buffer_q[i].empty()) {
|
|---|
| 99 | MPI_Isend(&send_buffer_x[i].front(), send_buffer_x[i].size(), MPI_DOUBLE, i, 0, comm_global, &Request());
|
|---|
| 100 | MPI_Isend(&send_buffer_q[i].front(), send_buffer_q[i].size(), MPI_DOUBLE, i, 1, comm_global, &Request());
|
|---|
| 101 | MPI_Isend(&send_buffer_ind[i].front(), send_buffer_ind[i].size(), MPI_INT, i, 2, comm_global, &Request());
|
|---|
| 102 | }
|
|---|
| 103 |
|
|---|
| 104 | #ifndef VMG_ONE_SIDED
|
|---|
| 105 | receiver[i] = send_buffer_q[i].size();
|
|---|
| 106 | #endif
|
|---|
| 107 | }
|
|---|
| 108 |
|
|---|
| 109 | /*
|
|---|
| 110 | * Receive particles
|
|---|
| 111 | */
|
|---|
| 112 | for (int i=0; i<size; ++i) {
|
|---|
| 113 |
|
|---|
| 114 | if (recv_sizes[i] > 0) {
|
|---|
| 115 |
|
|---|
| 116 | recv_buffer_x[i].resize(3*recv_sizes[i]);
|
|---|
| 117 | recv_buffer_q[i].resize(recv_sizes[i]);
|
|---|
| 118 | recv_buffer_ind[i].resize(recv_sizes[i]);
|
|---|
| 119 |
|
|---|
| 120 | MPI_Irecv(&recv_buffer_x[i].front(), 3*recv_sizes[i], MPI_DOUBLE, i, 0, comm_global, &Request());
|
|---|
| 121 | MPI_Irecv(&recv_buffer_q[i].front(), recv_sizes[i], MPI_DOUBLE, i, 1, comm_global, &Request());
|
|---|
| 122 | MPI_Irecv(&recv_buffer_ind[i].front(), recv_sizes[i], MPI_INT, i, 2, comm_global, &Request());
|
|---|
| 123 |
|
|---|
| 124 | }
|
|---|
| 125 |
|
|---|
| 126 | }
|
|---|
| 127 |
|
|---|
| 128 | WaitAll();
|
|---|
| 129 |
|
|---|
| 130 | particles.clear();
|
|---|
| 131 |
|
|---|
| 132 | for (int i=0; i<size; ++i)
|
|---|
| 133 | for (int j=0; j<recv_sizes[i]; ++j)
|
|---|
| 134 | particles.push_back(Particle(&recv_buffer_x[i][3*j], recv_buffer_q[i][j], 0.0, 0.0, i, recv_buffer_ind[i][j]));
|
|---|
| 135 | }
|
|---|
| 136 |
|
|---|
| 137 | void Particle::CommMPI::CommParticlesBack(std::list<Particle>& particles)
|
|---|
| 138 | {
|
|---|
| 139 | std::list<Particle>::iterator iter;
|
|---|
| 140 |
|
|---|
| 141 | #ifdef VMG_ONE_SIDED
|
|---|
| 142 | if (!win_created) {
|
|---|
| 143 | vmg_float* p = MG::GetFactory().GetObjectStorageArray<vmg_float>("PARTICLE_POTENTIAL_ARRAY");
|
|---|
| 144 | const vmg_int& num_particles_local = MG::GetFactory().GetObjectStorageVal<vmg_int>("PARTICLE_NUM_LOCAL");
|
|---|
| 145 | MPI_Win_create(p, num_particles_local*sizeof(vmg_float), sizeof(vmg_float), info, comm_global, &win);
|
|---|
| 146 | win_created = true;
|
|---|
| 147 | }
|
|---|
| 148 |
|
|---|
| 149 | MPI_Win_fence(MPI_MODE_NOPRECEDE, win);
|
|---|
| 150 |
|
|---|
| 151 | for (iter=particles.begin(); iter!=particles.end(); ++iter)
|
|---|
| 152 | MPI_Put(&iter->Pot(), 1, MPI_DOUBLE, iter->Rank(), iter->Index(), 1, MPI_DOUBLE, win);
|
|---|
| 153 |
|
|---|
| 154 | MPI_Win_fence(MPI_MODE_NOSTORE | MPI_MODE_NOSUCCEED, win);
|
|---|
| 155 | #else
|
|---|
| 156 | int rank, size;
|
|---|
| 157 | MPI_Comm_rank(comm_global, &rank);
|
|---|
| 158 | MPI_Comm_size(comm_global, &size);
|
|---|
| 159 |
|
|---|
| 160 | std::vector< std::vector<vmg_float> > send_buffer_float(size);
|
|---|
| 161 | std::vector< std::vector<vmg_float> > recv_buffer_float(size);
|
|---|
| 162 | std::vector< std::vector<vmg_int> > send_buffer_index(size);
|
|---|
| 163 | std::vector< std::vector<vmg_int> > recv_buffer_index(size);
|
|---|
| 164 |
|
|---|
| 165 | vmg_int* size_receive = MG::GetFactory().GetObjectStorageArray<vmg_int>("PARTICLE_RECEIVER_ARRAY");
|
|---|
| 166 | vmg_float* p = MG::GetFactory().GetObjectStorageArray<vmg_float>("PARTICLE_POTENTIAL_ARRAY");
|
|---|
| 167 | vmg_float* f = MG::GetFactory().GetObjectStorageArray<vmg_float>("PARTICLE_FIELD_ARRAY");
|
|---|
| 168 |
|
|---|
| 169 | // Build send buffer
|
|---|
| 170 | for (iter=particles.begin(); iter!=particles.end(); ++iter) {
|
|---|
| 171 | send_buffer_float[iter->Rank()].push_back(iter->Pot());
|
|---|
| 172 | send_buffer_float[iter->Rank()].push_back(iter->Field()[0]);
|
|---|
| 173 | send_buffer_float[iter->Rank()].push_back(iter->Field()[1]);
|
|---|
| 174 | send_buffer_float[iter->Rank()].push_back(iter->Field()[2]);
|
|---|
| 175 | send_buffer_index[iter->Rank()].push_back(iter->Index());
|
|---|
| 176 | }
|
|---|
| 177 |
|
|---|
| 178 | // Send potentials
|
|---|
| 179 | for (int i=0; i<size; ++i) {
|
|---|
| 180 | if (!send_buffer_float[i].empty()) {
|
|---|
| 181 | MPI_Isend(&send_buffer_float[i].front(), send_buffer_float[i].size(), MPI_DOUBLE, i, 699+rank, comm_global, &Request());
|
|---|
| 182 | MPI_Isend(&send_buffer_index[i].front(), send_buffer_index[i].size(), MPI_INT, i, 32111+rank, comm_global, &Request());
|
|---|
| 183 | }
|
|---|
| 184 | }
|
|---|
| 185 |
|
|---|
| 186 | //Receive potentials
|
|---|
| 187 | for (int i=0; i<size; ++i) {
|
|---|
| 188 | if (size_receive[i] > 0) {
|
|---|
| 189 | recv_buffer_float[i].resize(4*size_receive[i]);
|
|---|
| 190 | recv_buffer_index[i].resize(size_receive[i]);
|
|---|
| 191 | MPI_Irecv(&recv_buffer_float[i].front(), 4*size_receive[i], MPI_DOUBLE, i, 699+i, comm_global, &Request());
|
|---|
| 192 | MPI_Irecv(&recv_buffer_index[i].front(), size_receive[i], MPI_INT, i, 32111+i, comm_global, &Request());
|
|---|
| 193 | }
|
|---|
| 194 | }
|
|---|
| 195 |
|
|---|
| 196 | WaitAll();
|
|---|
| 197 |
|
|---|
| 198 | // Add potential values
|
|---|
| 199 | for (int i=0; i<size; ++i)
|
|---|
| 200 | for (vmg_int j=0; j<size_receive[i]; ++j) {
|
|---|
| 201 | p[recv_buffer_index[i][j]] = recv_buffer_float[i][4*j];
|
|---|
| 202 | std::memcpy(&f[recv_buffer_index[i][j]], &recv_buffer_float[i][4*j+1], 3*sizeof(vmg_float));
|
|---|
| 203 | }
|
|---|
| 204 | #endif
|
|---|
| 205 |
|
|---|
| 206 | }
|
|---|
| 207 |
|
|---|
| 208 | void Particle::CommMPI::CommLCListToGhosts(LinkedCellList& lc)
|
|---|
| 209 | {
|
|---|
| 210 | VMG::MPI::DatatypesLocal types(lc, comm_global, false);
|
|---|
| 211 | std::vector<int> send_size(types.NB().size());
|
|---|
| 212 | vmg_int recv_size;
|
|---|
| 213 | std::list<Particle*>::iterator iter;
|
|---|
| 214 | Index ind;
|
|---|
| 215 | Vector offset;
|
|---|
| 216 |
|
|---|
| 217 | const Vector halo_length = lc.Local().HaloSize1() * lc.Extent().MeshWidth();
|
|---|
| 218 |
|
|---|
| 219 | lc.ClearHalo();
|
|---|
| 220 |
|
|---|
| 221 | for (unsigned int i=0; i<types.NB().size(); ++i)
|
|---|
| 222 | if (types.NB()[i].Feasible()) {
|
|---|
| 223 |
|
|---|
| 224 | for (int j=0; j<3; ++j)
|
|---|
| 225 | if ((types.Offset()[i][j] < 0 && lc.Global().LocalBegin()[j] == 0) ||
|
|---|
| 226 | (types.Offset()[i][j] > 0 && lc.Global().LocalEnd()[j] == lc.Global().GlobalSize()[j]))
|
|---|
| 227 | offset[j] = -1.0 * types.Offset()[i][j] * lc.Extent().Size()[j];
|
|---|
| 228 | else
|
|---|
| 229 | offset[j] = 0.0;
|
|---|
| 230 |
|
|---|
| 231 | for (ind.X() = types.NB()[i].Starts().X(); ind.X() < types.NB()[i].Starts().X()+types.NB()[i].Subsizes().X(); ++ind.X())
|
|---|
| 232 | for (ind.Y() = types.NB()[i].Starts().Y(); ind.Y() < types.NB()[i].Starts().Y()+types.NB()[i].Subsizes().Y(); ++ind.Y())
|
|---|
| 233 | for (ind.Z() = types.NB()[i].Starts().Z(); ind.Z() < types.NB()[i].Starts().Z()+types.NB()[i].Subsizes().Z(); ++ind.Z())
|
|---|
| 234 | for (iter=lc(ind).begin(); iter!=lc(ind).end(); ++iter) {
|
|---|
| 235 |
|
|---|
| 236 | for (int j=0; j<3; ++j)
|
|---|
| 237 | types.NB()[i].Buffer().push_back((*iter)->Pos()[j] + offset[j]);
|
|---|
| 238 | types.NB()[i].Buffer().push_back((*iter)->Charge());
|
|---|
| 239 |
|
|---|
| 240 | assert(lc.Extent().Begin().IsComponentwiseLessOrEqual((*iter)->Pos()));
|
|---|
| 241 | assert(lc.Extent().End().IsComponentwiseGreaterOrEqual((*iter)->Pos()));
|
|---|
| 242 | assert(lc.Extent().Begin().IsComponentwiseLessOrEqual((*iter)->Pos() + offset + halo_length));
|
|---|
| 243 | assert(lc.Extent().End().IsComponentwiseGreaterOrEqual((*iter)->Pos() + offset - halo_length));
|
|---|
| 244 | }
|
|---|
| 245 |
|
|---|
| 246 | send_size[i] = types.NB()[i].Buffer().size();
|
|---|
| 247 | MPI_Isend(&send_size[i], 1, MPI_INT, types.NB()[i].Rank(), 2048+types.NB()[i].TagSend(), comm_global, &Request());
|
|---|
| 248 |
|
|---|
| 249 | if (send_size[i] > 0)
|
|---|
| 250 | MPI_Isend(&types.NB()[i].Buffer().front(), send_size[i], MPI_DOUBLE,
|
|---|
| 251 | types.NB()[i].Rank(), 4096+types.NB()[i].TagSend(),
|
|---|
| 252 | comm_global, &Request());
|
|---|
| 253 | }
|
|---|
| 254 |
|
|---|
| 255 | for (unsigned int i=0; i<types.Halo().size(); ++i)
|
|---|
| 256 | if (types.Halo()[i].Feasible()) {
|
|---|
| 257 | MPI_Recv(&recv_size, 1, MPI_INT, types.Halo()[i].Rank(), 2048+types.Halo()[i].TagReceive(), comm_global, MPI_STATUS_IGNORE);
|
|---|
| 258 | if (recv_size > 0) {
|
|---|
| 259 | types.Halo()[i].Buffer().resize(recv_size);
|
|---|
| 260 | MPI_Irecv(&types.Halo()[i].Buffer().front(), recv_size, MPI_DOUBLE,
|
|---|
| 261 | types.Halo()[i].Rank(), 4096+types.Halo()[i].TagReceive(),
|
|---|
| 262 | comm_global, &Request());
|
|---|
| 263 | }
|
|---|
| 264 | }
|
|---|
| 265 |
|
|---|
| 266 | WaitAll();
|
|---|
| 267 |
|
|---|
| 268 | for (unsigned int i=0; i<types.Halo().size(); ++i)
|
|---|
| 269 | for (unsigned int j=0; j<types.Halo()[i].Buffer().size(); j+=4)
|
|---|
| 270 | lc.AddParticleToHalo(&types.Halo()[i].Buffer()[j], types.Halo()[i].Buffer()[j+3]);
|
|---|
| 271 | }
|
|---|
| 272 |
|
|---|
| 273 | #endif /* HAVE_MPI */
|
|---|