source: src/units/particle/comm_mpi_particle.cpp@ f003a9

Last change on this file since f003a9 was f003a9, checked in by Julian Iseringhausen <isering@…>, 14 years ago

Refactored vmg in order to separate the core library and the particle simulation part properly.

git-svn-id: https://svn.version.fz-juelich.de/scafacos/trunk@1798 5161e1c8-67bf-11de-9fd5-51895aff932f

  • Property mode set to 100644
File size: 9.9 KB
Line 
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
26using namespace VMG;
27
28void 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 MPI_Datatype dt_temp;
60
61 std::memcpy(&global_extent[6*rank], grid.Global().LocalBegin().vec(), 3*sizeof(int));
62 std::memcpy(&global_extent[6*rank+3], grid.Global().LocalEnd().vec(), 3*sizeof(int));
63
64 MPI_Allgather(MPI_IN_PLACE, 6, MPI_INT, &global_extent.front(), 6, MPI_INT, comm_global);
65
66 for (int i=0; i<size; ++i) {
67 begin_remote[i] = static_cast<Index>(&global_extent[6*i]);
68 end_remote[i] = static_cast<Index>(&global_extent[6*i+3]);
69 }
70
71 for (int i=0; i<num_particles_local; ++i) {
72 index = static_cast<Index>((Vector(&x[3*i]) - grid.Extent().Begin()) / grid.Extent().MeshWidth());
73 for (int j=0; j<size; ++j)
74 if (index.IsInBounds(begin_remote[j], end_remote[j])) {
75 send_buffer_x[j].push_back(x[3*i+0]);
76 send_buffer_x[j].push_back(x[3*i+1]);
77 send_buffer_x[j].push_back(x[3*i+2]);
78 send_buffer_q[j].push_back(q[i]);
79 send_buffer_ind[j].push_back(i);
80 break;
81 }
82 }
83
84 /*
85 * Communicate which process gets how many particles
86 */
87 for (int i=0; i<size; ++i)
88 send_sizes[i] = send_buffer_q[i].size();
89
90 MPI_Alltoall(&send_sizes.front(), 1, MPI_INT, &recv_sizes.front(), 1, MPI_INT, comm_global);
91
92 assert(RequestsPending() == 0);
93
94 /*
95 * Send particles
96 */
97 for (int i=0; i<size; ++i) {
98
99 if (!send_buffer_q[i].empty()) {
100 MPI_Isend(&send_buffer_x[i].front(), send_buffer_x[i].size(), MPI_DOUBLE, i, 0, comm_global, &Request());
101 MPI_Isend(&send_buffer_q[i].front(), send_buffer_q[i].size(), MPI_DOUBLE, i, 1, comm_global, &Request());
102 MPI_Isend(&send_buffer_ind[i].front(), send_buffer_ind[i].size(), MPI_INT, i, 2, comm_global, &Request());
103 }
104
105#ifndef VMG_ONE_SIDED
106 receiver[i] = send_buffer_q[i].size();
107#endif
108 }
109
110 /*
111 * Receive particles
112 */
113 for (int i=0; i<size; ++i) {
114
115 if (recv_sizes[i] > 0) {
116
117 recv_buffer_x[i].resize(3*recv_sizes[i]);
118 recv_buffer_q[i].resize(recv_sizes[i]);
119 recv_buffer_ind[i].resize(recv_sizes[i]);
120
121 MPI_Irecv(&recv_buffer_x[i].front(), 3*recv_sizes[i], MPI_DOUBLE, i, 0, comm_global, &Request());
122 MPI_Irecv(&recv_buffer_q[i].front(), recv_sizes[i], MPI_DOUBLE, i, 1, comm_global, &Request());
123 MPI_Irecv(&recv_buffer_ind[i].front(), recv_sizes[i], MPI_INT, i, 2, comm_global, &Request());
124
125 }
126
127 }
128
129 WaitAll();
130
131 particles.clear();
132
133 for (int i=0; i<size; ++i)
134 for (int j=0; j<recv_sizes[i]; ++j)
135 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]));
136}
137
138void Particle::CommMPI::CommParticlesBack(std::list<Particle>& particles)
139{
140 std::list<Particle>::iterator iter;
141
142#ifdef VMG_ONE_SIDED
143 if (!win_created) {
144 vmg_float* p = MG::GetFactory().GetObjectStorageArray<vmg_float>("PARTICLE_POTENTIAL_ARRAY");
145 const vmg_int& num_particles_local = MG::GetFactory().GetObjectStorageVal<vmg_int>("PARTICLE_NUM_LOCAL");
146 MPI_Win_create(p, num_particles_local*sizeof(vmg_float), sizeof(vmg_float), info, comm_global, &win);
147 win_created = true;
148 }
149
150 MPI_Win_fence(MPI_MODE_NOPRECEDE, win);
151
152 for (iter=particles.begin(); iter!=particles.end(); ++iter)
153 MPI_Put(&iter->Pot(), 1, MPI_DOUBLE, iter->Rank(), iter->Index(), 1, MPI_DOUBLE, win);
154
155 MPI_Win_fence(MPI_MODE_NOSTORE | MPI_MODE_NOSUCCEED, win);
156#else
157 int rank, size;
158 MPI_Comm_rank(comm_global, &rank);
159 MPI_Comm_size(comm_global, &size);
160
161 std::vector< std::vector<vmg_float> > send_buffer_float(size);
162 std::vector< std::vector<vmg_float> > recv_buffer_float(size);
163 std::vector< std::vector<vmg_int> > send_buffer_index(size);
164 std::vector< std::vector<vmg_int> > recv_buffer_index(size);
165
166 vmg_int* size_receive = MG::GetFactory().GetObjectStorageArray<vmg_int>("PARTICLE_RECEIVER_ARRAY");
167 vmg_float* p = MG::GetFactory().GetObjectStorageArray<vmg_float>("PARTICLE_POTENTIAL_ARRAY");
168 vmg_float* f = MG::GetFactory().GetObjectStorageArray<vmg_float>("PARTICLE_FIELD_ARRAY");
169
170 // Build send buffer
171 for (iter=particles.begin(); iter!=particles.end(); ++iter) {
172 send_buffer_float[iter->Rank()].push_back(iter->Pot());
173 send_buffer_float[iter->Rank()].push_back(iter->Field()[0]);
174 send_buffer_float[iter->Rank()].push_back(iter->Field()[1]);
175 send_buffer_float[iter->Rank()].push_back(iter->Field()[2]);
176 send_buffer_index[iter->Rank()].push_back(iter->Index());
177 }
178
179 // Send potentials
180 for (int i=0; i<size; ++i) {
181 if (!send_buffer_float[i].empty()) {
182 MPI_Isend(&send_buffer_float[i].front(), send_buffer_float[i].size(), MPI_DOUBLE, i, 699+rank, comm_global, &Request());
183 MPI_Isend(&send_buffer_index[i].front(), send_buffer_index[i].size(), MPI_INT, i, 32111+rank, comm_global, &Request());
184 }
185 }
186
187 //Receive potentials
188 for (int i=0; i<size; ++i) {
189 if (size_receive[i] > 0) {
190 recv_buffer_float[i].resize(4*size_receive[i]);
191 recv_buffer_index[i].resize(size_receive[i]);
192 MPI_Irecv(&recv_buffer_float[i].front(), 4*size_receive[i], MPI_DOUBLE, i, 699+i, comm_global, &Request());
193 MPI_Irecv(&recv_buffer_index[i].front(), size_receive[i], MPI_INT, i, 32111+i, comm_global, &Request());
194 }
195 }
196
197 WaitAll();
198
199 // Add potential values
200 for (int i=0; i<size; ++i)
201 for (unsigned int j=0; j<size_receive[i]; ++j) {
202 p[recv_buffer_index[i][j]] = recv_buffer_float[i][4*j];
203 std::memcpy(&f[recv_buffer_index[i][j]], &recv_buffer_float[i][4*j+1], 3*sizeof(vmg_float));
204 }
205#endif
206
207}
208
209void Particle::CommMPI::CommLCListToGhosts(LinkedCellList& lc)
210{
211 VMG::MPI::DatatypesLocal types(lc, comm_global, false);
212 std::vector<int> send_size(types.NB().size());
213 vmg_int recv_size;
214 std::list<Particle*>::iterator iter;
215 Index ind;
216 Vector offset;
217
218 const Vector halo_length = lc.Local().HaloSize1() * lc.Extent().MeshWidth();
219
220 lc.ClearHalo();
221
222 for (unsigned int i=0; i<types.NB().size(); ++i)
223 if (types.NB()[i].Feasible()) {
224
225 for (int j=0; j<3; ++j)
226 if ((types.Offset()[i][j] < 0 && lc.Global().LocalBegin()[j] == 0) ||
227 (types.Offset()[i][j] > 0 && lc.Global().LocalEnd()[j] == lc.Global().GlobalSize()[j]))
228 offset[j] = -1.0 * types.Offset()[i][j] * lc.Extent().Size()[j];
229 else
230 offset[j] = 0.0;
231
232 for (ind.X() = types.NB()[i].Starts().X(); ind.X() < types.NB()[i].Starts().X()+types.NB()[i].Subsizes().X(); ++ind.X())
233 for (ind.Y() = types.NB()[i].Starts().Y(); ind.Y() < types.NB()[i].Starts().Y()+types.NB()[i].Subsizes().Y(); ++ind.Y())
234 for (ind.Z() = types.NB()[i].Starts().Z(); ind.Z() < types.NB()[i].Starts().Z()+types.NB()[i].Subsizes().Z(); ++ind.Z())
235 for (iter=lc(ind).begin(); iter!=lc(ind).end(); ++iter) {
236
237 for (int j=0; j<3; ++j)
238 types.NB()[i].Buffer().push_back((*iter)->Pos()[j] + offset[j]);
239 types.NB()[i].Buffer().push_back((*iter)->Charge());
240
241 assert(lc.Extent().Begin().IsComponentwiseLessOrEqual((*iter)->Pos()));
242 assert(lc.Extent().End().IsComponentwiseGreaterOrEqual((*iter)->Pos()));
243 assert(lc.Extent().Begin().IsComponentwiseLessOrEqual((*iter)->Pos() + offset + halo_length));
244 assert(lc.Extent().End().IsComponentwiseGreaterOrEqual((*iter)->Pos() + offset - halo_length));
245 }
246
247 send_size[i] = types.NB()[i].Buffer().size();
248 MPI_Isend(&send_size[i], 1, MPI_INT, types.NB()[i].Rank(), 2048+types.NB()[i].TagSend(), comm_global, &Request());
249
250 if (send_size[i] > 0)
251 MPI_Isend(&types.NB()[i].Buffer().front(), send_size[i], MPI_DOUBLE,
252 types.NB()[i].Rank(), 4096+types.NB()[i].TagSend(),
253 comm_global, &Request());
254 }
255
256 for (unsigned int i=0; i<types.Halo().size(); ++i)
257 if (types.Halo()[i].Feasible()) {
258 MPI_Recv(&recv_size, 1, MPI_INT, types.Halo()[i].Rank(), 2048+types.Halo()[i].TagReceive(), comm_global, MPI_STATUS_IGNORE);
259 if (recv_size > 0) {
260 types.Halo()[i].Buffer().resize(recv_size);
261 MPI_Irecv(&types.Halo()[i].Buffer().front(), recv_size, MPI_DOUBLE,
262 types.Halo()[i].Rank(), 4096+types.Halo()[i].TagReceive(),
263 comm_global, &Request());
264 }
265 }
266
267 WaitAll();
268
269 for (unsigned int i=0; i<types.Halo().size(); ++i)
270 for (unsigned int j=0; j<types.Halo()[i].Buffer().size(); j+=4)
271 lc.AddParticleToHalo(&types.Halo()[i].Buffer()[j], types.Halo()[i].Buffer()[j+3]);
272}
273
274#endif /* HAVE_MPI */
Note: See TracBrowser for help on using the repository browser.