source: src/comm/comm_mpi.cpp@ 1610dc

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

Did some work on the MPI communication.

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

  • Property mode set to 100644
File size: 39.0 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
18#ifdef HAVE_BOOST_FILESYSTEM
19#include <boost/filesystem.hpp>
20namespace fs = boost::filesystem;
21#endif
22
23#include <cstring>
24#include <sstream>
25
26#include "base/helper.hpp"
27#include "base/linked_cell_list.hpp"
28#include "base/timer.hpp"
29#include "comm/comm_mpi.hpp"
30#include "grid/grid.hpp"
31#include "grid/multigrid.hpp"
32#include "grid/tempgrid.hpp"
33#include "mg.hpp"
34
35static char print_buffer[512];
36
37using namespace VMG;
38
39void CommMPI::CommSubgrid(Grid& grid_old, Grid& grid_new)
40{
41 Timer::Start("CommSubgrid");
42
43 MPI_Comm comm = comm_info.GetUnionCommunicator(grid_old, grid_new);
44
45 if (comm != MPI_COMM_NULL) {
46
47 int comm_rank, comm_size;
48 MPI_Comm_rank(comm, &comm_rank);
49 MPI_Comm_size(comm, &comm_size);
50
51 Index begin_local, end_local;
52 Index begin_remote, end_remote;
53 Index sizes, subsizes, starts;
54 Index offset_new, offset_old;
55 int* buffer = new int[6*comm_size];
56
57 for (int i=0; i<3; ++i) {
58
59 offset_old[i] = (grid_old.Local().HaloEnd1()[i] > 0 ? grid_old.Local().Begin()[i] : 0);
60 offset_new[i] = (grid_new.Local().HaloEnd1()[i] > 0 ? grid_new.Local().Begin()[i] : 0);
61
62 buffer[6*comm_rank + i ] = grid_new.Global().BeginLocal()[i];
63 buffer[6*comm_rank + i+3] = grid_new.Global().EndLocal()[i];
64
65 }
66
67 //TODO: Create LUT to avoid global communication here
68 MPI_Allgather(MPI_IN_PLACE, 6, MPI_INT, buffer, 6, MPI_INT, comm);
69
70 for (int i=0; i<comm_size; ++i) {
71
72 for (int j=0; j<3; ++j) {
73
74 begin_remote[j] = buffer[6*i + j];
75 end_remote[j] = buffer[6*i + j+3];
76
77 }
78
79 begin_remote = begin_remote.Clamp(grid_old.Global().BeginLocal(), grid_old.Global().EndLocal());
80 end_remote = end_remote.Clamp(grid_old.Global().BeginLocal(), grid_old.Global().EndLocal());
81
82 if ((end_remote - begin_remote).Product() > 0) {
83
84 MPI_Datatype dt_temp;
85
86 sizes = grid_old.Local().SizeTotal();
87 subsizes = end_remote - begin_remote;
88 starts = begin_remote - grid_old.Global().BeginLocal() + offset_old;
89
90 MPI_Type_create_subarray(3, sizes.vec(), subsizes.vec(), starts.vec(), MPI_ORDER_C, MPI_DOUBLE, &dt_temp);
91
92 MPI_Type_commit(&dt_temp);
93
94 MPI_Isend(&grid_old(0), 1, dt_temp, i, 0, comm, &Request());
95
96 MPI_Type_free(&dt_temp);
97
98 }
99
100 }
101
102 for (int i=0; i<3; ++i) {
103
104 buffer[6*comm_rank + i ] = grid_old.Global().BeginLocal()[i];
105 buffer[6*comm_rank + i+3] = grid_old.Global().EndLocal()[i];
106
107 }
108
109 //TODO: Create LUT to avoid global communication here
110 MPI_Allgather(MPI_IN_PLACE, 6, MPI_INT, buffer, 6, MPI_INT, comm);
111
112 for (int i=0; i<comm_size; ++i) {
113
114 for (int j=0; j<3; ++j) {
115
116 begin_remote[j] = buffer[6*i + j];
117 end_remote[j] = buffer[6*i + j+3];
118
119 }
120
121 begin_remote = begin_remote.Clamp(grid_new.Global().BeginLocal(), grid_new.Global().EndLocal());
122 end_remote = end_remote.Clamp(grid_new.Global().BeginLocal(), grid_new.Global().EndLocal());
123
124 if ((end_remote - begin_remote).Product() > 0) {
125
126 MPI_Datatype dt_temp;
127
128 sizes = grid_new.Local().SizeTotal();
129 subsizes = end_remote - begin_remote;
130 starts = begin_remote - grid_new.Global().BeginLocal() + offset_new;
131
132 MPI_Type_create_subarray(3, sizes.vec(), subsizes.vec(), starts.vec(), MPI_ORDER_C, MPI_DOUBLE, &dt_temp);
133
134 MPI_Type_commit(&dt_temp);
135
136 MPI_Irecv(&grid_new(0), 1, dt_temp, i, 0, comm, &Request());
137
138 MPI_Type_free(&dt_temp);
139
140 }
141
142 }
143
144 ClearRequests();
145
146 delete [] buffer;
147
148 }
149
150 CommToGhosts(grid_new);
151
152 Timer::Stop("CommSubgrid");
153}
154
155void CommMPI::CommAddSubgrid(Grid& grid_old, Grid& grid_new)
156{
157 Timer::Start("CommAddSubgrid");
158
159 MPI_Comm comm = comm_info.GetUnionCommunicator(grid_old, grid_new);
160
161 if (comm != MPI_COMM_NULL) {
162
163 int comm_rank, comm_size;
164 MPI_Comm_rank(comm, &comm_rank);
165 MPI_Comm_size(comm, &comm_size);
166
167 Index begin_local, end_local;
168 Index begin_remote, end_remote;
169 Index sizes, subsizes, starts;
170 int* buffer = new int[6*comm_size];
171
172 for (int i=0; i<3; ++i) {
173
174 buffer[6*comm_rank + i ] = grid_new.Global().BeginLocal()[i];
175 buffer[6*comm_rank + i+3] = grid_new.Global().EndLocal()[i];
176
177 }
178
179 //TODO: Create LUT to avoid global communication here
180 MPI_Allgather(MPI_IN_PLACE, 6, MPI_INT, buffer, 6, MPI_INT, comm);
181
182 for (int i=0; i<comm_size; ++i) {
183 for (int j=0; j<3; ++j) {
184
185 begin_remote[j] = buffer[6*i + j];
186 end_remote[j] = buffer[6*i + j+3];
187
188 }
189
190 begin_remote = begin_remote.Clamp(grid_old.Global().BeginLocal(), grid_old.Global().EndLocal());
191 end_remote = end_remote.Clamp(grid_old.Global().BeginLocal(), grid_old.Global().EndLocal());
192
193 if ((end_remote - begin_remote).Product() > 0) {
194
195 MPI_Datatype dt_temp;
196
197 sizes = grid_old.Local().SizeTotal();
198 subsizes = end_remote - begin_remote;
199 starts = begin_remote - grid_old.Global().BeginLocal();
200
201 MPI_Type_create_subarray(3, sizes.vec(), subsizes.vec(), starts.vec(), MPI_ORDER_C, MPI_DOUBLE, &dt_temp);
202
203 MPI_Type_commit(&dt_temp);
204
205 MPI_Isend(&grid_old(0), 1, dt_temp, i, 0, comm, &Request());
206
207 MPI_Type_free(&dt_temp);
208
209 }
210
211 }
212
213 for (int i=0; i<3; ++i) {
214
215 buffer[6*comm_rank + i ] = grid_old.Global().BeginLocal()[i];
216 buffer[6*comm_rank + i+3] = grid_old.Global().EndLocal()[i];
217
218 }
219
220 //TODO: Create LUT to avoid global communication here
221 MPI_Allgather(MPI_IN_PLACE, 6, MPI_INT, buffer, 6, MPI_INT, comm);
222
223 for (int i=0; i<comm_size; ++i) {
224
225 for (int j=0; j<3; ++j) {
226
227 begin_remote[j] = buffer[6*i + j];
228 end_remote[j] = buffer[6*i + j+3];
229
230 }
231
232 begin_remote = begin_remote.Clamp(grid_new.Global().BeginLocal(), grid_new.Global().EndLocal());
233 end_remote = end_remote.Clamp(grid_new.Global().BeginLocal(), grid_new.Global().EndLocal());
234
235 if ((end_remote - begin_remote).Product() > 0) {
236
237 MPI_Datatype dt_temp;
238
239 sizes = grid_new.Local().SizeTotal();
240 subsizes = end_remote - begin_remote;
241 starts = begin_remote - grid_old.Global().BeginLocal();
242
243 MPI_Type_create_subarray(3, sizes.vec(), subsizes.vec(), starts.vec(), MPI_ORDER_C, MPI_DOUBLE, &dt_temp);
244
245 MPI_Type_commit(&dt_temp);
246
247 ReceiveAndAddSubgrid(grid_new, comm, dt_temp, i, 0);
248
249 MPI_Type_free(&dt_temp);
250
251 }
252
253 }
254
255 ClearRequests();
256
257 delete [] buffer;
258
259 }
260
261 CommToGhosts(grid_new);
262
263 Timer::Stop("CommAddSubgrid");
264}
265
266void CommMPI::CommFromGhosts(Grid& grid)
267{
268 Timer::Start("CommFromGhosts");
269
270 MPI_Comm comm = comm_info.GetCommunicator(grid);
271
272 if (comm != MPI_COMM_NULL) {
273
274 const Index has_halo_1 = grid.Local().HasHalo1();
275 const Index has_halo_2 = grid.Local().HasHalo2();
276 Tuple<int> neighbors;
277
278 for (int i=0; i<3; ++i) {
279
280 MPI_Cart_shift(comm, i, 1, &neighbors.Left(), &neighbors.Right());
281
282 if (has_halo_1[i]) {
283 MPI_Datatype dts1 = comm_info.GetDatatypeSubarray(grid, grid.Iterators().Halo1()[i]);
284 MPI_Isend(&grid(0), 1, dts1, neighbors.Left(), 1, comm, &Request());
285 }
286
287 if (has_halo_2[i]) {
288
289 MPI_Datatype dts2 = comm_info.GetDatatypeSubarray(grid, grid.Iterators().Halo2()[i]);
290 MPI_Isend(&grid(0), 1, dts2, neighbors.Right(), 2, comm, &Request());
291
292 MPI_Datatype dtr1 = comm_info.GetDatatypeSubarray(grid, grid.Iterators().NearBoundary2()[i]);
293 ReceiveAndAddSubgrid(grid, comm, dtr1, neighbors.Right(), 1);
294
295 }
296
297 if (has_halo_1[i]) {
298 MPI_Datatype dtr2 = comm_info.GetDatatypeSubarray(grid, grid.Iterators().NearBoundary1()[i]);
299 ReceiveAndAddSubgrid(grid, comm, dtr2, neighbors.Left(), 2);
300 }
301
302 }
303
304 ClearRequests();
305
306 }
307
308 Timer::Stop("CommFromGhosts");
309}
310
311void CommMPI::CommToGhosts(Grid& grid)
312{
313 Timer::Start("CommToGhosts");
314
315 MPI_Comm comm = comm_info.GetCommunicator(grid);
316
317 if (comm != MPI_COMM_NULL) {
318
319 const Index has_halo_1 = grid.Local().HasHalo1();
320 const Index has_halo_2 = grid.Local().HasHalo2();
321 Tuple<int> neighbors;
322
323 for (int i=2; i>=0; --i) {
324
325 MPI_Cart_shift(comm, i, 1, &neighbors.Left(), &neighbors.Right());
326
327 MPI_Datatype dts_left, dts_right;
328 MPI_Datatype dtr_left, dtr_right;
329 int num_left, num_right;
330
331 if (has_halo_1[i]) {
332 dts_left = comm_info.GetDatatypeSubarray(grid, grid.Iterators().NearBoundary1()[i]);
333 dtr_left = comm_info.GetDatatypeSubarray(grid, grid.Iterators().Halo1()[i]);
334 num_left = 1;
335 }else {
336 dts_left = MPI_INT;
337 dtr_left = MPI_INT;
338 num_left = 0;
339 }
340
341 if (has_halo_2[i]) {
342 dts_right = comm_info.GetDatatypeSubarray(grid, grid.Iterators().NearBoundary2()[i]);
343 dtr_right = comm_info.GetDatatypeSubarray(grid, grid.Iterators().Halo2()[i]);
344 num_right = 1;
345 }else {
346 dts_right = MPI_INT;
347 dtr_right = MPI_INT;
348 num_right = 0;
349 }
350
351 Timer::Start("MPI_Sendrecv");
352 MPI_Sendrecv(&(grid(0)), num_right, dts_right, neighbors.Right(), 3,
353 &(grid(0)), num_left, dtr_left, neighbors.Left(), 3,
354 comm, MPI_STATUS_IGNORE);
355 MPI_Sendrecv(&(grid(0)), num_left, dts_left, neighbors.Left(), 4,
356 &(grid(0)), num_right, dtr_right, neighbors.Right(), 4,
357 comm, MPI_STATUS_IGNORE);
358 Timer::Stop("MPI_Sendrecv");
359
360 }
361 }
362
363 Timer::Stop("CommToGhosts");
364}
365
366Grid& CommMPI::GetCoarserGrid(Multigrid& multigrid)
367{
368 std::map<const Grid*, TempGrid*>::iterator iter = coarser_grids.find(&multigrid());
369
370 if (iter == coarser_grids.end()) {
371 TempGrid* temp_grid = new TempGrid();
372 temp_grid->SetPropertiesToCoarser(multigrid(), BoundaryConditions());
373 iter = coarser_grids.insert(std::make_pair(&multigrid(), temp_grid)).first;
374 }
375
376 TempGrid& grid_coarser = *iter->second;
377
378 CommSubgrid(multigrid(multigrid.Level()-1), grid_coarser);
379
380 return grid_coarser;
381}
382
383Grid& CommMPI::GetFinerGrid(Multigrid& multigrid)
384{
385 std::map<const Grid*, TempGrid*>::iterator iter = finer_grids.find(&multigrid());
386
387 if (iter == finer_grids.end()) {
388 TempGrid* temp_grid = new TempGrid();
389 temp_grid->SetPropertiesToFiner(multigrid(), BoundaryConditions());
390 iter = finer_grids.insert(std::make_pair(&multigrid(), temp_grid)).first;
391 }
392
393 TempGrid& grid_finer = *iter->second;
394
395 CommSubgrid(multigrid(multigrid.Level()+1), grid_finer);
396
397 return grid_finer;
398}
399
400Grid& CommMPI::GetGlobalCoarseGrid(Multigrid& multigrid)
401{
402 if (global_coarse_grid == NULL) {
403
404 global_coarse_grid = new TempGrid();
405 global_coarse_grid->SetPropertiesToGlobalCoarseGrid();
406
407 }
408
409 CommSubgrid(multigrid(multigrid.MinLevel()), *global_coarse_grid);
410
411 return *global_coarse_grid;
412}
413
414//
415// Decodes a derived datatype created by MPI_Type_create_subarray
416// and adds values to grid.
417// I chose this approach since MPI neither supports a point-to-point reduce nor
418// a non-blocking global reduce.
419// For more information see Chapter 4.1 "DERIVED DATATYPES" of the MPI standard
420//
421void CommMPI::ReceiveAndAddSubgrid(Grid& grid, const MPI_Comm& comm, MPI_Datatype& type,
422 const int& rank, const int& tag)
423{
424 int counter = 0;
425 int ibuffer[11];
426 Index i;
427 MPI_Datatype base_type;
428 MPI_Aint address;
429
430 MPI_Type_get_contents(type, 11, 0, 1, ibuffer, &address, &base_type);
431
432 if (static_cast<int>(receive_buffer.size()) < ibuffer[4]*ibuffer[5]*ibuffer[6])
433 receive_buffer.resize(ibuffer[4]*ibuffer[5]*ibuffer[6]);
434
435 MPI_Recv(&(receive_buffer.front()), ibuffer[4]*ibuffer[5]*ibuffer[6], MPI_DOUBLE, rank, tag, comm, MPI_STATUS_IGNORE);
436
437 for (i.X()=ibuffer[7]; i.X()<ibuffer[7]+ibuffer[4]; ++i.X())
438 for (i.Y()=ibuffer[8]; i.Y()<ibuffer[8]+ibuffer[5]; ++i.Y())
439 for (i.Z()=ibuffer[9]; i.Z()<ibuffer[9]+ibuffer[6]; ++i.Z())
440 grid(i) += receive_buffer[counter++];
441}
442
443void CommMPI::CommParticles(const Grid& grid, std::list<Particle::Particle>& particles)
444{
445 Timer::Start("CommParticles");
446
447 Factory& factory = MG::GetFactory();
448
449 const vmg_int& num_particles_local = factory.GetObjectStorageVal<vmg_int>("PARTICLE_NUM_LOCAL");
450 vmg_float* x = factory.GetObjectStorageArray<vmg_float>("PARTICLE_POS_ARRAY");
451 vmg_float* q = factory.GetObjectStorageArray<vmg_float>("PARTICLE_CHARGE_ARRAY");
452
453 int rank, size;
454 MPI_Comm_rank(comm_global, &rank);
455 MPI_Comm_size(comm_global, &size);
456
457 Index index;
458 int global_extent[6*size], send_sizes[size], recv_sizes[size];
459 Index begin_remote[size], end_remote[size];
460 std::vector<int> pos_indices[size], charge_indices[size];
461 MPI_Datatype dt_temp;
462 std::vector<vmg_float> buffer_x, buffer_q;
463 std::vector<vmg_int> buffer_ind;
464
465 std::memcpy(&global_extent[6*rank], grid.Global().BeginLocal().vec(), 3*sizeof(vmg_float));
466 std::memcpy(&global_extent[6*rank+3], grid.Global().EndLocal().vec(), 3*sizeof(vmg_float));
467
468 MPI_Allgather(MPI_IN_PLACE, 6, MPI_INT, global_extent, 6, MPI_INT, MPI_COMM_WORLD);
469
470 for (int i=0; i<size; ++i) {
471 begin_remote[i] = static_cast<Vector>(&global_extent[6*i]);
472 end_remote[i] = static_cast<Vector>(&global_extent[6*i+3]);
473 }
474
475 for (int i=0; i<num_particles_local; ++i) {
476 index = static_cast<Index>((Vector(&x[3*i]) - grid.Extent().Begin()) / grid.Extent().MeshWidth());
477 for (int j=0; j<size; ++j)
478 if (index.IsInBounds(begin_remote[j], end_remote[j])) {
479 pos_indices[j].push_back(3*i);
480 charge_indices[j].push_back(i);
481 break;
482 }
483 }
484
485 /*
486 * Communicate which process gets how many particles
487 */
488 for (int i=0; i<size; ++i)
489 send_sizes[i] = charge_indices[i].size();
490
491 MPI_Alltoall(send_sizes, 1, MPI_INT, recv_sizes, 1, MPI_INT, comm_global);
492
493 /*
494 * Send particles
495 */
496 for (int i=0; i<size; ++i) {
497
498 if (pos_indices[i].size() > 0) {
499
500 MPI_Type_create_indexed_block(pos_indices[i].size(), 3, &pos_indices[i].front(), MPI_DOUBLE, &dt_temp);
501 MPI_Type_commit(&dt_temp);
502 MPI_Isend(x, 1, dt_temp, i, 0, MPI_COMM_WORLD, &Request());
503 MPI_Type_free(&dt_temp);
504
505 MPI_Type_create_indexed_block(charge_indices[i].size(), 1, &charge_indices[i].front(), MPI_DOUBLE, &dt_temp);
506 MPI_Type_commit(&dt_temp);
507 MPI_Isend(q, 1, dt_temp, i, 1, MPI_COMM_WORLD, &Request());
508 MPI_Type_free(&dt_temp);
509
510 MPI_Isend(&charge_indices[i].front(), charge_indices[i].size(), MPI_INT, i, 2, MPI_COMM_WORLD, &Request());
511
512 }
513
514 }
515
516 /*
517 * Receive particles
518 */
519 for (int i=0; i<size; ++i) {
520
521 if (recv_sizes[i] > 0) {
522
523 buffer_x.resize(3*recv_sizes[i]);
524 buffer_q.resize(recv_sizes[i]);
525 buffer_ind.resize(recv_sizes[i]);
526
527 MPI_Recv(&buffer_x.front(), 3*recv_sizes[i], MPI_DOUBLE, i, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
528 MPI_Recv(&buffer_q.front(), recv_sizes[i], MPI_DOUBLE, i, 1, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
529 MPI_Recv(&buffer_ind.front(), recv_sizes[i], MPI_INT, i, 2, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
530
531 particles.clear();
532
533 for (int j=0; j<recv_sizes[i]; ++j)
534 particles.push_back(Particle::Particle(&buffer_x[3*j], buffer_q[j], 0.0, 0.0, i, buffer_ind[j]));
535
536 }
537
538 }
539
540 ClearRequests();
541
542 Timer::Stop("CommParticles");
543}
544
545void CommMPI::CommLCListGhosts(const Grid& grid, Particle::LinkedCellList& lc)
546{
547 Timer::Start("CommLCListGhosts");
548
549 std::vector<vmg_float> send_buffer_1, send_buffer_2, recv_buffer;
550 vmg_int send_size_1, send_size_2, recv_size;
551 std::list<Particle::Particle>::iterator p_iter;
552 Grid::iterator g_iter;
553
554 Tuple<int> neighbors;
555
556 for (int i=2; i>=0; --i) {
557
558 /*
559 * Get ranks of neighbors
560 */
561 MPI_Cart_shift(comm_global, i, 1, &neighbors.Left(), &neighbors.Right());
562
563 /*
564 * Send to left
565 */
566 send_buffer_1.clear();
567
568 for (g_iter=lc.Iterators().Boundary1()[i].Begin(); g_iter!=lc.Iterators().Boundary1()[i].End(); ++g_iter)
569 lc(*g_iter).clear();
570
571 for (g_iter=lc.Iterators().NearBoundary1()[i].Begin(); g_iter!=lc.Iterators().NearBoundary1()[i].End(); ++g_iter)
572 for (p_iter=lc(*g_iter).begin(); p_iter!=lc(*g_iter).end(); ++p_iter) {
573
574 if (grid.Global().BeginLocal()[i] == 0)
575 for (int j=0; j<3; ++j)
576 send_buffer_1.push_back(p_iter->Pos()[j] + grid.Extent().Size()[j]);
577 else
578 for (int j=0; j<3; ++j)
579 send_buffer_1.push_back(p_iter->Pos()[j]);
580
581 send_buffer_1.push_back(p_iter->Charge());
582
583 }
584
585 send_size_1 = send_buffer_1.size();
586
587 MPI_Isend(&send_size_1, 1, MPI_INT, neighbors.Left(), 0, comm_global, &Request());
588
589 if (send_size_1 > 0)
590 MPI_Isend(&send_buffer_1.front(), send_size_1, MPI_DOUBLE, neighbors.Left(), 1, comm_global, &Request());
591
592 /*
593 * Send to right
594 */
595 send_buffer_2.clear();
596
597 for (g_iter=lc.Iterators().Boundary2()[i].Begin(); g_iter!=lc.Iterators().Boundary2()[i].End(); ++g_iter)
598 lc(*g_iter).clear();
599
600 for (g_iter=lc.Iterators().NearBoundary2()[i].Begin(); g_iter!=lc.Iterators().NearBoundary2()[i].End(); ++g_iter)
601 for (p_iter=lc(*g_iter).begin(); p_iter!=lc(*g_iter).end(); ++p_iter) {
602
603 if (grid.Global().EndLocal()[i] == grid.Global().SizeGlobal()[i])
604 for (int j=0; j<3; ++j)
605 send_buffer_2.push_back(p_iter->Pos()[j] - grid.Extent().Size()[j]);
606 else
607 for (int j=0; j<3; ++j)
608 send_buffer_2.push_back(p_iter->Pos()[j]);
609
610 send_buffer_2.push_back(p_iter->Charge());
611
612 }
613
614 send_size_2 = send_buffer_2.size();
615
616 MPI_Isend(&send_size_2, 1, MPI_INT, neighbors.Right(), 2, comm_global, &Request());
617
618 if (send_size_2 > 0)
619 MPI_Isend(&send_buffer_2.front(), send_size_2, MPI_DOUBLE, neighbors.Right(), 3, comm_global, &Request());
620
621 /*
622 * Receive from right
623 */
624 MPI_Recv(&recv_size, 1, MPI_INT, neighbors.Right(), 0, comm_global, MPI_STATUS_IGNORE);
625
626 if (recv_size > 0) {
627
628 recv_buffer.resize(recv_size);
629 MPI_Recv(&recv_buffer.front(), recv_size, MPI_DOUBLE, neighbors.Right(), 1, comm_global, MPI_STATUS_IGNORE);
630
631 for (vmg_int i=0; i<recv_size; i+=4)
632 lc.AddParticleToComplete(&recv_buffer[i], recv_buffer[i+3]);
633
634 }
635
636 /*
637 * Receive from left
638 */
639 MPI_Recv(&recv_size, 1, MPI_INT, neighbors.Left(), 2, comm_global, MPI_STATUS_IGNORE);
640
641 if (recv_size > 0) {
642
643 recv_buffer.resize(recv_size);
644 MPI_Recv(&recv_buffer.front(), recv_size, MPI_DOUBLE, neighbors.Left(), 3, comm_global, MPI_STATUS_IGNORE);
645
646 for (vmg_int i=0; i<recv_size; i+=4)
647 lc.AddParticleToComplete(&recv_buffer[i], recv_buffer[i+3]);
648
649 }
650
651
652 ClearRequests();
653 }
654
655 Timer::Stop("CommLCListGhosts");
656}
657
658void CommMPI::CommAddPotential(std::list<Particle::Particle>& particles)
659{
660 Timer::Start("CommAddPotential");
661
662 std::list<Particle::Particle>::iterator iter;
663
664 if (!win_created) {
665 vmg_float* p = MG::GetFactory().GetObjectStorageArray<vmg_float>("PARTICLE_POTENTIAL_ARRAY");
666 const vmg_int& num_particles_local = MG::GetFactory().GetObjectStorageVal<vmg_int>("PARTICLE_NUM_LOCAL");
667 MPI_Win_create(p, num_particles_local*sizeof(vmg_float), sizeof(vmg_float), info, comm_global, &win);
668 win_created = true;
669 }
670
671 MPI_Win_fence(MPI_MODE_NOPRECEDE, win);
672
673 for (iter=particles.begin(); iter!=particles.end(); ++iter)
674 MPI_Accumulate(&iter->Pot(), 1, MPI_DOUBLE, iter->Rank(), iter->Index(), 1, MPI_DOUBLE, MPI_SUM, win);
675
676 MPI_Win_fence(MPI_MODE_NOSTORE | MPI_MODE_NOSUCCEED, win);
677
678 Timer::Stop("CommAddPotential");
679}
680
681void CommMPI::CommAddPotential(Particle::LinkedCellList& lc)
682{
683 Timer::Start("CommAddPotential");
684
685 Grid::iterator g_iter;
686 std::list<Particle::Particle>::iterator p_iter;
687
688 if (!win_created) {
689 vmg_float* p = MG::GetFactory().GetObjectStorageArray<vmg_float>("PARTICLE_POTENTIAL_ARRAY");
690 const vmg_int& num_particles_local = MG::GetFactory().GetObjectStorageVal<vmg_int>("PARTICLE_NUM_LOCAL");
691 MPI_Win_create(p, num_particles_local*sizeof(vmg_float), sizeof(vmg_float), info, comm_global, &win);
692 win_created = true;
693 }
694
695 MPI_Win_fence(MPI_MODE_NOPRECEDE, win);
696
697 for (g_iter=lc.Iterators().Local().Begin(); g_iter!=lc.Iterators().Local().End(); ++g_iter)
698 for (p_iter=lc(*g_iter).begin(); p_iter!=lc(*g_iter).end(); ++p_iter)
699 MPI_Accumulate(&p_iter->Pot(), 1, MPI_DOUBLE, p_iter->Rank(), p_iter->Index(), 1, MPI_DOUBLE, MPI_SUM, win);
700
701 MPI_Win_fence(MPI_MODE_NOSTORE | MPI_MODE_NOSUCCEED, win);
702
703 Timer::Stop("CommAddPotential");
704}
705
706vmg_float CommMPI::GlobalSum(const vmg_float& value)
707{
708 vmg_float result = value;
709
710 MPI_Allreduce(MPI_IN_PLACE, &result, 1, MPI_DOUBLE, MPI_SUM, comm_global);
711
712 return result;
713}
714
715vmg_float CommMPI::GlobalSumRoot(const vmg_float& value)
716{
717 vmg_float recv_buffer;
718 vmg_float send_buffer = value;
719
720 MPI_Reduce(&send_buffer, &recv_buffer, 1, MPI_DOUBLE, MPI_SUM, 0, comm_global);
721
722 return recv_buffer;
723}
724
725void CommMPI::GlobalSumArray(vmg_float* array, const vmg_int& size)
726{
727 MPI_Allreduce(MPI_IN_PLACE, array, size, MPI_DOUBLE, MPI_SUM, comm_global);
728}
729
730void CommMPI::PrintString(const char* format, ...)
731{
732 va_list args;
733 va_start(args, format);
734 vsprintf(print_buffer, format, args);
735 printf("VMG: Rank %d: %s\n", GlobalRank(), print_buffer);
736 va_end(args);
737}
738
739void CommMPI::PrintStringOnce(const char* format, ...)
740{
741 if (GlobalRank() == 0) {
742 va_list args;
743 va_start(args, format);
744 vsprintf(print_buffer, format, args);
745 printf("VMG: Rank %d: %s\n", GlobalRank(), print_buffer);
746 va_end(args);
747 }
748}
749
750void CommMPI::PrintXML(const std::string& filename, const std::string& xml_data)
751{
752 MPI_File file;
753 std::stringstream path, xml_header;
754
755 pugi::xml_document doc;
756 pugi::xml_node node_data = doc.append_child("Global").append_child("NumProcesses").append_child(pugi::node_pcdata);
757 node_data.set_value(Helper::ToString(GlobalProcs().Product()).c_str());
758 doc.save(xml_header);
759
760 path << OutputPath() << filename;
761
762 char* filename_array = Helper::GetCharArray(path.str());
763 char* xml_header_array = Helper::GetCharArray(xml_header.str());
764 char* str_array = Helper::GetCharArray(xml_data);
765
766 MPI_File_open(MPI_COMM_SELF, filename_array, MPI_MODE_WRONLY|MPI_MODE_CREATE, MPI_INFO_NULL, &file);
767 MPI_File_set_size(file, 0);
768 MPI_File_write(file, xml_header_array, xml_header.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
769 MPI_File_write(file, str_array, xml_data.size(), MPI_CHAR, MPI_STATUS_IGNORE);
770 MPI_File_close(&file);
771
772 delete [] filename_array;
773 delete [] xml_header_array;
774 delete [] str_array;
775}
776
777void CommMPI::PrintXMLAll(const std::string& filename, const std::string& xml_data)
778{
779 MPI_File file;
780 std::stringstream path;
781
782 path << OutputPath() << filename;
783
784 char* filename_array = Helper::GetCharArray(path.str());
785 char* str_array = Helper::GetCharArray(xml_data);
786
787 MPI_File_open(comm_global, filename_array, MPI_MODE_WRONLY|MPI_MODE_CREATE, MPI_INFO_NULL, &file);
788 MPI_File_set_size(file, 0);
789
790 if (GlobalRank() == 0) {
791 std::stringstream xml_header;
792 pugi::xml_document doc;
793 pugi::xml_node node_data = doc.append_child("Global").append_child("NumProcesses").append_child(pugi::node_pcdata);
794 node_data.set_value(Helper::ToString(GlobalProcs().Product()).c_str());
795 doc.save(xml_header);
796
797 char* xml_header_array = Helper::GetCharArray(xml_header.str());
798
799 MPI_File_write_shared(file, xml_header_array, xml_header.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
800
801 delete [] xml_header_array;
802 }
803
804 MPI_File_write_ordered(file, str_array, xml_data.size(), MPI_CHAR, MPI_STATUS_IGNORE);
805 MPI_File_close(&file);
806
807 delete [] filename_array;
808 delete [] str_array;
809}
810
811void CommMPI::PrintAllSettings()
812{
813 std::stringstream buf;
814 MPI_File file;
815 MPI_Status status;
816
817 Multigrid* mg = MG::GetFactory().Get("SOL")->Cast<Multigrid>();
818
819 buf << OutputPath() << "settings.txt";
820 char *filename = Helper::GetCharArray(buf.str());
821
822 MPI_File_open(comm_global, filename, MPI_MODE_WRONLY|MPI_MODE_CREATE, MPI_INFO_NULL, &file);
823 MPI_File_set_size(file, 0);
824 MPI_File_close(&file);
825
826 for (int i=mg->MinLevel(); i<=mg->MaxLevel(); ++i) {
827
828 MPI_Comm comm = comm_info.GetCommunicator((*mg)(i));
829
830 if (comm != MPI_COMM_NULL) {
831
832 MPI_File_open(comm, filename, MPI_MODE_WRONLY|MPI_MODE_CREATE|MPI_MODE_APPEND, MPI_INFO_NULL, &file);
833
834 int comm_rank;
835 MPI_Comm_rank(comm, &comm_rank);
836
837 if (comm_rank == 0) {
838
839 buf.str("");
840 buf << "PROCESS INDEPENDENT INFORMATION:" << std::endl
841 << "GRID LEVEL " << i << std::endl
842 << "BEGINFINEST: " << (*mg)(i).Global().BeginFinest() << std::endl
843 << "ENDFINEST: " << (*mg)(i).Global().EndFinest() << std::endl
844 << "SIZEFINEST: " << (*mg)(i).Global().SizeFinest() << std::endl
845 << "SIZEGLOBAL: " << (*mg)(i).Global().SizeGlobal() << std::endl
846 << "BOUNDARYTYPE: " << (*mg)(i).Global().BoundaryType() << std::endl
847 << "EXTENT" << std::endl
848 << "SIZE: " << (*mg)(i).Extent().Size() << std::endl
849 << "SIZEFACTOR: " << (*mg)(i).Extent().SizeFactor() << std::endl
850 << "BEGIN: " << (*mg)(i).Extent().Begin() << std::endl
851 << "END: " << (*mg)(i).Extent().End() << std::endl
852 << "MESHWIDTH: " << (*mg)(i).Extent().MeshWidth() << std::endl
853 << std::endl
854 << "PROCESS DEPENDENT INFORMATION:" << std::endl;
855
856 char* char_buf = Helper::GetCharArray(buf.str());
857
858 MPI_File_write_shared(file, char_buf, buf.str().size(), MPI_CHAR, &status);
859
860 delete [] char_buf;
861
862 }
863
864 buf.str("");
865 buf << "PROCESS " << comm_rank << ":" << std::endl
866 << "GLOBAL" << std::endl
867 << "BEGINLOCAL: " << (*mg)(i).Global().BeginLocal() << std::endl
868 << "ENDLOCAL: " << (*mg)(i).Global().EndLocal() << std::endl
869 << "SIZELOCAL: " << (*mg)(i).Global().SizeLocal() << std::endl
870 << "LOCAL" << std::endl
871 << "SIZE: " << (*mg)(i).Local().Size() << std::endl
872 << "SIZETOTAL: " << (*mg)(i).Local().SizeTotal() << std::endl
873 << "BEGIN: " << (*mg)(i).Local().Begin() << std::endl
874 << "END: " << (*mg)(i).Local().End() << std::endl
875 << "HALOBEGIN1: " << (*mg)(i).Local().HaloBegin1() << std::endl
876 << "HALOEND1: " << (*mg)(i).Local().HaloEnd1() << std::endl
877 << "HALOBEGIN2: " << (*mg)(i).Local().HaloBegin2() << std::endl
878 << "HALOEND2: " << (*mg)(i).Local().HaloEnd2() << std::endl
879 << "BOUNDARYBEGIN1: " << (*mg)(i).Local().BoundaryBegin1() << std::endl
880 << "BOUNDARYEND1: " << (*mg)(i).Local().BoundaryEnd1() << std::endl
881 << "BOUNDARYBEGIN2: " << (*mg)(i).Local().BoundaryBegin2() << std::endl
882 << "BOUNDARYEND2: " << (*mg)(i).Local().BoundaryEnd2() << std::endl
883 << "ALIGNMENTBEGIN: " << (*mg)(i).Local().AlignmentBegin() << std::endl
884 << "ALIGNMENTEND: " << (*mg)(i).Local().AlignmentEnd() << std::endl
885 << std::endl;
886
887 char* char_buf = Helper::GetCharArray(buf.str());
888
889 MPI_File_write_ordered(file, char_buf, buf.str().size(), MPI_CHAR, &status);
890
891 delete [] char_buf;
892
893 MPI_File_close(&file);
894
895 }
896 }
897
898 delete [] filename;
899
900}
901
902void CommMPI::PrintDefect(Grid& sol, Grid& rhs, const char* information)
903{
904 TempGrid *temp = MG::GetTempGrid();
905 temp->SetProperties(sol);
906 temp->ImportFromResidual(sol, rhs);
907
908 PrintInnerGrid(*temp, information);
909}
910
911void CommMPI::PrintGrid(Grid& grid, const char* information)
912{
913 Index i;
914 std::stringstream buf;
915
916 Index begin = 0;
917 Index end = grid.Local().End();
918 Index begin_local, end_local, begin_global, end_global;
919
920 CommToGhosts(grid);
921
922 for (int j=0; j<3; ++j)
923 if (grid.Global().EndLocal()[j] == grid.Global().SizeGlobal()[j])
924 end[j] = grid.Local().SizeTotal()[j];
925
926 for (i.Z()=begin.Z(); i.Z()<end.Z(); ++i.Z())
927 for (i.Y()=begin.Y(); i.Y()<end.Y(); ++i.Y())
928 for (i.X()=begin.X(); i.X()<end.X(); ++i.X())
929 buf << std::scientific << grid.GetVal(i) << " ";
930
931 begin_local = grid.Global().BeginLocal();
932 end_local = grid.Global().EndLocal();
933 begin_global = 0;
934
935 for (int i=0; i<3; ++i) {
936
937 if (BoundaryConditions()[i] == Dirichlet ||
938 BoundaryConditions()[i] == Open) {
939
940 if (begin_local[i] != 0)
941 begin_local[i] -= grid.Local().HaloEnd1()[i] - grid.Local().HaloBegin1()[i];
942
943 end_global[i] = grid.Global().SizeGlobal()[i];
944
945 }else if (BoundaryConditions()[i] == Periodic) {
946
947 if (end_local[i] == grid.Global().SizeGlobal()[i])
948 end_local[i] += grid.Local().HaloEnd1()[i] - grid.Local().HaloBegin1()[i] + grid.Local().HaloEnd2()[i] - grid.Local().HaloBegin2()[i];
949 else
950 end_local[i] += grid.Local().HaloEnd1()[i] - grid.Local().HaloBegin1()[i];
951
952 end_global[i] = grid.Global().SizeGlobal()[i] +
953 grid.Local().HaloEnd1()[i] - grid.Local().HaloBegin1()[i] +
954 grid.Local().HaloEnd2()[i] - grid.Local().HaloBegin2()[i];
955
956 }
957
958 }
959
960 CreateOutputFiles(grid, buf, information,
961 begin_global, end_global,
962 begin_local, end_local);
963}
964
965void CommMPI::PrintCorrectedGrid(Grid& grid, const char* information)
966{
967 Index i;
968
969 std::stringstream buf;
970
971 Index begin = 0;
972 Index end = grid.Local().End();
973 Index begin_local, end_local, begin_global, end_global;
974
975 CommToGhosts(grid);
976
977 for (int j=0; j<3; ++j)
978 if (grid.Global().EndLocal()[j] == grid.Global().SizeGlobal()[j])
979 end[j] = grid.Local().SizeTotal()[j];
980
981 for (i.Z()=begin.Z(); i.Z()<end.Z(); ++i.Z())
982 for (i.Y()=begin.Y(); i.Y()<end.Y(); ++i.Y())
983 for (i.X()=begin.X(); i.X()<end.X(); ++i.X())
984 buf << std::scientific << grid.GetVal(i) - Grid::Correction() << " ";
985
986 begin_local = grid.Global().BeginLocal();
987 end_local = grid.Global().EndLocal();
988 begin_global = 0;
989
990 for (int i=0; i<3; ++i) {
991
992 if (BoundaryConditions()[i] == Dirichlet ||
993 BoundaryConditions()[i] == Open) {
994
995 if (begin_local[i] != 0)
996 --begin_local[i];
997
998 end_global[i] = grid.Global().SizeGlobal()[i];
999
1000 }else if (BoundaryConditions()[i] == Periodic) {
1001
1002 if (begin_local[i] != 0)
1003 --begin_local[i];
1004
1005 if (end_local[i] == grid.Global().SizeGlobal()[i])
1006 end_local[i] += 2;
1007
1008 end_global[i] = grid.Global().SizeGlobal()[i] + 2;
1009
1010 }
1011
1012 }
1013
1014 CreateOutputFiles(grid, buf, information,
1015 0, grid.Global().SizeGlobal()+2,
1016 begin, end);
1017}
1018
1019void CommMPI::PrintInnerGrid(Grid& grid, const char* information)
1020{
1021 Index i;
1022 Index begin_local, end_local, begin_global, end_global;
1023 std::stringstream buf;
1024
1025 CommToGhosts(grid);
1026
1027 Index begin = 0;
1028 Index end = grid.Local().End();
1029
1030 /* VTK parallel ImageData files seem to require some overlap. */
1031 for (int j=0; j<3; ++j)
1032 if (grid.Global().BeginLocal()[j] == 0)
1033 ++begin[j];
1034
1035 /* Write grid values into buffer */
1036 for (i.Z()=begin.Z(); i.Z()<end.Z(); ++i.Z())
1037 for (i.Y()=begin.Y(); i.Y()<end.Y(); ++i.Y())
1038 for (i.X()=begin.X(); i.X()<end.X(); ++i.X())
1039 buf << std::scientific << grid.GetVal(i) << " ";
1040
1041 for (int j=0; j<3; ++j) {
1042
1043 if (BoundaryConditions()[j] == Dirichlet ||
1044 BoundaryConditions()[j] == Open) {
1045
1046 begin_global[j] = 1;
1047 end_global[j] = grid.Global().SizeGlobal()[j] - 1;
1048
1049 if (grid.Global().BeginLocal()[j] == 0)
1050 begin_local[j] = 1;
1051 else
1052 begin_local[j] = grid.Global().BeginLocal()[j] - 1;
1053
1054 if (grid.Global().EndLocal()[j] == grid.Global().SizeGlobal()[j])
1055 end_local[j] = grid.Global().EndLocal()[j] - 1;
1056 else
1057 end_local[j] = grid.Global().EndLocal()[j];
1058
1059 }else if (BoundaryConditions()[j] == Periodic) {
1060
1061 begin_global[j] = 0;
1062 end_global[j] = grid.Global().SizeGlobal()[j];
1063
1064 if (grid.Global().BeginLocal()[j] == 0)
1065 begin_local[j] = grid.Global().BeginLocal()[j];
1066 else
1067 begin_local[j] = grid.Global().BeginLocal()[j] - 1;
1068
1069 end_local[j] = grid.Global().EndLocal()[j];
1070
1071 }else {
1072
1073 assert(0 == "Boundary condition not handled properly.");
1074
1075 }
1076
1077 }
1078
1079 CreateOutputFiles(grid, buf, information,
1080 begin_global, end_global,
1081 begin_local, end_local);
1082}
1083
1084void CommMPI::PrintInnerCorrectedGrid(Grid& grid, const char* information)
1085{
1086 Index i;
1087 const Index begin = grid.Local().Begin();
1088 const Index end = grid.Local().End();
1089 std::stringstream buf;
1090
1091 for (i.Z()=begin.Z(); i.Z()<end.Z(); ++i.Z())
1092 for (i.Y()=begin.Y(); i.Y()<end.Y(); ++i.Y())
1093 for (i.X()=begin.X(); i.X()<end.X(); ++i.X())
1094 buf << std::scientific << grid.GetVal(i) - Grid::Correction() << " ";
1095
1096 CreateOutputFiles(grid, buf, information,
1097 0, grid.Global().SizeGlobal(),
1098 grid.Global().BeginLocal(), grid.Global().EndLocal());
1099}
1100
1101void CommMPI::CreateOutputFiles(const Grid& grid, const std::stringstream& serial_data, const char* information,
1102 const Index& begin_global, const Index& end_global,
1103 const Index& begin_local, const Index& end_local)
1104{
1105 int output_count = OutputCount();
1106
1107 MPI_Comm comm = comm_info.GetCommunicator(grid);
1108
1109 if (comm != MPI_COMM_NULL) {
1110
1111 MPI_File file;
1112 std::string conv_information = Helper::ReplaceWhitespaces(information, "_");
1113
1114 CreateParallelOutputFile(grid, comm, output_count, conv_information.c_str(),
1115 begin_global, end_global, begin_local, end_local);
1116
1117 file = CreateSerialOutputFile(grid, comm, output_count, conv_information.c_str(),
1118 begin_global, end_global, begin_local, end_local);
1119
1120 char *char_buf = Helper::GetCharArray(serial_data.str());
1121 MPI_File_write(file, char_buf, serial_data.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
1122 delete [] char_buf;
1123
1124 FinalizeSerialOutputFile(file);
1125
1126 }
1127}
1128
1129void CommMPI::CreateParallelOutputFile(const Grid& grid, MPI_Comm& comm,
1130 const int& output_count, const char* information,
1131 const Index& begin_global, const Index& end_global,
1132 const Index& begin_local, const Index& end_local)
1133{
1134 int rank;
1135 MPI_File file;
1136 char parallel_filename[513], serial_filename[513];
1137 std::stringstream buf;
1138
1139 MPI_Comm_rank(comm, &rank);
1140
1141 sprintf(parallel_filename, "%s%04d.pvti", OutputPath().c_str(), output_count);
1142 sprintf(serial_filename, "%04d_%d.vti", output_count, rank);
1143
1144 MPI_File_open(comm, parallel_filename, MPI_MODE_WRONLY|MPI_MODE_CREATE, MPI_INFO_NULL, &file);
1145 MPI_File_set_size(file, 0);
1146
1147 if (rank == 0) {
1148
1149 buf << "<?xml version=\"1.0\"?>" << std::endl
1150 << "<VTKFile type=\"PImageData\" version=\"0.1\" byte_order=\"LittleEndian\">" << std::endl
1151 << " <PImageData WholeExtent=\"";
1152
1153 for (int i=0; i<3; ++i)
1154 buf << begin_global[i] << " " << end_global[i] - 1 << " ";
1155
1156 buf << "\"" << std::endl
1157 << " GhostLevel=\"0\" Origin=\"0 0 0\" Spacing=\"";
1158
1159 for (int i=0; i<3; ++i)
1160 buf << grid.Extent().MeshWidth()[i] << " ";
1161
1162 buf << "\">" << std::endl
1163 << " <PPointData Scalars=\"" << information << "\">" << std::endl
1164 << " <PDataArray type=\"Float32\" Name=\"" << information << "\"/>" << std::endl
1165 << " </PPointData>" << std::endl;
1166
1167 char* char_buf = Helper::GetCharArray(buf.str());
1168
1169 MPI_File_write_shared(file, char_buf, buf.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
1170
1171 delete [] char_buf;
1172 }
1173
1174 buf.str("");
1175 buf << " <Piece Extent=\"";
1176
1177 for (int i=0; i<3; ++i)
1178 buf << begin_local[i] << " " << end_local[i] - 1 << " ";
1179
1180 buf << "\" Source=\"" << serial_filename << "\"/>" << std::endl;
1181
1182 char* char_buf = Helper::GetCharArray(buf.str());
1183
1184 MPI_File_write_ordered(file, char_buf, buf.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
1185
1186 delete [] char_buf;
1187
1188 if (rank == 0) {
1189
1190 buf.str("");
1191
1192 buf << " </PImageData>" << std::endl
1193 << "</VTKFile>" << std::endl;
1194
1195 char* char_buf = Helper::GetCharArray(buf.str());
1196
1197 MPI_File_write_shared(file, char_buf, buf.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
1198
1199 delete [] char_buf;
1200 }
1201
1202 MPI_File_close(&file);
1203}
1204
1205MPI_File CommMPI::CreateSerialOutputFile(const Grid& grid, MPI_Comm& comm,
1206 const int& output_count, const char* information,
1207 const Index& begin_global, const Index& end_global,
1208 const Index& begin_local, const Index& end_local)
1209{
1210 char serial_filename[513];
1211 int rank;
1212 MPI_File file;
1213 std::stringstream buf;
1214
1215 MPI_Comm_rank(comm_global, &rank);
1216
1217 sprintf(serial_filename, "%s%04d_%d.vti", OutputPath().c_str(), output_count, rank);
1218
1219 MPI_File_open(MPI_COMM_SELF, serial_filename, MPI_MODE_WRONLY|MPI_MODE_CREATE, MPI_INFO_NULL, &file);
1220
1221 buf << "<?xml version=\"1.0\"?>" << std::endl
1222 << "<VTKFile type=\"ImageData\" version=\"0.1\" byte_order=\"LittleEndian\">" << std::endl
1223 << " <ImageData WholeExtent=\"";
1224
1225 for (int i=0; i<3; ++i)
1226 buf << begin_global[i] << " " << end_global[i] - 1 << " ";
1227
1228 buf << "\"" << std::endl
1229 << " Origin=\"0 0 0\" Spacing=\"";
1230
1231 for (int i=0; i<3; ++i)
1232 buf << grid.Extent().MeshWidth()[i] << " ";
1233
1234 buf << "\">" << std::endl
1235 << " <Piece Extent=\"";
1236
1237 for (int i=0; i<3; ++i)
1238 buf << begin_local[i] << " " << end_local[i] - 1 << " ";
1239
1240 buf << "\">" << std::endl
1241 << " <PointData Scalars=\"" << information << "\">" << std::endl
1242 << " <DataArray type=\"Float32\" Name=\"" << information << "\" format=\"ascii\">" << std::endl
1243 << " ";
1244
1245 char* char_buf = Helper::GetCharArray(buf.str());
1246 MPI_File_write(file, char_buf, buf.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
1247 delete [] char_buf;
1248
1249 return file;
1250}
1251
1252void CommMPI::FinalizeSerialOutputFile(MPI_File& file)
1253{
1254 std::stringstream buf;
1255
1256 buf << std::endl
1257 << " </DataArray>" << std::endl
1258 << " </PointData>" << std::endl
1259 << " </Piece>" << std::endl
1260 << " </ImageData>" << std::endl
1261 << "</VTKFile>" << std::endl;
1262
1263 char* char_buf = Helper::GetCharArray(buf.str());
1264 MPI_File_write(file, char_buf, buf.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
1265 delete [] char_buf;
1266
1267 MPI_File_close(&file);
1268}
1269
1270int CommMPI::GlobalRank() const
1271{
1272 int rank;
1273
1274 MPI_Comm_rank(comm_global, &rank);
1275
1276 return rank;
1277}
1278
1279Index CommMPI::GlobalPos() const
1280{
1281 Index dims, periods, coords;
1282
1283 MPI_Cart_get(comm_global, 3, dims.vec(), periods.vec(), coords.vec());
1284
1285 return coords;
1286}
1287
1288Index CommMPI::GlobalProcs() const
1289{
1290 Index dims, periods, coords;
1291
1292 MPI_Cart_get(comm_global, 3, dims.vec(), periods.vec(), coords.vec());
1293
1294 return dims;
1295}
1296
1297void CommMPI::InitCommMPI(const MPI_Comm& comm)
1298{
1299 int status, size, rank;
1300 int dims[3] = {0, 0, 0};
1301 int periods[3];
1302
1303 global_coarse_grid = NULL;
1304
1305 for (int i=0; i<3; ++i)
1306 periods[i] = (BoundaryConditions()[i] == Periodic ? 1 : 0);
1307
1308 MPI_Comm_size(comm, &size);
1309 MPI_Comm_rank(comm, &rank);
1310
1311 MPI_Topo_test(comm, &status);
1312
1313 if (status == MPI_CART) {
1314
1315 comm_global = comm;
1316
1317 }else {
1318
1319 MPI_Dims_create(size, 3, dims);
1320
1321 MPI_Cart_create(comm, 3, dims, periods, 1, &comm_global);
1322
1323 }
1324
1325 comm_info.SetCommGlobal(comm_global);
1326
1327 MPI_Info_create(&info);
1328 char key[] = "no_locks";
1329 char val[] = "true";
1330 MPI_Info_set(info, key, val);
1331
1332}
1333
1334CommMPI::~CommMPI()
1335{
1336 std::map<const Grid*, TempGrid*>::iterator iter;
1337
1338 for (iter=finer_grids.begin(); iter!=finer_grids.end(); ++iter)
1339 delete iter->second;
1340
1341 for (iter=coarser_grids.begin(); iter!=coarser_grids.end(); ++iter)
1342 delete iter->second;
1343
1344 MPI_Comm_free(&comm_global);
1345
1346 if (win_created)
1347 MPI_Win_free(&win);
1348
1349 MPI_Info_free(&info);
1350
1351 delete global_coarse_grid;
1352}
1353
1354std::string CommMPI::CreateOutputDirectory()
1355{
1356#ifdef HAVE_BOOST_FILESYSTEM
1357 std::string path, unique_path;
1358 std::stringstream unique_suffix;
1359 int suffix_counter = 0;
1360 char buffer[129];
1361 time_t rawtime;
1362 struct tm *timeinfo;
1363 int path_size;
1364 char* path_array;
1365
1366 if (GlobalRank() == 0) {
1367
1368 time(&rawtime);
1369 timeinfo = localtime(&rawtime);
1370 strftime(buffer, 128, "./output/%Y_%m_%d_%H_%M_%S/", timeinfo);
1371 path = buffer;
1372
1373 if (!fs::exists("output"))
1374 fs::create_directory("output");
1375
1376 unique_path = path;
1377
1378 while (fs::exists(unique_path.c_str())) {
1379
1380 unique_suffix.str("");
1381 unique_suffix << "_" << suffix_counter++ << "/";
1382
1383 unique_path = path;
1384 unique_path.replace(unique_path.size()-1, 1, unique_suffix.str());
1385
1386 }
1387
1388 fs::create_directory(unique_path.c_str());
1389
1390 path_size = unique_path.size() + 1;
1391 path_array = Helper::GetCharArray(unique_path);
1392
1393 MPI_Bcast(&path_size, 1, MPI_INT, 0, comm_global);
1394
1395 }else {
1396
1397 MPI_Bcast(&path_size, 1, MPI_INT, 0, comm_global);
1398 path_array = new char[path_size];
1399
1400 }
1401
1402 MPI_Bcast(path_array, path_size, MPI_CHAR, 0, comm_global);
1403
1404 unique_path = path_array;
1405
1406 delete [] path_array;
1407
1408 return unique_path;
1409
1410#else
1411
1412 return "./";
1413
1414#endif
1415}
1416
1417
1418#endif /* HAVE_MPI */
Note: See TracBrowser for help on using the repository browser.