source: src/comm/comm_mpi.cpp@ 716da7

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

Fix energy calculation.

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

  • Property mode set to 100644
File size: 39.4 KB
Line 
1/**
2 * @file comm_mpi.cpp
3 * @author Julian Iseringhausen <isering@ins.uni-bonn.de>
4 * @date Wed Jun 16 13:21:06 2010
5 *
6 * @brief Class for MPI-based 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#ifdef HAVE_BOOST_FILESYSTEM
23#include <boost/filesystem.hpp>
24namespace fs = boost::filesystem;
25#endif
26
27#ifdef HAVE_VTK
28#include <vtkAbstractArray.h>
29#include <vtkImageData.h>
30#include <vtkPointData.h>
31#include <vtkSmartPointer.h>
32#include <vtkXMLImageDataWriter.h>
33#endif
34
35#include <cstring>
36#include <sstream>
37
38#include "base/helper.hpp"
39#include "base/linked_cell_list.hpp"
40#include "base/tuple.hpp"
41#include "comm/comm_mpi.hpp"
42#include "comm/mpi/datatypes_local.hpp"
43#include "grid/grid.hpp"
44#include "grid/multigrid.hpp"
45#include "grid/tempgrid.hpp"
46#include "mg.hpp"
47#include "base/timer.hpp"
48static char print_buffer[512];
49
50using namespace VMG;
51
52void CommMPI::IsendAll(Grid& grid, std::vector<VMG::MPI::Datatype>& types, const MPI_Comm& comm, const int& tag_start)
53{
54 for (std::vector<VMG::MPI::Datatype>::const_iterator iter=types.begin(); iter!=types.end(); ++iter)
55 iter->Isend(grid, tag_start, comm, Request());
56}
57
58void CommMPI::IrecvAll(Grid& grid, std::vector<VMG::MPI::Datatype>& types, const MPI_Comm& comm, const int& tag_start)
59{
60 for (std::vector<VMG::MPI::Datatype>::const_iterator iter=types.begin(); iter!=types.end(); ++iter)
61 iter->Irecv(grid, tag_start, comm, Request());
62}
63
64void CommMPI::IsendAllBuffered(const Grid& grid, std::vector<VMG::MPI::Datatype>& types, const MPI_Comm& comm, const int& tag_start)
65{
66 for (std::vector<VMG::MPI::Datatype>::iterator iter=types.begin(); iter!=types.end(); ++iter)
67 iter->IsendBuffered(grid, tag_start, comm, Request());
68}
69
70void CommMPI::IrecvAllBuffered(std::vector<VMG::MPI::Datatype>& types, const MPI_Comm& comm, const int& tag_start)
71{
72 for (std::vector<VMG::MPI::Datatype>::iterator iter=types.begin(); iter!=types.end(); ++iter)
73 iter->IrecvBuffered(tag_start, comm, Request());
74}
75
76void CommMPI::ReplaceBufferAll(Grid& grid, const std::vector<VMG::MPI::Datatype>& types)
77{
78 for (std::vector<VMG::MPI::Datatype>::const_iterator iter=types.begin(); iter!= types.end(); ++iter)
79 iter->GridReplace(grid);
80}
81
82void CommMPI::AddBufferAll(Grid& grid, const std::vector<VMG::MPI::Datatype>& types)
83{
84 for (std::vector<VMG::MPI::Datatype>::const_iterator iter=types.begin(); iter!= types.end(); ++iter)
85 iter->GridSum(grid);
86}
87
88void CommMPI::CommSubgrid(Grid& grid_old, Grid& grid_new, const int& direction)
89{
90 MPI_Comm comm = settings.CommunicatorGlobal(grid_old);
91 if (comm != MPI_COMM_NULL) {
92 VMG::MPI::DatatypesGlobal& datatypes = settings.DatatypesGlobal(grid_old, grid_new, direction);
93 IrecvAllBuffered(datatypes.Receive(), comm, 0411);
94 IsendAllBuffered(grid_old, datatypes.Send(), comm, 0411);
95 WaitAll();
96 ReplaceBufferAll(grid_new, datatypes.Receive());
97 }
98}
99
100void CommMPI::CommAddSubgrid(Grid& grid_old, Grid& grid_new, const int& direction)
101{
102 MPI_Comm comm = settings.CommunicatorGlobal(grid_old);
103 if (comm != MPI_COMM_NULL) {
104 VMG::MPI::DatatypesGlobal& datatypes = settings.DatatypesGlobal(grid_old, grid_new, direction);
105 IrecvAllBuffered(datatypes.Receive(), comm, 1806);
106 IsendAllBuffered(grid_old, datatypes.Send(), comm, 1806);
107 WaitAll();
108 AddBufferAll(grid_new, datatypes.Receive());
109 }
110}
111
112void CommMPI::CommToGhosts(Grid& grid)
113{
114 MPI_Comm comm = settings.CommunicatorLocal(grid);
115 if (comm != MPI_COMM_NULL) {
116 VMG::MPI::DatatypesLocal& types = settings.DatatypesLocal(grid);
117 IrecvAllBuffered(types.Halo(), comm, 2310);
118 IsendAllBuffered(grid, types.NB(), comm, 2310);
119 WaitAll();
120 ReplaceBufferAll(grid, types.Halo());
121 }
122}
123
124void CommMPI::CommToGhostsAsyncStart(Grid& grid)
125{
126 MPI_Comm comm = settings.CommunicatorLocal(grid);
127 if (comm != MPI_COMM_NULL) {
128 VMG::MPI::DatatypesLocal& types = settings.DatatypesLocal(grid);
129 IrecvAllBuffered(types.Halo(), comm, 2412);
130 IsendAllBuffered(grid, types.NB(), comm, 2412);
131 TestAll();
132 }
133}
134
135void CommMPI::CommToGhostsAsyncFinish(Grid& grid)
136{
137 WaitAll();
138 ReplaceBufferAll(grid, settings.DatatypesLocal(grid).Halo());
139}
140
141void CommMPI::CommFromGhosts(Grid& grid)
142{
143 MPI_Comm comm = settings.CommunicatorLocal(grid);
144 if (comm != MPI_COMM_NULL) {
145 VMG::MPI::DatatypesLocal& types = settings.DatatypesLocal(grid);
146 IrecvAllBuffered(types.NB(), comm, 1337);
147 IsendAllBuffered(grid, types.Halo(), comm, 1337);
148 WaitAll();
149 AddBufferAll(grid, types.NB());
150 }
151}
152
153void CommMPI::CommFromGhostsAsyncStart(Grid& grid)
154{
155 MPI_Comm comm = settings.CommunicatorLocal(grid);
156 if (comm != MPI_COMM_NULL) {
157 VMG::MPI::DatatypesLocal& types = settings.DatatypesLocal(grid);
158 IrecvAllBuffered(types.NB(), comm, 0xc0ffee);
159 IsendAllBuffered(grid, types.Halo(), comm, 0xc0ffee);
160 TestAll();
161 }
162}
163
164void CommMPI::CommFromGhostsAsyncFinish(Grid& grid)
165{
166 WaitAll();
167 AddBufferAll(grid, settings.DatatypesLocal(grid).NB());
168}
169
170void CommMPI::CommParticles(const Grid& grid, std::list<Particle::Particle>& particles)
171{
172 Factory& factory = MG::GetFactory();
173
174 const vmg_int& num_particles_local = factory.GetObjectStorageVal<vmg_int>("PARTICLE_NUM_LOCAL");
175 vmg_float* x = factory.GetObjectStorageArray<vmg_float>("PARTICLE_POS_ARRAY");
176 vmg_float* q = factory.GetObjectStorageArray<vmg_float>("PARTICLE_CHARGE_ARRAY");
177
178 int rank, size;
179 MPI_Comm_rank(comm_global, &rank);
180 MPI_Comm_size(comm_global, &size);
181
182#ifndef VMG_ONE_SIDED
183 vmg_int* receiver;
184 if (!factory.TestObject("PARTICLE_RECEIVER_ARRAY"))
185 new ObjectStorageArray<vmg_int>("PARTICLE_RECEIVER_ARRAY", size);
186 receiver = factory.GetObjectStorageArray<vmg_int>("PARTICLE_RECEIVER_ARRAY");
187#endif
188
189 Index index;
190 std::vector<int> global_extent(6*size);
191 std::vector<int> send_sizes(size);
192 std::vector<int> recv_sizes(size);
193 std::vector<Index> begin_remote(size);
194 std::vector<Index> end_remote(size);
195 std::vector< std::vector<vmg_float> > send_buffer_x(size);
196 std::vector< std::vector<vmg_float> > send_buffer_q(size);
197 std::vector< std::vector<vmg_int> > send_buffer_ind(size);
198 std::vector< std::vector<vmg_float> > recv_buffer_x(size);
199 std::vector< std::vector<vmg_float> > recv_buffer_q(size);
200 std::vector< std::vector<vmg_int> > recv_buffer_ind(size);
201 MPI_Datatype dt_temp;
202
203 std::memcpy(&global_extent[6*rank], grid.Global().LocalBegin().vec(), 3*sizeof(int));
204 std::memcpy(&global_extent[6*rank+3], grid.Global().LocalEnd().vec(), 3*sizeof(int));
205
206 MPI_Allgather(MPI_IN_PLACE, 6, MPI_INT, &global_extent.front(), 6, MPI_INT, comm_global);
207
208 for (int i=0; i<size; ++i) {
209 begin_remote[i] = static_cast<Index>(&global_extent[6*i]);
210 end_remote[i] = static_cast<Index>(&global_extent[6*i+3]);
211 }
212
213 for (int i=0; i<num_particles_local; ++i) {
214 index = static_cast<Index>((Vector(&x[3*i]) - grid.Extent().Begin()) / grid.Extent().MeshWidth());
215 for (int j=0; j<size; ++j)
216 if (index.IsInBounds(begin_remote[j], end_remote[j])) {
217 send_buffer_x[j].push_back(x[3*i+0]);
218 send_buffer_x[j].push_back(x[3*i+1]);
219 send_buffer_x[j].push_back(x[3*i+2]);
220 send_buffer_q[j].push_back(q[i]);
221 send_buffer_ind[j].push_back(i);
222 break;
223 }
224 }
225
226 /*
227 * Communicate which process gets how many particles
228 */
229 for (int i=0; i<size; ++i)
230 send_sizes[i] = send_buffer_q[i].size();
231
232 MPI_Alltoall(&send_sizes.front(), 1, MPI_INT, &recv_sizes.front(), 1, MPI_INT, comm_global);
233
234 assert(RequestsPending() == 0);
235
236 /*
237 * Send particles
238 */
239 for (int i=0; i<size; ++i) {
240
241 if (!send_buffer_q[i].empty()) {
242 MPI_Isend(&send_buffer_x[i].front(), send_buffer_x[i].size(), MPI_DOUBLE, i, 0, comm_global, &Request());
243 MPI_Isend(&send_buffer_q[i].front(), send_buffer_q[i].size(), MPI_DOUBLE, i, 1, comm_global, &Request());
244 MPI_Isend(&send_buffer_ind[i].front(), send_buffer_ind[i].size(), MPI_INT, i, 2, comm_global, &Request());
245 }
246
247#ifndef VMG_ONE_SIDED
248 receiver[i] = send_buffer_q[i].size();
249#endif
250 }
251
252 /*
253 * Receive particles
254 */
255 for (int i=0; i<size; ++i) {
256
257 if (recv_sizes[i] > 0) {
258
259 recv_buffer_x[i].resize(3*recv_sizes[i]);
260 recv_buffer_q[i].resize(recv_sizes[i]);
261 recv_buffer_ind[i].resize(recv_sizes[i]);
262
263 MPI_Irecv(&recv_buffer_x[i].front(), 3*recv_sizes[i], MPI_DOUBLE, i, 0, comm_global, &Request());
264 MPI_Irecv(&recv_buffer_q[i].front(), recv_sizes[i], MPI_DOUBLE, i, 1, comm_global, &Request());
265 MPI_Irecv(&recv_buffer_ind[i].front(), recv_sizes[i], MPI_INT, i, 2, comm_global, &Request());
266
267 }
268
269 }
270
271 WaitAll();
272
273 particles.clear();
274
275 for (int i=0; i<size; ++i)
276 for (int j=0; j<recv_sizes[i]; ++j)
277 particles.push_back(Particle::Particle(&recv_buffer_x[i][3*j], recv_buffer_q[i][j], 0.0, 0.0, i, recv_buffer_ind[i][j]));
278}
279
280void CommMPI::CommParticlesBack(std::list<Particle::Particle>& particles)
281{
282 std::list<Particle::Particle>::iterator iter;
283
284#ifdef VMG_ONE_SIDED
285 if (!win_created) {
286 vmg_float* p = MG::GetFactory().GetObjectStorageArray<vmg_float>("PARTICLE_POTENTIAL_ARRAY");
287 const vmg_int& num_particles_local = MG::GetFactory().GetObjectStorageVal<vmg_int>("PARTICLE_NUM_LOCAL");
288 MPI_Win_create(p, num_particles_local*sizeof(vmg_float), sizeof(vmg_float), info, comm_global, &win);
289 win_created = true;
290 }
291
292 MPI_Win_fence(MPI_MODE_NOPRECEDE, win);
293
294 for (iter=particles.begin(); iter!=particles.end(); ++iter)
295 MPI_Put(&iter->Pot(), 1, MPI_DOUBLE, iter->Rank(), iter->Index(), 1, MPI_DOUBLE, win);
296
297 MPI_Win_fence(MPI_MODE_NOSTORE | MPI_MODE_NOSUCCEED, win);
298#else
299 int rank, size;
300 MPI_Comm_rank(comm_global, &rank);
301 MPI_Comm_size(comm_global, &size);
302
303 std::vector< std::vector<vmg_float> > send_buffer_float(size);
304 std::vector< std::vector<vmg_float> > recv_buffer_float(size);
305 std::vector< std::vector<vmg_int> > send_buffer_index(size);
306 std::vector< std::vector<vmg_int> > recv_buffer_index(size);
307
308 vmg_int* size_receive = MG::GetFactory().GetObjectStorageArray<vmg_int>("PARTICLE_RECEIVER_ARRAY");
309 vmg_float* p = MG::GetFactory().GetObjectStorageArray<vmg_float>("PARTICLE_POTENTIAL_ARRAY");
310 vmg_float* f = MG::GetFactory().GetObjectStorageArray<vmg_float>("PARTICLE_FIELD_ARRAY");
311
312 // Build send buffer
313 for (iter=particles.begin(); iter!=particles.end(); ++iter) {
314 send_buffer_float[iter->Rank()].push_back(iter->Pot());
315 send_buffer_float[iter->Rank()].push_back(iter->Field()[0]);
316 send_buffer_float[iter->Rank()].push_back(iter->Field()[1]);
317 send_buffer_float[iter->Rank()].push_back(iter->Field()[2]);
318 send_buffer_index[iter->Rank()].push_back(iter->Index());
319 }
320
321 // Send potentials
322 for (int i=0; i<size; ++i) {
323 if (!send_buffer_float[i].empty()) {
324 MPI_Isend(&send_buffer_float[i].front(), send_buffer_float[i].size(), MPI_DOUBLE, i, 699+rank, comm_global, &Request());
325 MPI_Isend(&send_buffer_index[i].front(), send_buffer_index[i].size(), MPI_INT, i, 32111+rank, comm_global, &Request());
326 }
327 }
328
329 //Receive potentials
330 for (int i=0; i<size; ++i) {
331 if (size_receive[i] > 0) {
332 recv_buffer_float[i].resize(4*size_receive[i]);
333 recv_buffer_index[i].resize(size_receive[i]);
334 MPI_Irecv(&recv_buffer_float[i].front(), 4*size_receive[i], MPI_DOUBLE, i, 699+i, comm_global, &Request());
335 MPI_Irecv(&recv_buffer_index[i].front(), size_receive[i], MPI_INT, i, 32111+i, comm_global, &Request());
336 }
337 }
338
339 WaitAll();
340
341 // Add potential values
342 for (int i=0; i<size; ++i)
343 for (unsigned int j=0; j<size_receive[i]; ++j) {
344 p[recv_buffer_index[i][j]] = recv_buffer_float[i][4*j];
345 std::memcpy(&f[recv_buffer_index[i][j]], &recv_buffer_float[i][4*j+1], 3*sizeof(vmg_float));
346 }
347#endif
348
349}
350
351void CommMPI::CommLCListToGhosts(Particle::LinkedCellList& lc)
352{
353 VMG::MPI::DatatypesLocal types(lc, comm_global);
354 std::vector<int> send_size(types.NB().size());
355 vmg_int recv_size;
356 std::list<Particle::Particle*>::iterator iter;
357 Index ind;
358 Vector offset;
359
360 const Vector halo_length = lc.Local().HaloSize1() * lc.Extent().MeshWidth();
361
362 lc.ClearHalo();
363
364 for (unsigned int i=0; i<types.NB().size(); ++i)
365 if (types.NB()[i].Feasible()) {
366
367 for (int j=0; j<3; ++j)
368 if ((types.Offset()[i][j] < 0 && lc.Global().LocalBegin()[j] == 0) ||
369 (types.Offset()[i][j] > 0 && lc.Global().LocalEnd()[j] == lc.Global().GlobalSize()[j]))
370 offset[j] = -1 * types.Offset()[i][j] * lc.Extent().Size()[j];
371 else
372 offset[j] = 0;
373
374 for (ind.X() = types.NB()[i].Starts().X(); ind.X() < types.NB()[i].Starts().X()+types.NB()[i].Subsizes().X(); ++ind.X())
375 for (ind.Y() = types.NB()[i].Starts().Y(); ind.Y() < types.NB()[i].Starts().Y()+types.NB()[i].Subsizes().Y(); ++ind.Y())
376 for (ind.Z() = types.NB()[i].Starts().Z(); ind.Z() < types.NB()[i].Starts().Z()+types.NB()[i].Subsizes().Z(); ++ind.Z())
377 for (iter=lc(ind).begin(); iter!=lc(ind).end(); ++iter) {
378 for (int j=0; j<3; ++j)
379 types.NB()[i].Buffer().push_back((*iter)->Pos()[j] + offset[j]);
380 types.NB()[i].Buffer().push_back((*iter)->Charge());
381
382 assert(lc.Extent().Begin().IsComponentwiseLessOrEqual((*iter)->Pos()));
383 assert(lc.Extent().End().IsComponentwiseGreaterOrEqual((*iter)->Pos()));
384 assert(lc.Extent().Begin().IsComponentwiseLessOrEqual((*iter)->Pos() + offset + halo_length));
385 assert(lc.Extent().End().IsComponentwiseGreaterOrEqual((*iter)->Pos() + offset - halo_length));
386 }
387
388 send_size[i] = types.NB()[i].Buffer().size();
389 MPI_Isend(&send_size[i], 1, MPI_INT, types.NB()[i].Rank(), 2048+types.NB()[i].TagSend(), comm_global, &Request());
390
391 if (send_size[i] > 0)
392 MPI_Isend(&types.NB()[i].Buffer().front(), send_size[i], MPI_DOUBLE,
393 types.NB()[i].Rank(), 4096+types.NB()[i].TagSend(),
394 comm_global, &Request());
395 }
396
397 for (unsigned int i=0; i<types.Halo().size(); ++i)
398 if (types.Halo()[i].Feasible()) {
399 MPI_Recv(&recv_size, 1, MPI_INT, types.Halo()[i].Rank(), 2048+types.Halo()[i].TagReceive(), comm_global, MPI_STATUS_IGNORE);
400 if (recv_size > 0) {
401 types.Halo()[i].Buffer().resize(recv_size);
402 MPI_Irecv(&types.Halo()[i].Buffer().front(), recv_size, MPI_DOUBLE,
403 types.Halo()[i].Rank(), 4096+types.Halo()[i].TagReceive(),
404 comm_global, &Request());
405 }
406 }
407
408 WaitAll();
409
410 for (unsigned int i=0; i<types.Halo().size(); ++i)
411 for (unsigned int j=0; j<types.Halo()[i].Buffer().size(); j+=4)
412 lc.AddParticleToHalo(&types.Halo()[i].Buffer()[j], types.Halo()[i].Buffer()[j+3]);
413}
414
415void CommMPI::CommLCListFromGhosts(const Grid& grid, Particle::LinkedCellList& lc)
416{
417 VMG::MPI::DatatypesLocal types(lc, comm_global);
418 std::vector< std::vector<vmg_float> > send_buffer(types.NB().size());
419 std::vector< std::vector<vmg_float> > recv_buffer(types.Halo().size());
420 std::vector<vmg_int> send_size(types.NB().size());
421 vmg_int recv_size;
422 std::list<Particle::Particle*>::iterator iter;
423 Index ind;
424
425 for (unsigned int i=0; i<types.Halo().size(); ++i) {
426 assert(types.Halo()[i].Subsizes().X() == types.NB()[i].Subsizes().X());
427 assert(types.Halo()[i].Subsizes().Y() == types.NB()[i].Subsizes().Y());
428 assert(types.Halo()[i].Subsizes().Z() == types.NB()[i].Subsizes().Z());
429
430 for (ind.X() = 0; ind.X() < types.Halo()[i].Subsizes().X(); ++ind.X())
431 for (ind.Y() = 0; ind.Y() < types.Halo()[i].Subsizes().Y(); ++ind.Y())
432 for (ind.Z() = 0; ind.Z() < types.Halo()[i].Subsizes().Z(); ++ind.Z())
433 assert(lc(ind+types.Halo()[i].Starts()).size() == lc(ind+types.NB()[i].Starts()).size());
434 }
435
436 for (unsigned int i=0; i<types.Halo().size(); ++i)
437 if (types.Halo()[i].Feasible()) {
438
439 for (ind.X() = types.Halo()[i].Starts().X(); ind.X() < types.Halo()[i].Starts().X()+types.Halo()[i].Subsizes().X(); ++ind.X())
440 for (ind.Y() = types.Halo()[i].Starts().X(); ind.Y() < types.Halo()[i].Starts().Y()+types.Halo()[i].Subsizes().Y(); ++ind.Y())
441 for (ind.Z() = types.Halo()[i].Starts().X(); ind.Z() < types.Halo()[i].Starts().Z()+types.Halo()[i].Subsizes().Z(); ++ind.Z())
442 for (iter=lc(ind).begin(); iter!=lc(ind).end(); ++iter)
443 send_buffer[i].push_back((*iter)->Pot());
444
445 MPI_Isend(&send_buffer[i].front(), send_buffer[i].size(), MPI_DOUBLE, types.Halo()[i].Rank(), 342523+types.Halo()[i].TagSend(), comm_global, &Request());
446 }
447
448 for (unsigned int i=0; i<types.NB().size(); ++i)
449 if (types.NB()[i].Feasible()) {
450
451 recv_size = 0;
452 for (ind.X() = types.NB()[i].Starts().X(); ind.X() < types.NB()[i].Starts().X()+types.NB()[i].Subsizes().X(); ++ind.X())
453 for (ind.Y() = types.NB()[i].Starts().X(); ind.Y() < types.NB()[i].Starts().Y()+types.NB()[i].Subsizes().Y(); ++ind.Y())
454 for (ind.Z() = types.NB()[i].Starts().X(); ind.Z() < types.NB()[i].Starts().Z()+types.NB()[i].Subsizes().Z(); ++ind.Z())
455 recv_size += lc(ind).size();
456
457 if (recv_size > 0) {
458 recv_buffer[i].resize(recv_size);
459 MPI_Irecv(&recv_buffer[i].front(), recv_size, MPI_DOUBLE, types.NB()[i].Rank(), 342523+types.NB()[i].TagReceive(), comm_global, &Request());
460 }
461 }
462
463 WaitAll();
464
465 for (unsigned int i=0; i<types.NB().size(); ++i) {
466 unsigned int c=0;
467 if (types.NB()[i].Feasible()) {
468 for (ind.X() = types.NB()[i].Starts().X(); ind.X() < types.NB()[i].Starts().X()+types.NB()[i].Subsizes().X(); ++ind.X())
469 for (ind.Y() = types.NB()[i].Starts().X(); ind.Y() < types.NB()[i].Starts().Y()+types.NB()[i].Subsizes().Y(); ++ind.Y())
470 for (ind.Z() = types.NB()[i].Starts().X(); ind.Z() < types.NB()[i].Starts().Z()+types.NB()[i].Subsizes().Z(); ++ind.Z())
471 for (iter=lc(ind).begin(); iter!=lc(ind).end(); ++iter)
472 (*iter)->Pot() += recv_buffer[i][c++];
473 }
474 }
475}
476
477vmg_float CommMPI::GlobalSum(const vmg_float& value)
478{
479 vmg_float result = value;
480 MPI_Allreduce(MPI_IN_PLACE, &result, 1, MPI_DOUBLE, MPI_SUM, comm_global);
481 return result;
482}
483
484vmg_float CommMPI::GlobalSumRoot(const vmg_float& value)
485{
486 vmg_float recv_buffer = value;
487 vmg_float send_buffer = value;
488 MPI_Reduce(&send_buffer, &recv_buffer, 1, MPI_DOUBLE, MPI_SUM, 0, comm_global);
489 return recv_buffer;
490}
491
492void CommMPI::GlobalSumArray(vmg_float* array, const vmg_int& size)
493{
494 MPI_Allreduce(MPI_IN_PLACE, array, size, MPI_DOUBLE, MPI_SUM, comm_global);
495}
496
497vmg_int CommMPI::GlobalSum(const vmg_int& value)
498{
499 vmg_int result = value;
500 MPI_Allreduce(MPI_IN_PLACE, &result, 1, MPI_INT, MPI_SUM, comm_global);
501 return result;
502}
503
504vmg_int CommMPI::GlobalSumRoot(const vmg_int& value)
505{
506 vmg_int recv_buffer = value;
507 vmg_int send_buffer = value;
508 MPI_Reduce(&send_buffer, &recv_buffer, 1, MPI_INT, MPI_SUM, 0, comm_global);
509 return recv_buffer;
510}
511
512void CommMPI::GlobalSumArray(vmg_int* array, const vmg_int& size)
513{
514 MPI_Allreduce(MPI_IN_PLACE, array, size, MPI_INT, MPI_SUM, comm_global);
515}
516
517vmg_float CommMPI::LevelSum(const Grid& grid, const vmg_float& value)
518{
519 vmg_float result = value;
520 MPI_Comm comm = settings.CommunicatorLocal(grid);
521 assert(comm != MPI_COMM_NULL);
522 MPI_Allreduce(MPI_IN_PLACE, &result, 1, MPI_DOUBLE, MPI_SUM, comm);
523 return result;
524}
525
526vmg_float CommMPI::LevelSumRoot(const Grid& grid, const vmg_float& value)
527{
528 vmg_float recv_buffer = value;
529 vmg_float send_buffer = value;
530 MPI_Comm comm = settings.CommunicatorLocal(grid);
531 assert(comm != MPI_COMM_NULL);
532 MPI_Reduce(&send_buffer, &recv_buffer, 1, MPI_DOUBLE, MPI_SUM, 0, comm);
533 return recv_buffer;
534}
535
536void CommMPI::LevelSumArray(const Grid& grid, vmg_float* array, const vmg_int& size)
537{
538 MPI_Comm comm = settings.CommunicatorLocal(grid);
539 assert(comm != MPI_COMM_NULL);
540 MPI_Allreduce(MPI_IN_PLACE, array, size, MPI_DOUBLE, MPI_SUM, comm);
541}
542
543vmg_int CommMPI::LevelSum(const Grid& grid, const vmg_int& value)
544{
545 vmg_int result = value;
546 MPI_Comm comm = settings.CommunicatorLocal(grid);
547 assert(comm != MPI_COMM_NULL);
548 MPI_Allreduce(MPI_IN_PLACE, &result, 1, MPI_INT, MPI_SUM, comm);
549 return result;
550}
551
552vmg_int CommMPI::LevelSumRoot(const Grid& grid, const vmg_int& value)
553{
554 vmg_int recv_buffer = value;
555 vmg_int send_buffer = value;
556 MPI_Comm comm = settings.CommunicatorLocal(grid);
557 assert(comm != MPI_COMM_NULL);
558 MPI_Reduce(&send_buffer, &recv_buffer, 1, MPI_INT, MPI_SUM, 0, comm);
559 return recv_buffer;
560}
561
562void CommMPI::LevelSumArray(const Grid& grid, vmg_int* array, const vmg_int& size)
563{
564 MPI_Comm comm = settings.CommunicatorLocal(grid);
565 assert(comm != MPI_COMM_NULL);
566 MPI_Allreduce(MPI_IN_PLACE, array, size, MPI_INT, MPI_SUM, comm);
567}
568
569void CommMPI::PrintString(const char* format, ...)
570{
571 va_list args;
572 va_start(args, format);
573 vsprintf(print_buffer, format, args);
574 printf("VMG: Rank %d: %s\n", GlobalRank(), print_buffer);
575 va_end(args);
576}
577
578void CommMPI::PrintStringOnce(const char* format, ...)
579{
580 if (GlobalRank() == 0) {
581 va_list args;
582 va_start(args, format);
583 vsprintf(print_buffer, format, args);
584 printf("VMG: Rank %d: %s\n", GlobalRank(), print_buffer);
585 va_end(args);
586 }
587}
588
589void CommMPI::PrintXML(const std::string& filename, const std::string& xml_data)
590{
591 MPI_File file;
592 std::stringstream path, xml_header;
593
594 pugi::xml_document doc;
595 pugi::xml_node node_data = doc.append_child("Global").append_child("NumProcesses").append_child(pugi::node_pcdata);
596 node_data.set_value(Helper::ToString(GlobalProcs().Product()).c_str());
597 doc.save(xml_header);
598
599 path << OutputPath() << filename;
600
601 char* filename_array = Helper::GetCharArray(path.str());
602 char* xml_header_array = Helper::GetCharArray(xml_header.str());
603 char* str_array = Helper::GetCharArray(xml_data);
604
605 MPI_File_open(MPI_COMM_SELF, filename_array, MPI_MODE_WRONLY|MPI_MODE_CREATE, MPI_INFO_NULL, &file);
606 MPI_File_set_size(file, 0);
607 MPI_File_write(file, xml_header_array, xml_header.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
608 MPI_File_write(file, str_array, xml_data.size(), MPI_CHAR, MPI_STATUS_IGNORE);
609 MPI_File_close(&file);
610
611 delete [] filename_array;
612 delete [] xml_header_array;
613 delete [] str_array;
614}
615
616void CommMPI::PrintXMLAll(const std::string& filename, const std::string& xml_data)
617{
618 MPI_File file;
619 std::stringstream path;
620
621 path << OutputPath() << filename;
622
623 char* filename_array = Helper::GetCharArray(path.str());
624 char* str_array = Helper::GetCharArray(xml_data);
625
626 MPI_File_open(comm_global, filename_array, MPI_MODE_WRONLY|MPI_MODE_CREATE, MPI_INFO_NULL, &file);
627 MPI_File_set_size(file, 0);
628
629 if (GlobalRank() == 0) {
630 std::stringstream xml_header;
631 pugi::xml_document doc;
632 pugi::xml_node node_data = doc.append_child("Global").append_child("NumProcesses").append_child(pugi::node_pcdata);
633 node_data.set_value(Helper::ToString(GlobalProcs().Product()).c_str());
634 doc.save(xml_header);
635
636 char* xml_header_array = Helper::GetCharArray(xml_header.str());
637
638 MPI_File_write_shared(file, xml_header_array, xml_header.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
639
640 delete [] xml_header_array;
641 }
642
643 MPI_File_write_ordered(file, str_array, xml_data.size(), MPI_CHAR, MPI_STATUS_IGNORE);
644 MPI_File_close(&file);
645
646 delete [] filename_array;
647 delete [] str_array;
648}
649
650void CommMPI::PrintGridInformation(const Grid& grid, char* filename, const std::string& name)
651{
652 std::stringstream buf;
653 MPI_File file;
654 int rank, size;
655 int size_local, size_local_max;
656
657 MPI_Comm comm = settings.CommunicatorGlobal(grid);
658 MPI_Comm comm_local = settings.CommunicatorLocal(grid);
659
660 if (comm_local != MPI_COMM_NULL)
661 MPI_Comm_size(comm_local, &size_local);
662 else
663 size_local = 0;
664
665 if (comm != MPI_COMM_NULL) {
666
667 MPI_Reduce(&size_local, &size_local_max, 1, MPI_INT, MPI_MAX, 0, comm);
668
669 MPI_File_open(comm, filename, MPI_MODE_WRONLY|MPI_MODE_CREATE|MPI_MODE_APPEND, MPI_INFO_NULL, &file);
670
671 MPI_Comm_rank(comm, &rank);
672 MPI_Comm_size(comm, &size);
673
674 if (rank == 0) {
675
676 buf << "###########################################################" << std::endl
677 << "GLOBAL INFORMATION:" << std::endl
678 << " NAME: " << name << std::endl
679 << " LEVEL: " << grid.Level() << std::endl
680 << " COMM_SIZE_GLOBAL: " << size << std::endl
681 << " COMM_SIZE_LOCAL: " << size_local_max << std::endl
682 << " GLOBAL:" << std::endl
683 << " GLOBAL_FINER_BEGIN: " << grid.Global().GlobalFinerBegin() << std::endl
684 << " GLOBAL_FINER_END: " << grid.Global().GlobalFinerEnd() << std::endl
685 << " GLOBAL_FINER_SIZE: " << grid.Global().GlobalFinerSize() << std::endl
686 << " FINEST_ABS_BEGIN: " << grid.Global().FinestAbsBegin() << std::endl
687 << " FINEST_ABS_END: " << grid.Global().FinestAbsEnd() << std::endl
688 << " FINEST_ABS_SIZE: " << grid.Global().FinestAbsSize() << std::endl
689 << " GLOBAL_SIZE: " << grid.Global().GlobalSize() << std::endl
690 << " EXTENT:" << std::endl
691 << " BEGIN: " << grid.Extent().Begin() << std::endl
692 << " END: " << grid.Extent().End() << std::endl
693 << " SIZE: " << grid.Extent().Size() << std::endl
694 << " MESH_WIDTH: " << grid.Extent().MeshWidth() << std::endl
695 << std::endl
696 << "LOCAL INFORMATION:" << std::endl;
697 }
698
699 buf << "RANK " << rank << ":" << std::endl
700 << " GLOBAL:" << std::endl
701 << " LOCAL_BEGIN: " << grid.Global().LocalBegin() << std::endl
702 << " LOCAL_END: " << grid.Global().LocalEnd() << std::endl
703 << " LOCAL_SIZE: " << grid.Global().LocalSize() << std::endl
704 << " LOCAL_FINER_BEGIN: " << grid.Global().LocalFinerBegin() << std::endl
705 << " LOCAL_FINER_END: " << grid.Global().LocalFinerEnd() << std::endl
706 << " LOCAL_FINER_SIZE: " << grid.Global().LocalFinerSize() << std::endl
707 << " BOUNDARY_TYPE: " << grid.Global().BoundaryType() << std::endl
708 << " LOCAL:" << std::endl
709 << " BEGIN: " << grid.Local().Begin() << std::endl
710 << " END: " << grid.Local().End() << std::endl
711 << " SIZE: " << grid.Local().Size() << std::endl
712 << " SIZE_TOTAL: " << grid.Local().SizeTotal() << std::endl
713 << " HALO_BEGIN_1: " << grid.Local().HaloBegin1() << std::endl
714 << " HALO_END_1: " << grid.Local().HaloEnd1() << std::endl
715 << " HALO_SIZE_1: " << grid.Local().HaloSize1() << std::endl
716 << " HALO_BEGIN_2: " << grid.Local().HaloBegin2() << std::endl
717 << " HALO_END_2: " << grid.Local().HaloEnd2() << std::endl
718 << " HALO_SIZE_2: " << grid.Local().HaloSize2() << std::endl
719 << " BOUNDARY_BEGIN_1: " << grid.Local().BoundaryBegin1() << std::endl
720 << " BOUNDARY_END_1: " << grid.Local().BoundaryEnd1() << std::endl
721 << " BOUNDARY_SIZE_1: " << grid.Local().BoundarySize1() << std::endl
722 << " BOUNDARY_BEGIN_2: " << grid.Local().BoundaryBegin2() << std::endl
723 << " BOUNDARY_END_2: " << grid.Local().BoundaryEnd2() << std::endl
724 << " BOUNDARY_SIZE_2: " << grid.Local().BoundarySize2() << std::endl
725 << " FINER_BEGIN: " << grid.Local().FinerBegin() << std::endl
726 << " FINER_END: " << grid.Local().FinerEnd() << std::endl
727 << " FINER_SIZE: " << grid.Local().FinerSize() << std::endl;
728
729 if (rank == size-1)
730 buf << "###########################################################" << std::endl;
731
732 char* char_buf = Helper::GetCharArray(buf.str());
733 MPI_File_write_ordered(file, char_buf, buf.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
734 delete [] char_buf;
735
736 MPI_File_close(&file);
737
738 }
739}
740
741void CommMPI::PrintAllSettings()
742{
743 std::stringstream buf;
744 MPI_File file;
745
746 const Multigrid& mg = *MG::GetSol();
747
748 buf << OutputPath() << "settings.txt";
749 char *filename = Helper::GetCharArray(buf.str());
750
751 MPI_File_open(comm_global, filename, MPI_MODE_WRONLY|MPI_MODE_CREATE, MPI_INFO_NULL, &file);
752 MPI_File_set_size(file, 0);
753 MPI_File_close(&file);
754
755 for (int i=mg.MinLevel(); i<=mg.MaxLevel(); ++i)
756 PrintGridInformation(mg(i), filename, "MULTIGRID");
757
758 for (int i=mg.MinLevel()+1; i<=mg.MaxLevel(); ++i)
759 PrintGridInformation(settings.CoarserGrid(mg(i)), filename, "COARSER_GRID");
760
761 for (int i=mg.MinLevel(); i<mg.MaxLevel(); ++i)
762 PrintGridInformation(settings.FinerGrid(mg(i)), filename, "FINER_GRID");
763
764 delete [] filename;
765
766}
767
768void CommMPI::PrintGrid(Grid& grid, const char* information)
769{
770 int output_count = OutputCount();
771
772#ifdef HAVE_VTK
773
774 if (settings.CommunicatorLocal(grid) != MPI_COMM_NULL) {
775
776 Index end, end_global;
777
778 for (int i=0; i<3; ++i) {
779 end[i] = grid.Local().End()[i];
780 end_global[i] = grid.Global().LocalEnd()[i];
781 }
782
783 vtkSmartPointer<vtkImageData> image = vtkSmartPointer<vtkImageData>::New();
784 image->SetExtent(grid.Global().LocalBegin().X(), end_global.X()-1,
785 grid.Global().LocalBegin().Y(), end_global.Y()-1,
786 grid.Global().LocalBegin().Z(), end_global.Z()-1);
787 image->SetSpacing(grid.Extent().MeshWidth().vec());
788 image->SetOrigin(grid.Extent().Begin().vec());
789 image->SetScalarTypeToDouble();
790 image->SetNumberOfScalarComponents(1);
791 image->AllocateScalars();
792 image->GetPointData()->GetScalars()->SetName(information);
793
794 Index i;
795 for (i.X()=grid.Local().Begin().X(); i.X()<end.X(); ++i.X())
796 for (i.Y()=grid.Local().Begin().Y(); i.Y()<end.Y(); ++i.Y())
797 for (i.Z()=grid.Local().Begin().Z(); i.Z()<end.Z(); ++i.Z())
798 image->SetScalarComponentFromDouble(i.X() - grid.Local().Begin().X() + grid.Global().LocalBegin().X(),
799 i.Y() - grid.Local().Begin().Y() + grid.Global().LocalBegin().Y(),
800 i.Z() - grid.Local().Begin().Z() + grid.Global().LocalBegin().Z(),
801 0, grid.GetVal(i));
802
803 image->Update();
804
805 int rank, size;
806 MPI_Comm_rank(comm_global, &rank);
807 MPI_Comm_size(comm_global, &size);
808
809 std::stringstream filename;
810 filename << OutputPath() << std::setw(4) << std::setfill('0') << output_count << "_" << rank << ".vti";
811
812 vtkSmartPointer<vtkXMLImageDataWriter> writer = vtkSmartPointer<vtkXMLImageDataWriter>::New();
813 writer->SetFileName(filename.str().c_str());
814 writer->SetCompressorTypeToNone();
815 writer->SetDataModeToAscii();
816 writer->SetInput(image);
817 writer->Update();
818 writer->Write();
819
820 }
821
822#else /* HAVE_VTK */
823 Index i;
824 std::stringstream buf;
825
826 Index begin, end;
827 Index begin_local, end_local, begin_global, end_global;
828
829 CommToGhosts(grid);
830
831 for (int i=0; i<3; ++i) {
832 end[i] = grid.Local().End()[i] + (grid.Global().LocalEnd()[i] == grid.Global().GlobalSize()[i] ? 0 : grid.Local().HaloSize1()[i]);
833 end_local[i] = grid.Global().LocalEnd()[i] - (grid.Global().LocalEnd()[i] == grid.Global().GlobalSize()[i] ? 1 : 0);
834 }
835
836 begin = grid.Local().Begin();
837 begin_local = grid.Global().LocalBegin();
838 begin_global = 0;
839 end_global = grid.Global().GlobalSize()-1;
840
841 for (i.Z()=begin.Z(); i.Z()<end.Z(); ++i.Z())
842 for (i.Y()=begin.Y(); i.Y()<end.Y(); ++i.Y())
843 for (i.X()=begin.X(); i.X()<end.X(); ++i.X())
844 buf << std::scientific << grid.GetVal(i) << " ";
845
846 CreateOutputFiles(grid, buf, information,
847 begin_global, end_global,
848 begin_local, end_local,
849 output_count);
850#endif /* HAVE_VTK */
851}
852
853void CommMPI::PrintDefect(Grid& sol, Grid& rhs, const char* information)
854{
855 TempGrid *temp = MG::GetTempGrid();
856 temp->SetProperties(sol);
857 temp->ImportFromResidual(sol, rhs);
858 PrintGrid(*temp, information);
859}
860
861void CommMPI::CreateOutputFiles(const Grid& grid, const std::stringstream& serial_data, const char* information,
862 const Index& begin_global, const Index& end_global,
863 const Index& begin_local, const Index& end_local,
864 const int& output_count)
865{
866 MPI_Comm comm = settings.CommunicatorGlobal(grid);
867
868 if (comm != MPI_COMM_NULL) {
869
870 MPI_File file;
871 std::string conv_information = Helper::ReplaceWhitespaces(information, "_");
872
873 CreateParallelOutputFile(grid, comm, output_count, conv_information.c_str(),
874 begin_global, end_global, begin_local, end_local);
875
876 file = CreateSerialOutputFile(grid, comm, output_count, conv_information.c_str(),
877 begin_global, end_global, begin_local, end_local);
878
879 char *char_buf = Helper::GetCharArray(serial_data.str());
880 MPI_File_write(file, char_buf, serial_data.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
881 delete [] char_buf;
882
883 FinalizeSerialOutputFile(file);
884
885 }
886}
887
888void CommMPI::CreateParallelOutputFile(const Grid& grid, MPI_Comm& comm,
889 const int& output_count, const char* information,
890 const Index& begin_global, const Index& end_global,
891 const Index& begin_local, const Index& end_local)
892{
893 int rank;
894 MPI_File file;
895 char parallel_filename[513], serial_filename[513];
896 std::stringstream buf;
897
898 MPI_Comm_rank(comm, &rank);
899
900 sprintf(parallel_filename, "%s%04d.pvti", OutputPath().c_str(), output_count);
901 sprintf(serial_filename, "%04d_%d.vti", output_count, rank);
902
903 MPI_File_open(comm, parallel_filename, MPI_MODE_WRONLY|MPI_MODE_CREATE, MPI_INFO_NULL, &file);
904 MPI_File_set_size(file, 0);
905
906 if (rank == 0) {
907
908 buf << "<?xml version=\"1.0\"?>" << std::endl
909 << "<VTKFile type=\"PImageData\" version=\"0.1\" byte_order=\"LittleEndian\">" << std::endl
910 << " <PImageData WholeExtent=\"";
911
912 for (int i=0; i<3; ++i)
913 buf << begin_global[i] << " " << end_global[i] << " ";
914
915 buf << "\"" << std::endl
916 << " GhostLevel=\"0\" Origin=\"0 0 0\" Spacing=\"";
917
918 for (int i=0; i<3; ++i)
919 buf << grid.Extent().MeshWidth()[i] << " ";
920
921 buf << "\">" << std::endl
922 << " <PPointData Scalars=\"" << information << "\">" << std::endl
923 << " <PDataArray type=\"Float32\" Name=\"" << information << "\"/>" << std::endl
924 << " </PPointData>" << std::endl;
925
926 char* char_buf = Helper::GetCharArray(buf.str());
927
928 MPI_File_write_shared(file, char_buf, buf.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
929
930 delete [] char_buf;
931 }
932
933 buf.str("");
934
935 if ((end_local-begin_local).Product() > 0) {
936 buf << " <Piece Extent=\"";
937
938 for (int i=0; i<3; ++i)
939 buf << begin_local[i] << " " << end_local[i] << " ";
940
941 buf << "\" Source=\"" << serial_filename << "\"/>" << std::endl;
942 }
943
944 char* char_buf = Helper::GetCharArray(buf.str());
945
946 MPI_File_write_ordered(file, char_buf, buf.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
947
948 delete [] char_buf;
949
950 if (rank == 0) {
951
952 buf.str("");
953
954 buf << " </PImageData>" << std::endl
955 << "</VTKFile>" << std::endl;
956
957 char* char_buf = Helper::GetCharArray(buf.str());
958
959 MPI_File_write_shared(file, char_buf, buf.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
960
961 delete [] char_buf;
962 }
963
964 MPI_File_close(&file);
965}
966
967MPI_File CommMPI::CreateSerialOutputFile(const Grid& grid, MPI_Comm& comm,
968 const int& output_count, const char* information,
969 const Index& begin_global, const Index& end_global,
970 const Index& begin_local, const Index& end_local)
971{
972 char serial_filename[513];
973 int rank;
974 MPI_File file;
975 std::stringstream buf;
976
977 MPI_Comm_rank(comm_global, &rank);
978
979 sprintf(serial_filename, "%s%04d_%d.vti", OutputPath().c_str(), output_count, rank);
980
981 MPI_File_open(MPI_COMM_SELF, serial_filename, MPI_MODE_WRONLY|MPI_MODE_CREATE, MPI_INFO_NULL, &file);
982
983 buf << "<?xml version=\"1.0\"?>" << std::endl
984 << "<VTKFile type=\"ImageData\" version=\"0.1\" byte_order=\"LittleEndian\">" << std::endl
985 << " <ImageData WholeExtent=\"";
986
987 for (int i=0; i<3; ++i)
988 buf << begin_global[i] << " " << end_global[i] << " ";
989
990 buf << "\"" << std::endl
991 << " Origin=\"0 0 0\" Spacing=\"";
992
993 for (int i=0; i<3; ++i)
994 buf << grid.Extent().MeshWidth()[i] << " ";
995
996 buf << "\">" << std::endl
997 << " <Piece Extent=\"";
998
999 for (int i=0; i<3; ++i)
1000 buf << begin_local[i] << " " << end_local[i] << " ";
1001
1002 buf << "\">" << std::endl
1003 << " <PointData Scalars=\"" << information << "\">" << std::endl
1004 << " <DataArray type=\"Float32\" Name=\"" << information << "\" format=\"ascii\">" << std::endl
1005 << " ";
1006
1007 char* char_buf = Helper::GetCharArray(buf.str());
1008 MPI_File_write(file, char_buf, buf.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
1009 delete [] char_buf;
1010
1011 return file;
1012}
1013
1014void CommMPI::FinalizeSerialOutputFile(MPI_File& file)
1015{
1016 std::stringstream buf;
1017
1018 buf << std::endl
1019 << " </DataArray>" << std::endl
1020 << " </PointData>" << std::endl
1021 << " </Piece>" << std::endl
1022 << " </ImageData>" << std::endl
1023 << "</VTKFile>" << std::endl;
1024
1025 char* char_buf = Helper::GetCharArray(buf.str());
1026 MPI_File_write(file, char_buf, buf.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
1027 delete [] char_buf;
1028
1029 MPI_File_close(&file);
1030}
1031
1032int CommMPI::GlobalRank() const
1033{
1034 int rank;
1035 MPI_Comm_rank(comm_global, &rank);
1036 return rank;
1037}
1038
1039int CommMPI::GlobalSize() const
1040{
1041 int size;
1042 MPI_Comm_size(comm_global, &size);
1043 return size;
1044}
1045
1046Index CommMPI::GlobalPos() const
1047{
1048 Index dims, periods, coords;
1049 MPI_Cart_get(comm_global, 3, dims.vec(), periods.vec(), coords.vec());
1050 return coords;
1051}
1052
1053Index CommMPI::GlobalProcs() const
1054{
1055 Index dims, periods, coords;
1056 MPI_Cart_get(comm_global, 3, dims.vec(), periods.vec(), coords.vec());
1057 return dims;
1058}
1059
1060int CommMPI::Rank(const Grid& grid) const
1061{
1062 int rank;
1063 MPI_Comm comm = settings.CommunicatorLocal(grid);
1064 assert(comm != MPI_COMM_NULL);
1065 MPI_Comm_rank(comm, &rank);
1066 return rank;
1067}
1068
1069int CommMPI::Size(const Grid& grid) const
1070{
1071 int size;
1072 MPI_Comm comm = settings.CommunicatorLocal(grid);
1073 assert(comm != MPI_COMM_NULL);
1074 MPI_Comm_size(comm, &size);
1075 return size;
1076}
1077
1078Index CommMPI::Pos(const Grid& grid) const
1079{
1080 Index dims, periods, coords;
1081 MPI_Comm comm = settings.CommunicatorLocal(grid);
1082 assert(comm != MPI_COMM_NULL);
1083 MPI_Cart_get(comm, 3, dims.vec(), periods.vec(), coords.vec());
1084 return coords;
1085}
1086
1087Index CommMPI::Procs(const Grid& grid) const
1088{
1089 Index dims, periods, coords;
1090 MPI_Comm comm = settings.CommunicatorLocal(grid);
1091 assert(comm != MPI_COMM_NULL);
1092 MPI_Cart_get(comm, 3, dims.vec(), periods.vec(), coords.vec());
1093 return dims;
1094}
1095
1096void CommMPI::InitCommMPI(const MPI_Comm& comm)
1097{
1098 int status, size, rank;
1099 int dims[3] = {0, 0, 0};
1100 int periods[3];
1101
1102 for (int i=0; i<3; ++i)
1103 periods[i] = (BoundaryConditions()[i] == Periodic ? 1 : 0);
1104
1105 MPI_Comm_size(comm, &size);
1106 MPI_Comm_rank(comm, &rank);
1107
1108 MPI_Topo_test(comm, &status);
1109
1110 if (status == MPI_CART) {
1111
1112 comm_global = comm;
1113
1114 }else {
1115
1116 MPI_Dims_create(size, 3, dims);
1117
1118 MPI_Cart_create(comm, 3, dims, periods, 1, &comm_global);
1119
1120 }
1121
1122 MPI_Info_create(&info);
1123 char key[] = "no_locks";
1124 char val[] = "true";
1125 MPI_Info_set(info, key, val);
1126
1127}
1128
1129CommMPI::~CommMPI()
1130{
1131 MPI_Comm_free(&comm_global);
1132#ifdef VMG_ONE_SIDED
1133 if (win_created)
1134 MPI_Win_free(&win);
1135#endif
1136 MPI_Info_free(&info);
1137}
1138
1139Grid& CommMPI::GetCoarserGrid(Multigrid& multigrid)
1140{
1141 return settings.CoarserGrid(multigrid(multigrid.Level()));
1142}
1143
1144Grid& CommMPI::GetFinerGrid(Multigrid& multigrid)
1145{
1146 return settings.FinerGrid(multigrid(multigrid.Level()));
1147}
1148
1149std::string CommMPI::CreateOutputDirectory()
1150{
1151#ifdef HAVE_BOOST_FILESYSTEM
1152 std::string path, unique_path;
1153 std::stringstream unique_suffix;
1154 int suffix_counter = 0;
1155 char buffer[129];
1156 time_t rawtime;
1157 struct tm *timeinfo;
1158 int path_size;
1159 char* path_array;
1160
1161 if (GlobalRank() == 0) {
1162
1163 time(&rawtime);
1164 timeinfo = localtime(&rawtime);
1165 strftime(buffer, 128, "./output/%Y_%m_%d_%H_%M_%S/", timeinfo);
1166 path = buffer;
1167
1168 if (!fs::exists("output"))
1169 fs::create_directory("output");
1170
1171 unique_path = path;
1172
1173 while (fs::exists(unique_path.c_str())) {
1174
1175 unique_suffix.str("");
1176 unique_suffix << "_" << suffix_counter++ << "/";
1177
1178 unique_path = path;
1179 unique_path.replace(unique_path.size()-1, 1, unique_suffix.str());
1180
1181 }
1182
1183 fs::create_directory(unique_path.c_str());
1184
1185 path_size = unique_path.size() + 1;
1186 path_array = Helper::GetCharArray(unique_path);
1187
1188 MPI_Bcast(&path_size, 1, MPI_INT, 0, comm_global);
1189
1190 }else {
1191
1192 MPI_Bcast(&path_size, 1, MPI_INT, 0, comm_global);
1193 path_array = new char[path_size];
1194
1195 }
1196
1197 MPI_Bcast(path_array, path_size, MPI_CHAR, 0, comm_global);
1198
1199 unique_path = path_array;
1200
1201 delete [] path_array;
1202
1203 return unique_path;
1204
1205#else
1206
1207 return "./";
1208
1209#endif
1210}
1211
1212
1213#endif /* HAVE_MPI */
Note: See TracBrowser for help on using the repository browser.