| 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>
|
|---|
| 20 | namespace 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 |
|
|---|
| 35 | static char print_buffer[512];
|
|---|
| 36 |
|
|---|
| 37 | using namespace VMG;
|
|---|
| 38 |
|
|---|
| 39 | void 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 |
|
|---|
| 155 | void 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 |
|
|---|
| 266 | void 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 |
|
|---|
| 311 | void CommMPI::CommToGhosts(Grid& grid)
|
|---|
| 312 | {
|
|---|
| 313 | Timer::Start("CommToGhosts");
|
|---|
| 314 |
|
|---|
| 315 | Timer::Start("CommToGhostsGetCommunicator");
|
|---|
| 316 | MPI_Comm comm = comm_info.GetCommunicator(grid);
|
|---|
| 317 | Timer::Stop("CommToGhostsGetCommunicator");
|
|---|
| 318 |
|
|---|
| 319 | if (comm != MPI_COMM_NULL) {
|
|---|
| 320 |
|
|---|
| 321 | const Index has_halo_1 = grid.Local().HasHalo1();
|
|---|
| 322 | const Index has_halo_2 = grid.Local().HasHalo2();
|
|---|
| 323 | Tuple<int> neighbors;
|
|---|
| 324 |
|
|---|
| 325 | for (int i=2; i>=0; --i) {
|
|---|
| 326 |
|
|---|
| 327 | Timer::Start("CommToGhostsCartShift");
|
|---|
| 328 | MPI_Cart_shift(comm, i, 1, &neighbors.Left(), &neighbors.Right());
|
|---|
| 329 | Timer::Stop("CommToGhostsCartShift");
|
|---|
| 330 |
|
|---|
| 331 | MPI_Datatype dts_left, dts_right;
|
|---|
| 332 | MPI_Datatype dtr_left, dtr_right;
|
|---|
| 333 | int num_left, num_right;
|
|---|
| 334 |
|
|---|
| 335 | if (has_halo_1[i]) {
|
|---|
| 336 | Timer::Start("CommToGhostsGetDatatypeSubarray");
|
|---|
| 337 | dts_left = comm_info.GetDatatypeSubarray(grid, grid.Iterators().NearBoundary1()[i]);
|
|---|
| 338 | dtr_left = comm_info.GetDatatypeSubarray(grid, grid.Iterators().Halo1()[i]);
|
|---|
| 339 | Timer::Stop("CommToGhostsGetDatatypeSubarray");
|
|---|
| 340 | num_left = 1;
|
|---|
| 341 | }else {
|
|---|
| 342 | dts_left = MPI_INT;
|
|---|
| 343 | dtr_left = MPI_INT;
|
|---|
| 344 | num_left = 0;
|
|---|
| 345 | }
|
|---|
| 346 |
|
|---|
| 347 | if (has_halo_2[i]) {
|
|---|
| 348 | Timer::Start("CommToGhostsGetDatatypeSubarray");
|
|---|
| 349 | dts_right = comm_info.GetDatatypeSubarray(grid, grid.Iterators().NearBoundary2()[i]);
|
|---|
| 350 | dtr_right = comm_info.GetDatatypeSubarray(grid, grid.Iterators().Halo2()[i]);
|
|---|
| 351 | Timer::Stop("CommToGhostsGetDatatypeSubarray");
|
|---|
| 352 | num_right = 1;
|
|---|
| 353 | }else {
|
|---|
| 354 | dts_right = MPI_INT;
|
|---|
| 355 | dtr_right = MPI_INT;
|
|---|
| 356 | num_right = 0;
|
|---|
| 357 | }
|
|---|
| 358 |
|
|---|
| 359 | Timer::Start("CommToGhostsCommunication");
|
|---|
| 360 | MPI_Sendrecv(&(grid(0)), num_right, dts_right, neighbors.Right(), 3,
|
|---|
| 361 | &(grid(0)), num_left, dtr_left, neighbors.Left(), 3,
|
|---|
| 362 | comm, MPI_STATUS_IGNORE);
|
|---|
| 363 | MPI_Sendrecv(&(grid(0)), num_left, dts_left, neighbors.Left(), 4,
|
|---|
| 364 | &(grid(0)), num_right, dtr_right, neighbors.Right(), 4,
|
|---|
| 365 | comm, MPI_STATUS_IGNORE);
|
|---|
| 366 | Timer::Stop("CommToGhostsCommunication");
|
|---|
| 367 |
|
|---|
| 368 | }
|
|---|
| 369 | }
|
|---|
| 370 |
|
|---|
| 371 | Timer::Stop("CommToGhosts");
|
|---|
| 372 | }
|
|---|
| 373 |
|
|---|
| 374 | void CommMPI::CommToGhostsRB(Grid& grid, const int& offset)
|
|---|
| 375 | {
|
|---|
| 376 | MPI_Comm comm = comm_info.GetCommunicator(grid);
|
|---|
| 377 |
|
|---|
| 378 | if (comm != MPI_COMM_NULL) {
|
|---|
| 379 |
|
|---|
| 380 | Index i;
|
|---|
| 381 | int num_left, num_right;
|
|---|
| 382 | MPI_Datatype dts_left, dts_right;
|
|---|
| 383 | MPI_Datatype dtr_left, dtr_right;
|
|---|
| 384 | Tuple<int> neighbors;
|
|---|
| 385 | std::vector<vmg_float> buffer_send_left, buffer_send_right;
|
|---|
| 386 | std::vector<vmg_float> buffer_recv_left, buffer_recv_right;
|
|---|
| 387 | std::vector<vmg_float>::const_iterator iter;
|
|---|
| 388 |
|
|---|
| 389 | const Index halo_size_1 = grid.Local().HaloEnd1() - grid.Local().HaloBegin1();
|
|---|
| 390 | const Index halo_size_2 = grid.Local().HaloEnd2() - grid.Local().HaloBegin2();
|
|---|
| 391 |
|
|---|
| 392 | buffer_send_left.clear();
|
|---|
| 393 | buffer_send_right.clear();
|
|---|
| 394 | buffer_recv_left.clear();
|
|---|
| 395 | buffer_recv_right.clear();
|
|---|
| 396 |
|
|---|
| 397 | MPI_Cart_shift(comm, 2, 1, &neighbors.Left(), &neighbors.Right());
|
|---|
| 398 |
|
|---|
| 399 | if (halo_size_1.Z()) {
|
|---|
| 400 |
|
|---|
| 401 | buffer_send_left.reserve(grid.Local().Size().X() *
|
|---|
| 402 | grid.Local().Size().Y() *
|
|---|
| 403 | halo_size_1.Z() / 2);
|
|---|
| 404 |
|
|---|
| 405 | for (i.X()=grid.Local().Begin().X(); i.X()<grid.Local().End().X(); ++i.X())
|
|---|
| 406 | for (i.Y()=grid.Local().Begin().Y(); i.Y()<grid.Local().End().Y(); ++i.Y())
|
|---|
| 407 | for (i.Z()=grid.Local().Begin().Z() + (i.X()+i.Y()+offset)%2;
|
|---|
| 408 | i.Z()<grid.Local().Begin().Z() + halo_size_1.Z();
|
|---|
| 409 | i.Z() += 2)
|
|---|
| 410 | buffer_send_left.push_back(grid.GetVal(i));
|
|---|
| 411 | }
|
|---|
| 412 |
|
|---|
| 413 | if (halo_size_2.Z()) {
|
|---|
| 414 |
|
|---|
| 415 | buffer_send_right.reserve(grid.Local().Size().X() *
|
|---|
| 416 | grid.Local().Size().Y() *
|
|---|
| 417 | halo_size_2.Z() / 2);
|
|---|
| 418 |
|
|---|
| 419 | for (i.X()=grid.Local().Begin().X(); i.X()<grid.Local().End().X(); ++i.X())
|
|---|
| 420 | for (i.Y()=grid.Local().Begin().Y(); i.Y()<grid.Local().End().Y(); ++i.Y())
|
|---|
| 421 | for (i.Z()=grid.Local().End().Z() - halo_size_2.Z() + (i.X()+i.Y()+offset)%2;
|
|---|
| 422 | i.Z()<grid.Local().End().Z();
|
|---|
| 423 | i.Z() += 2)
|
|---|
| 424 | buffer_send_right.push_back(grid.GetVal(i));
|
|---|
| 425 |
|
|---|
| 426 | }
|
|---|
| 427 |
|
|---|
| 428 | buffer_recv_left.resize(buffer_send_left.size());
|
|---|
| 429 | buffer_recv_right.resize(buffer_send_right.size());
|
|---|
| 430 |
|
|---|
| 431 | MPI_Sendrecv(&buffer_send_right.front(), buffer_send_right.size(), MPI_DOUBLE, neighbors.Right(), 7,
|
|---|
| 432 | &buffer_recv_left.front(), buffer_recv_left.size(), MPI_DOUBLE, neighbors.Left(), 7,
|
|---|
| 433 | comm, MPI_STATUS_IGNORE);
|
|---|
| 434 |
|
|---|
| 435 | MPI_Sendrecv(&buffer_send_left.front(), buffer_send_left.size(), MPI_DOUBLE, neighbors.Left(), 8,
|
|---|
| 436 | & buffer_recv_right.front(), buffer_recv_right.size(), MPI_DOUBLE, neighbors.Right(), 8,
|
|---|
| 437 | comm, MPI_STATUS_IGNORE);
|
|---|
| 438 |
|
|---|
| 439 | if (halo_size_1.Z()) {
|
|---|
| 440 |
|
|---|
| 441 | iter = buffer_recv_left.begin();
|
|---|
| 442 | for (i.X()=grid.Local().Begin().X(); i.X()<grid.Local().End().X(); ++i.X())
|
|---|
| 443 | for (i.Y()=grid.Local().Begin().Y(); i.Y()<grid.Local().End().Y(); ++i.Y())
|
|---|
| 444 | for (i.Z() = (i.X()+i.Y()+offset)%2; i.Z() < grid.Local().Begin().Z(); i.Z() += 2)
|
|---|
| 445 | grid(i) = *iter++;
|
|---|
| 446 |
|
|---|
| 447 | }
|
|---|
| 448 |
|
|---|
| 449 | if (halo_size_2.Z()) {
|
|---|
| 450 |
|
|---|
| 451 | iter = buffer_recv_right.begin();
|
|---|
| 452 | for (i.X()=grid.Local().Begin().X(); i.X()<grid.Local().End().X(); ++i.X())
|
|---|
| 453 | for (i.Y()=grid.Local().Begin().Y(); i.Y()<grid.Local().End().Y(); ++i.Y())
|
|---|
| 454 | for (i.Z()=grid.Local().End().Z() + (i.X()+i.Y()+offset)%2; i.Z()<grid.Local().SizeTotal().Z(); i.Z() += 2)
|
|---|
| 455 | grid(i) = *iter++;
|
|---|
| 456 |
|
|---|
| 457 | }
|
|---|
| 458 |
|
|---|
| 459 | buffer_send_left.clear();
|
|---|
| 460 | buffer_send_right.clear();
|
|---|
| 461 | buffer_recv_left.clear();
|
|---|
| 462 | buffer_recv_right.clear();
|
|---|
| 463 |
|
|---|
| 464 | MPI_Cart_shift(comm, 1, 1, &neighbors.Left(), &neighbors.Right());
|
|---|
| 465 |
|
|---|
| 466 | if (halo_size_1.Y()) {
|
|---|
| 467 |
|
|---|
| 468 | buffer_send_left.reserve(grid.Local().Size().X() *
|
|---|
| 469 | halo_size_1.Y() *
|
|---|
| 470 | grid.Local().SizeTotal().Z() / 2);
|
|---|
| 471 |
|
|---|
| 472 | for (i.X()=grid.Local().Begin().X(); i.X()<grid.Local().End().X(); ++i.X())
|
|---|
| 473 | for (i.Y()=grid.Local().Begin().Y(); i.Y()<grid.Local().Begin().Y() + halo_size_1.Y(); ++i.Y())
|
|---|
| 474 | for (i.Z()=(i.X()+i.Y()+offset)%2; i.Z()<grid.Local().SizeTotal().Z(); i.Z() += 2)
|
|---|
| 475 | buffer_send_left.push_back(grid.GetVal(i));
|
|---|
| 476 | }
|
|---|
| 477 |
|
|---|
| 478 | if (halo_size_2.Y()) {
|
|---|
| 479 |
|
|---|
| 480 | buffer_send_right.reserve(grid.Local().Size().X() *
|
|---|
| 481 | halo_size_2.Y() *
|
|---|
| 482 | grid.Local().SizeTotal().Z() / 2);
|
|---|
| 483 |
|
|---|
| 484 | for (i.X()=grid.Local().Begin().X(); i.X()<grid.Local().End().X(); ++i.X())
|
|---|
| 485 | for (i.Y()=grid.Local().End().Y() - halo_size_2.Y(); i.Y()<grid.Local().End().Y(); ++i.Y())
|
|---|
| 486 | for (i.Z()=(i.X()+i.Y()+offset)%2; i.Z()<grid.Local().SizeTotal().Z(); i.Z() += 2)
|
|---|
| 487 | buffer_send_right.push_back(grid.GetVal(i));
|
|---|
| 488 |
|
|---|
| 489 | }
|
|---|
| 490 |
|
|---|
| 491 | buffer_recv_left.resize(buffer_send_left.size());
|
|---|
| 492 | buffer_recv_right.resize(buffer_send_right.size());
|
|---|
| 493 |
|
|---|
| 494 | MPI_Sendrecv(&buffer_send_right.front(), buffer_send_right.size(), MPI_DOUBLE, neighbors.Right(), 7,
|
|---|
| 495 | &buffer_recv_left.front(), buffer_recv_left.size(), MPI_DOUBLE, neighbors.Left(), 7,
|
|---|
| 496 | comm, MPI_STATUS_IGNORE);
|
|---|
| 497 |
|
|---|
| 498 | MPI_Sendrecv(&buffer_send_left.front(), buffer_send_left.size(), MPI_DOUBLE, neighbors.Left(), 8,
|
|---|
| 499 | & buffer_recv_right.front(), buffer_recv_right.size(), MPI_DOUBLE, neighbors.Right(), 8,
|
|---|
| 500 | comm, MPI_STATUS_IGNORE);
|
|---|
| 501 |
|
|---|
| 502 | if (halo_size_1.Y()) {
|
|---|
| 503 |
|
|---|
| 504 | iter = buffer_recv_left.begin();
|
|---|
| 505 | for (i.X()=grid.Local().Begin().X(); i.X()<grid.Local().End().X(); ++i.X())
|
|---|
| 506 | for (i.Y()=0; i.Y()<grid.Local().Begin().Y(); ++i.Y())
|
|---|
| 507 | for (i.Z()=(i.X()+i.Y()+offset)%2; i.Z()<grid.Local().SizeTotal().Z(); i.Z() += 2)
|
|---|
| 508 | grid(i) = *iter++;
|
|---|
| 509 |
|
|---|
| 510 | }
|
|---|
| 511 |
|
|---|
| 512 | if (halo_size_2.Y()) {
|
|---|
| 513 |
|
|---|
| 514 | iter = buffer_recv_right.begin();
|
|---|
| 515 | for (i.X()=grid.Local().Begin().X(); i.X()<grid.Local().End().X(); ++i.X())
|
|---|
| 516 | for (i.Y()=grid.Local().End().Y(); i.Y()<grid.Local().SizeTotal().Y(); ++i.Y())
|
|---|
| 517 | for (i.Z()=(i.X()+i.Y()+offset)%2; i.Z()<grid.Local().SizeTotal().Z(); i.Z() += 2)
|
|---|
| 518 | grid(i) = *iter++;
|
|---|
| 519 |
|
|---|
| 520 | }
|
|---|
| 521 |
|
|---|
| 522 | buffer_send_left.clear();
|
|---|
| 523 | buffer_send_right.clear();
|
|---|
| 524 | buffer_recv_left.clear();
|
|---|
| 525 | buffer_recv_right.clear();
|
|---|
| 526 |
|
|---|
| 527 | MPI_Cart_shift(comm, 0, 1, &neighbors.Left(), &neighbors.Right());
|
|---|
| 528 |
|
|---|
| 529 | if (halo_size_1.X()) {
|
|---|
| 530 |
|
|---|
| 531 | buffer_send_left.reserve(halo_size_1.X() *
|
|---|
| 532 | grid.Local().SizeTotal().Y() *
|
|---|
| 533 | grid.Local().SizeTotal().Z() / 2);
|
|---|
| 534 |
|
|---|
| 535 | for (i.X()=grid.Local().Begin().X(); i.X()<grid.Local().Begin().X() + halo_size_1.X(); ++i.X())
|
|---|
| 536 | for (i.Y()=0; i.Y()<grid.Local().SizeTotal().Y(); ++i.Y())
|
|---|
| 537 | for (i.Z()=(i.X()+i.Y()+offset)%2; i.Z()<grid.Local().SizeTotal().Z(); i.Z() += 2)
|
|---|
| 538 | buffer_send_left.push_back(grid.GetVal(i));
|
|---|
| 539 | }
|
|---|
| 540 |
|
|---|
| 541 | if (halo_size_2.X()) {
|
|---|
| 542 |
|
|---|
| 543 | buffer_send_right.reserve(halo_size_2.X() *
|
|---|
| 544 | grid.Local().SizeTotal().Y() *
|
|---|
| 545 | grid.Local().SizeTotal().Z() / 2);
|
|---|
| 546 |
|
|---|
| 547 | for (i.X()=grid.Local().End().X() - halo_size_2.X(); i.X()<grid.Local().End().X(); ++i.X())
|
|---|
| 548 | for (i.Y()=0; i.Y()<grid.Local().SizeTotal().Y(); ++i.Y())
|
|---|
| 549 | for (i.Z()=(i.X()+i.Y()+offset)%2; i.Z()<grid.Local().SizeTotal().Z(); i.Z() += 2)
|
|---|
| 550 | buffer_send_right.push_back(grid.GetVal(i));
|
|---|
| 551 |
|
|---|
| 552 | }
|
|---|
| 553 |
|
|---|
| 554 | buffer_recv_left.resize(buffer_send_left.size());
|
|---|
| 555 | buffer_recv_right.resize(buffer_send_right.size());
|
|---|
| 556 |
|
|---|
| 557 | MPI_Sendrecv(&buffer_send_right.front(), buffer_send_right.size(), MPI_DOUBLE, neighbors.Right(), 7,
|
|---|
| 558 | &buffer_recv_left.front(), buffer_recv_left.size(), MPI_DOUBLE, neighbors.Left(), 7,
|
|---|
| 559 | comm, MPI_STATUS_IGNORE);
|
|---|
| 560 |
|
|---|
| 561 | MPI_Sendrecv(&buffer_send_left.front(), buffer_send_left.size(), MPI_DOUBLE, neighbors.Left(), 8,
|
|---|
| 562 | & buffer_recv_right.front(), buffer_recv_right.size(), MPI_DOUBLE, neighbors.Right(), 8,
|
|---|
| 563 | comm, MPI_STATUS_IGNORE);
|
|---|
| 564 |
|
|---|
| 565 | if (halo_size_1.X()) {
|
|---|
| 566 |
|
|---|
| 567 | iter = buffer_recv_left.begin();
|
|---|
| 568 | for (i.X()=0; i.X()<grid.Local().Begin().X(); ++i.X())
|
|---|
| 569 | for (i.Y()=0; i.Y()<grid.Local().SizeTotal().Y(); ++i.Y())
|
|---|
| 570 | for (i.Z()=(i.X()+i.Y()+offset)%2; i.Z()<grid.Local().SizeTotal().Z(); i.Z() += 2)
|
|---|
| 571 | grid(i) = *iter++;
|
|---|
| 572 | }
|
|---|
| 573 |
|
|---|
| 574 | if (halo_size_2.X()) {
|
|---|
| 575 |
|
|---|
| 576 | iter = buffer_recv_right.begin();
|
|---|
| 577 | for (i.X()=grid.Local().End().X(); i.X()<grid.Local().SizeTotal().X(); ++i.X())
|
|---|
| 578 | for (i.Y()=0; i.Y()<grid.Local().SizeTotal().Y(); ++i.Y())
|
|---|
| 579 | for (i.Z()=(i.X()+i.Y()+offset)%2; i.Z()<grid.Local().SizeTotal().Z(); i.Z() += 2)
|
|---|
| 580 | grid(i) = *iter++;
|
|---|
| 581 |
|
|---|
| 582 | }
|
|---|
| 583 |
|
|---|
| 584 | }
|
|---|
| 585 | }
|
|---|
| 586 |
|
|---|
| 587 | Grid& CommMPI::GetCoarserGrid(Multigrid& multigrid)
|
|---|
| 588 | {
|
|---|
| 589 | std::map<const Grid*, TempGrid*>::iterator iter = coarser_grids.find(&multigrid());
|
|---|
| 590 |
|
|---|
| 591 | if (iter == coarser_grids.end()) {
|
|---|
| 592 | TempGrid* temp_grid = new TempGrid();
|
|---|
| 593 | temp_grid->SetPropertiesToCoarser(multigrid(), BoundaryConditions());
|
|---|
| 594 | iter = coarser_grids.insert(std::make_pair(&multigrid(), temp_grid)).first;
|
|---|
| 595 | }
|
|---|
| 596 |
|
|---|
| 597 | TempGrid& grid_coarser = *iter->second;
|
|---|
| 598 |
|
|---|
| 599 | CommSubgrid(multigrid(multigrid.Level()-1), grid_coarser);
|
|---|
| 600 |
|
|---|
| 601 | return grid_coarser;
|
|---|
| 602 | }
|
|---|
| 603 |
|
|---|
| 604 | Grid& CommMPI::GetFinerGrid(Multigrid& multigrid)
|
|---|
| 605 | {
|
|---|
| 606 | std::map<const Grid*, TempGrid*>::iterator iter = finer_grids.find(&multigrid());
|
|---|
| 607 |
|
|---|
| 608 | if (iter == finer_grids.end()) {
|
|---|
| 609 | TempGrid* temp_grid = new TempGrid();
|
|---|
| 610 | temp_grid->SetPropertiesToFiner(multigrid(), BoundaryConditions());
|
|---|
| 611 | iter = finer_grids.insert(std::make_pair(&multigrid(), temp_grid)).first;
|
|---|
| 612 | }
|
|---|
| 613 |
|
|---|
| 614 | TempGrid& grid_finer = *iter->second;
|
|---|
| 615 |
|
|---|
| 616 | CommSubgrid(multigrid(multigrid.Level()+1), grid_finer);
|
|---|
| 617 |
|
|---|
| 618 | return grid_finer;
|
|---|
| 619 | }
|
|---|
| 620 |
|
|---|
| 621 | Grid& CommMPI::GetGlobalCoarseGrid(Multigrid& multigrid)
|
|---|
| 622 | {
|
|---|
| 623 | if (global_coarse_grid == NULL) {
|
|---|
| 624 |
|
|---|
| 625 | global_coarse_grid = new TempGrid();
|
|---|
| 626 | global_coarse_grid->SetPropertiesToGlobalCoarseGrid();
|
|---|
| 627 |
|
|---|
| 628 | }
|
|---|
| 629 |
|
|---|
| 630 | CommSubgrid(multigrid(multigrid.MinLevel()), *global_coarse_grid);
|
|---|
| 631 |
|
|---|
| 632 | return *global_coarse_grid;
|
|---|
| 633 | }
|
|---|
| 634 |
|
|---|
| 635 | //
|
|---|
| 636 | // Decodes a derived datatype created by MPI_Type_create_subarray
|
|---|
| 637 | // and adds values to grid.
|
|---|
| 638 | // I chose this approach since MPI neither supports a point-to-point reduce nor
|
|---|
| 639 | // a non-blocking global reduce.
|
|---|
| 640 | // For more information see Chapter 4.1 "DERIVED DATATYPES" of the MPI standard
|
|---|
| 641 | //
|
|---|
| 642 | void CommMPI::ReceiveAndAddSubgrid(Grid& grid, const MPI_Comm& comm, MPI_Datatype& type,
|
|---|
| 643 | const int& rank, const int& tag)
|
|---|
| 644 | {
|
|---|
| 645 | int counter = 0;
|
|---|
| 646 | int ibuffer[11];
|
|---|
| 647 | Index i;
|
|---|
| 648 | MPI_Datatype base_type;
|
|---|
| 649 | MPI_Aint address;
|
|---|
| 650 |
|
|---|
| 651 | MPI_Type_get_contents(type, 11, 0, 1, ibuffer, &address, &base_type);
|
|---|
| 652 |
|
|---|
| 653 | if (static_cast<int>(receive_buffer.size()) < ibuffer[4]*ibuffer[5]*ibuffer[6])
|
|---|
| 654 | receive_buffer.resize(ibuffer[4]*ibuffer[5]*ibuffer[6]);
|
|---|
| 655 |
|
|---|
| 656 | MPI_Recv(&(receive_buffer.front()), ibuffer[4]*ibuffer[5]*ibuffer[6], MPI_DOUBLE, rank, tag, comm, MPI_STATUS_IGNORE);
|
|---|
| 657 |
|
|---|
| 658 | for (i.X()=ibuffer[7]; i.X()<ibuffer[7]+ibuffer[4]; ++i.X())
|
|---|
| 659 | for (i.Y()=ibuffer[8]; i.Y()<ibuffer[8]+ibuffer[5]; ++i.Y())
|
|---|
| 660 | for (i.Z()=ibuffer[9]; i.Z()<ibuffer[9]+ibuffer[6]; ++i.Z())
|
|---|
| 661 | grid(i) += receive_buffer[counter++];
|
|---|
| 662 | }
|
|---|
| 663 |
|
|---|
| 664 | void CommMPI::CommParticles(const Grid& grid, std::list<Particle::Particle>& particles)
|
|---|
| 665 | {
|
|---|
| 666 | Timer::Start("CommParticles");
|
|---|
| 667 |
|
|---|
| 668 | Factory& factory = MG::GetFactory();
|
|---|
| 669 |
|
|---|
| 670 | const vmg_int& num_particles_local = factory.GetObjectStorageVal<vmg_int>("PARTICLE_NUM_LOCAL");
|
|---|
| 671 | vmg_float* x = factory.GetObjectStorageArray<vmg_float>("PARTICLE_POS_ARRAY");
|
|---|
| 672 | vmg_float* q = factory.GetObjectStorageArray<vmg_float>("PARTICLE_CHARGE_ARRAY");
|
|---|
| 673 |
|
|---|
| 674 | int rank, size;
|
|---|
| 675 | MPI_Comm_rank(comm_global, &rank);
|
|---|
| 676 | MPI_Comm_size(comm_global, &size);
|
|---|
| 677 |
|
|---|
| 678 | Index index;
|
|---|
| 679 | int global_extent[6*size], send_sizes[size], recv_sizes[size];
|
|---|
| 680 | Index begin_remote[size], end_remote[size];
|
|---|
| 681 | std::vector<int> pos_indices[size], charge_indices[size];
|
|---|
| 682 | MPI_Datatype dt_temp;
|
|---|
| 683 | std::vector<vmg_float> buffer_x, buffer_q;
|
|---|
| 684 | std::vector<vmg_int> buffer_ind;
|
|---|
| 685 |
|
|---|
| 686 | std::memcpy(&global_extent[6*rank], grid.Global().BeginLocal().vec(), 3*sizeof(int));
|
|---|
| 687 | std::memcpy(&global_extent[6*rank+3], grid.Global().EndLocal().vec(), 3*sizeof(int));
|
|---|
| 688 |
|
|---|
| 689 | MPI_Allgather(MPI_IN_PLACE, 6, MPI_INT, global_extent, 6, MPI_INT, MPI_COMM_WORLD);
|
|---|
| 690 |
|
|---|
| 691 | for (int i=0; i<size; ++i) {
|
|---|
| 692 | begin_remote[i] = static_cast<Vector>(&global_extent[6*i]);
|
|---|
| 693 | end_remote[i] = static_cast<Vector>(&global_extent[6*i+3]);
|
|---|
| 694 | }
|
|---|
| 695 |
|
|---|
| 696 | for (int i=0; i<num_particles_local; ++i) {
|
|---|
| 697 | index = static_cast<Index>((Vector(&x[3*i]) - grid.Extent().Begin()) / grid.Extent().MeshWidth());
|
|---|
| 698 | for (int j=0; j<size; ++j)
|
|---|
| 699 | if (index.IsInBounds(begin_remote[j], end_remote[j])) {
|
|---|
| 700 | pos_indices[j].push_back(3*i);
|
|---|
| 701 | charge_indices[j].push_back(i);
|
|---|
| 702 | break;
|
|---|
| 703 | }
|
|---|
| 704 | }
|
|---|
| 705 |
|
|---|
| 706 | /*
|
|---|
| 707 | * Communicate which process gets how many particles
|
|---|
| 708 | */
|
|---|
| 709 | for (int i=0; i<size; ++i)
|
|---|
| 710 | send_sizes[i] = charge_indices[i].size();
|
|---|
| 711 |
|
|---|
| 712 | MPI_Alltoall(send_sizes, 1, MPI_INT, recv_sizes, 1, MPI_INT, comm_global);
|
|---|
| 713 |
|
|---|
| 714 | /*
|
|---|
| 715 | * Send particles
|
|---|
| 716 | */
|
|---|
| 717 | for (int i=0; i<size; ++i) {
|
|---|
| 718 |
|
|---|
| 719 | if (pos_indices[i].size() > 0) {
|
|---|
| 720 |
|
|---|
| 721 | MPI_Type_create_indexed_block(pos_indices[i].size(), 3, &pos_indices[i].front(), MPI_DOUBLE, &dt_temp);
|
|---|
| 722 | MPI_Type_commit(&dt_temp);
|
|---|
| 723 | MPI_Isend(x, 1, dt_temp, i, 0, MPI_COMM_WORLD, &Request());
|
|---|
| 724 | MPI_Type_free(&dt_temp);
|
|---|
| 725 |
|
|---|
| 726 | MPI_Type_create_indexed_block(charge_indices[i].size(), 1, &charge_indices[i].front(), MPI_DOUBLE, &dt_temp);
|
|---|
| 727 | MPI_Type_commit(&dt_temp);
|
|---|
| 728 | MPI_Isend(q, 1, dt_temp, i, 1, MPI_COMM_WORLD, &Request());
|
|---|
| 729 | MPI_Type_free(&dt_temp);
|
|---|
| 730 |
|
|---|
| 731 | MPI_Isend(&charge_indices[i].front(), charge_indices[i].size(), MPI_INT, i, 2, MPI_COMM_WORLD, &Request());
|
|---|
| 732 |
|
|---|
| 733 | }
|
|---|
| 734 |
|
|---|
| 735 | }
|
|---|
| 736 |
|
|---|
| 737 | /*
|
|---|
| 738 | * Receive particles
|
|---|
| 739 | */
|
|---|
| 740 | for (int i=0; i<size; ++i) {
|
|---|
| 741 |
|
|---|
| 742 | if (recv_sizes[i] > 0) {
|
|---|
| 743 |
|
|---|
| 744 | buffer_x.resize(3*recv_sizes[i]);
|
|---|
| 745 | buffer_q.resize(recv_sizes[i]);
|
|---|
| 746 | buffer_ind.resize(recv_sizes[i]);
|
|---|
| 747 |
|
|---|
| 748 | MPI_Recv(&buffer_x.front(), 3*recv_sizes[i], MPI_DOUBLE, i, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
|
|---|
| 749 | MPI_Recv(&buffer_q.front(), recv_sizes[i], MPI_DOUBLE, i, 1, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
|
|---|
| 750 | MPI_Recv(&buffer_ind.front(), recv_sizes[i], MPI_INT, i, 2, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
|
|---|
| 751 |
|
|---|
| 752 | particles.clear();
|
|---|
| 753 |
|
|---|
| 754 | for (int j=0; j<recv_sizes[i]; ++j)
|
|---|
| 755 | particles.push_back(Particle::Particle(&buffer_x[3*j], buffer_q[j], 0.0, 0.0, i, buffer_ind[j]));
|
|---|
| 756 |
|
|---|
| 757 | }
|
|---|
| 758 |
|
|---|
| 759 | }
|
|---|
| 760 |
|
|---|
| 761 | ClearRequests();
|
|---|
| 762 |
|
|---|
| 763 | Timer::Stop("CommParticles");
|
|---|
| 764 | }
|
|---|
| 765 |
|
|---|
| 766 | void CommMPI::CommLCListGhosts(const Grid& grid, Particle::LinkedCellList& lc)
|
|---|
| 767 | {
|
|---|
| 768 | Timer::Start("CommLCListGhosts");
|
|---|
| 769 |
|
|---|
| 770 | std::vector<vmg_float> send_buffer_1, send_buffer_2, recv_buffer;
|
|---|
| 771 | vmg_int send_size_1, send_size_2, recv_size;
|
|---|
| 772 | std::list<Particle::Particle>::iterator p_iter;
|
|---|
| 773 | Grid::iterator g_iter;
|
|---|
| 774 |
|
|---|
| 775 | Tuple<int> neighbors;
|
|---|
| 776 |
|
|---|
| 777 | for (int i=2; i>=0; --i) {
|
|---|
| 778 |
|
|---|
| 779 | /*
|
|---|
| 780 | * Get ranks of neighbors
|
|---|
| 781 | */
|
|---|
| 782 | MPI_Cart_shift(comm_global, i, 1, &neighbors.Left(), &neighbors.Right());
|
|---|
| 783 |
|
|---|
| 784 | /*
|
|---|
| 785 | * Send to left
|
|---|
| 786 | */
|
|---|
| 787 | send_buffer_1.clear();
|
|---|
| 788 |
|
|---|
| 789 | for (g_iter=lc.Iterators().Boundary1()[i].Begin(); g_iter!=lc.Iterators().Boundary1()[i].End(); ++g_iter)
|
|---|
| 790 | lc(*g_iter).clear();
|
|---|
| 791 |
|
|---|
| 792 | for (g_iter=lc.Iterators().NearBoundary1()[i].Begin(); g_iter!=lc.Iterators().NearBoundary1()[i].End(); ++g_iter)
|
|---|
| 793 | for (p_iter=lc(*g_iter).begin(); p_iter!=lc(*g_iter).end(); ++p_iter) {
|
|---|
| 794 |
|
|---|
| 795 | if (grid.Global().BeginLocal()[i] == 0)
|
|---|
| 796 | for (int j=0; j<3; ++j)
|
|---|
| 797 | send_buffer_1.push_back(p_iter->Pos()[j] + grid.Extent().Size()[j]);
|
|---|
| 798 | else
|
|---|
| 799 | for (int j=0; j<3; ++j)
|
|---|
| 800 | send_buffer_1.push_back(p_iter->Pos()[j]);
|
|---|
| 801 |
|
|---|
| 802 | send_buffer_1.push_back(p_iter->Charge());
|
|---|
| 803 |
|
|---|
| 804 | }
|
|---|
| 805 |
|
|---|
| 806 | send_size_1 = send_buffer_1.size();
|
|---|
| 807 |
|
|---|
| 808 | MPI_Isend(&send_size_1, 1, MPI_INT, neighbors.Left(), 0, comm_global, &Request());
|
|---|
| 809 |
|
|---|
| 810 | if (send_size_1 > 0)
|
|---|
| 811 | MPI_Isend(&send_buffer_1.front(), send_size_1, MPI_DOUBLE, neighbors.Left(), 1, comm_global, &Request());
|
|---|
| 812 |
|
|---|
| 813 | /*
|
|---|
| 814 | * Send to right
|
|---|
| 815 | */
|
|---|
| 816 | send_buffer_2.clear();
|
|---|
| 817 |
|
|---|
| 818 | for (g_iter=lc.Iterators().Boundary2()[i].Begin(); g_iter!=lc.Iterators().Boundary2()[i].End(); ++g_iter)
|
|---|
| 819 | lc(*g_iter).clear();
|
|---|
| 820 |
|
|---|
| 821 | for (g_iter=lc.Iterators().NearBoundary2()[i].Begin(); g_iter!=lc.Iterators().NearBoundary2()[i].End(); ++g_iter)
|
|---|
| 822 | for (p_iter=lc(*g_iter).begin(); p_iter!=lc(*g_iter).end(); ++p_iter) {
|
|---|
| 823 |
|
|---|
| 824 | if (grid.Global().EndLocal()[i] == grid.Global().SizeGlobal()[i])
|
|---|
| 825 | for (int j=0; j<3; ++j)
|
|---|
| 826 | send_buffer_2.push_back(p_iter->Pos()[j] - grid.Extent().Size()[j]);
|
|---|
| 827 | else
|
|---|
| 828 | for (int j=0; j<3; ++j)
|
|---|
| 829 | send_buffer_2.push_back(p_iter->Pos()[j]);
|
|---|
| 830 |
|
|---|
| 831 | send_buffer_2.push_back(p_iter->Charge());
|
|---|
| 832 |
|
|---|
| 833 | }
|
|---|
| 834 |
|
|---|
| 835 | send_size_2 = send_buffer_2.size();
|
|---|
| 836 |
|
|---|
| 837 | MPI_Isend(&send_size_2, 1, MPI_INT, neighbors.Right(), 2, comm_global, &Request());
|
|---|
| 838 |
|
|---|
| 839 | if (send_size_2 > 0)
|
|---|
| 840 | MPI_Isend(&send_buffer_2.front(), send_size_2, MPI_DOUBLE, neighbors.Right(), 3, comm_global, &Request());
|
|---|
| 841 |
|
|---|
| 842 | /*
|
|---|
| 843 | * Receive from right
|
|---|
| 844 | */
|
|---|
| 845 | MPI_Recv(&recv_size, 1, MPI_INT, neighbors.Right(), 0, comm_global, MPI_STATUS_IGNORE);
|
|---|
| 846 |
|
|---|
| 847 | if (recv_size > 0) {
|
|---|
| 848 |
|
|---|
| 849 | recv_buffer.resize(recv_size);
|
|---|
| 850 | MPI_Recv(&recv_buffer.front(), recv_size, MPI_DOUBLE, neighbors.Right(), 1, comm_global, MPI_STATUS_IGNORE);
|
|---|
| 851 |
|
|---|
| 852 | for (vmg_int i=0; i<recv_size; i+=4)
|
|---|
| 853 | lc.AddParticleToComplete(&recv_buffer[i], recv_buffer[i+3]);
|
|---|
| 854 |
|
|---|
| 855 | }
|
|---|
| 856 |
|
|---|
| 857 | /*
|
|---|
| 858 | * Receive from left
|
|---|
| 859 | */
|
|---|
| 860 | MPI_Recv(&recv_size, 1, MPI_INT, neighbors.Left(), 2, comm_global, MPI_STATUS_IGNORE);
|
|---|
| 861 |
|
|---|
| 862 | if (recv_size > 0) {
|
|---|
| 863 |
|
|---|
| 864 | recv_buffer.resize(recv_size);
|
|---|
| 865 | MPI_Recv(&recv_buffer.front(), recv_size, MPI_DOUBLE, neighbors.Left(), 3, comm_global, MPI_STATUS_IGNORE);
|
|---|
| 866 |
|
|---|
| 867 | for (vmg_int i=0; i<recv_size; i+=4)
|
|---|
| 868 | lc.AddParticleToComplete(&recv_buffer[i], recv_buffer[i+3]);
|
|---|
| 869 |
|
|---|
| 870 | }
|
|---|
| 871 |
|
|---|
| 872 |
|
|---|
| 873 | ClearRequests();
|
|---|
| 874 | }
|
|---|
| 875 |
|
|---|
| 876 | Timer::Stop("CommLCListGhosts");
|
|---|
| 877 | }
|
|---|
| 878 |
|
|---|
| 879 | void CommMPI::CommAddPotential(std::list<Particle::Particle>& particles)
|
|---|
| 880 | {
|
|---|
| 881 | Timer::Start("CommAddPotential");
|
|---|
| 882 |
|
|---|
| 883 | std::list<Particle::Particle>::iterator iter;
|
|---|
| 884 |
|
|---|
| 885 | if (!win_created) {
|
|---|
| 886 | vmg_float* p = MG::GetFactory().GetObjectStorageArray<vmg_float>("PARTICLE_POTENTIAL_ARRAY");
|
|---|
| 887 | const vmg_int& num_particles_local = MG::GetFactory().GetObjectStorageVal<vmg_int>("PARTICLE_NUM_LOCAL");
|
|---|
| 888 | MPI_Win_create(p, num_particles_local*sizeof(vmg_float), sizeof(vmg_float), info, comm_global, &win);
|
|---|
| 889 | win_created = true;
|
|---|
| 890 | }
|
|---|
| 891 |
|
|---|
| 892 | MPI_Win_fence(MPI_MODE_NOPRECEDE, win);
|
|---|
| 893 |
|
|---|
| 894 | for (iter=particles.begin(); iter!=particles.end(); ++iter)
|
|---|
| 895 | MPI_Accumulate(&iter->Pot(), 1, MPI_DOUBLE, iter->Rank(), iter->Index(), 1, MPI_DOUBLE, MPI_SUM, win);
|
|---|
| 896 |
|
|---|
| 897 | MPI_Win_fence(MPI_MODE_NOSTORE | MPI_MODE_NOSUCCEED, win);
|
|---|
| 898 |
|
|---|
| 899 | Timer::Stop("CommAddPotential");
|
|---|
| 900 | }
|
|---|
| 901 |
|
|---|
| 902 | void CommMPI::CommAddPotential(Particle::LinkedCellList& lc)
|
|---|
| 903 | {
|
|---|
| 904 | Timer::Start("CommAddPotential");
|
|---|
| 905 |
|
|---|
| 906 | Grid::iterator g_iter;
|
|---|
| 907 | std::list<Particle::Particle>::iterator p_iter;
|
|---|
| 908 |
|
|---|
| 909 | if (!win_created) {
|
|---|
| 910 | vmg_float* p = MG::GetFactory().GetObjectStorageArray<vmg_float>("PARTICLE_POTENTIAL_ARRAY");
|
|---|
| 911 | const vmg_int& num_particles_local = MG::GetFactory().GetObjectStorageVal<vmg_int>("PARTICLE_NUM_LOCAL");
|
|---|
| 912 | MPI_Win_create(p, num_particles_local*sizeof(vmg_float), sizeof(vmg_float), info, comm_global, &win);
|
|---|
| 913 | win_created = true;
|
|---|
| 914 | }
|
|---|
| 915 |
|
|---|
| 916 | MPI_Win_fence(MPI_MODE_NOPRECEDE, win);
|
|---|
| 917 |
|
|---|
| 918 | for (g_iter=lc.Iterators().Local().Begin(); g_iter!=lc.Iterators().Local().End(); ++g_iter)
|
|---|
| 919 | for (p_iter=lc(*g_iter).begin(); p_iter!=lc(*g_iter).end(); ++p_iter)
|
|---|
| 920 | MPI_Accumulate(&p_iter->Pot(), 1, MPI_DOUBLE, p_iter->Rank(), p_iter->Index(), 1, MPI_DOUBLE, MPI_SUM, win);
|
|---|
| 921 |
|
|---|
| 922 | MPI_Win_fence(MPI_MODE_NOSTORE | MPI_MODE_NOSUCCEED, win);
|
|---|
| 923 |
|
|---|
| 924 | Timer::Stop("CommAddPotential");
|
|---|
| 925 | }
|
|---|
| 926 |
|
|---|
| 927 | vmg_float CommMPI::GlobalSum(const vmg_float& value)
|
|---|
| 928 | {
|
|---|
| 929 | vmg_float result = value;
|
|---|
| 930 |
|
|---|
| 931 | MPI_Allreduce(MPI_IN_PLACE, &result, 1, MPI_DOUBLE, MPI_SUM, comm_global);
|
|---|
| 932 |
|
|---|
| 933 | return result;
|
|---|
| 934 | }
|
|---|
| 935 |
|
|---|
| 936 | vmg_float CommMPI::GlobalSumRoot(const vmg_float& value)
|
|---|
| 937 | {
|
|---|
| 938 | vmg_float recv_buffer;
|
|---|
| 939 | vmg_float send_buffer = value;
|
|---|
| 940 |
|
|---|
| 941 | MPI_Reduce(&send_buffer, &recv_buffer, 1, MPI_DOUBLE, MPI_SUM, 0, comm_global);
|
|---|
| 942 |
|
|---|
| 943 | return recv_buffer;
|
|---|
| 944 | }
|
|---|
| 945 |
|
|---|
| 946 | void CommMPI::GlobalSumArray(vmg_float* array, const vmg_int& size)
|
|---|
| 947 | {
|
|---|
| 948 | MPI_Allreduce(MPI_IN_PLACE, array, size, MPI_DOUBLE, MPI_SUM, comm_global);
|
|---|
| 949 | }
|
|---|
| 950 |
|
|---|
| 951 | void CommMPI::PrintString(const char* format, ...)
|
|---|
| 952 | {
|
|---|
| 953 | va_list args;
|
|---|
| 954 | va_start(args, format);
|
|---|
| 955 | vsprintf(print_buffer, format, args);
|
|---|
| 956 | printf("VMG: Rank %d: %s\n", GlobalRank(), print_buffer);
|
|---|
| 957 | va_end(args);
|
|---|
| 958 | }
|
|---|
| 959 |
|
|---|
| 960 | void CommMPI::PrintStringOnce(const char* format, ...)
|
|---|
| 961 | {
|
|---|
| 962 | if (GlobalRank() == 0) {
|
|---|
| 963 | va_list args;
|
|---|
| 964 | va_start(args, format);
|
|---|
| 965 | vsprintf(print_buffer, format, args);
|
|---|
| 966 | printf("VMG: Rank %d: %s\n", GlobalRank(), print_buffer);
|
|---|
| 967 | va_end(args);
|
|---|
| 968 | }
|
|---|
| 969 | }
|
|---|
| 970 |
|
|---|
| 971 | void CommMPI::PrintXML(const std::string& filename, const std::string& xml_data)
|
|---|
| 972 | {
|
|---|
| 973 | MPI_File file;
|
|---|
| 974 | std::stringstream path, xml_header;
|
|---|
| 975 |
|
|---|
| 976 | pugi::xml_document doc;
|
|---|
| 977 | pugi::xml_node node_data = doc.append_child("Global").append_child("NumProcesses").append_child(pugi::node_pcdata);
|
|---|
| 978 | node_data.set_value(Helper::ToString(GlobalProcs().Product()).c_str());
|
|---|
| 979 | doc.save(xml_header);
|
|---|
| 980 |
|
|---|
| 981 | path << OutputPath() << filename;
|
|---|
| 982 |
|
|---|
| 983 | char* filename_array = Helper::GetCharArray(path.str());
|
|---|
| 984 | char* xml_header_array = Helper::GetCharArray(xml_header.str());
|
|---|
| 985 | char* str_array = Helper::GetCharArray(xml_data);
|
|---|
| 986 |
|
|---|
| 987 | MPI_File_open(MPI_COMM_SELF, filename_array, MPI_MODE_WRONLY|MPI_MODE_CREATE, MPI_INFO_NULL, &file);
|
|---|
| 988 | MPI_File_set_size(file, 0);
|
|---|
| 989 | MPI_File_write(file, xml_header_array, xml_header.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
|
|---|
| 990 | MPI_File_write(file, str_array, xml_data.size(), MPI_CHAR, MPI_STATUS_IGNORE);
|
|---|
| 991 | MPI_File_close(&file);
|
|---|
| 992 |
|
|---|
| 993 | delete [] filename_array;
|
|---|
| 994 | delete [] xml_header_array;
|
|---|
| 995 | delete [] str_array;
|
|---|
| 996 | }
|
|---|
| 997 |
|
|---|
| 998 | void CommMPI::PrintXMLAll(const std::string& filename, const std::string& xml_data)
|
|---|
| 999 | {
|
|---|
| 1000 | MPI_File file;
|
|---|
| 1001 | std::stringstream path;
|
|---|
| 1002 |
|
|---|
| 1003 | path << OutputPath() << filename;
|
|---|
| 1004 |
|
|---|
| 1005 | char* filename_array = Helper::GetCharArray(path.str());
|
|---|
| 1006 | char* str_array = Helper::GetCharArray(xml_data);
|
|---|
| 1007 |
|
|---|
| 1008 | MPI_File_open(comm_global, filename_array, MPI_MODE_WRONLY|MPI_MODE_CREATE, MPI_INFO_NULL, &file);
|
|---|
| 1009 | MPI_File_set_size(file, 0);
|
|---|
| 1010 |
|
|---|
| 1011 | if (GlobalRank() == 0) {
|
|---|
| 1012 | std::stringstream xml_header;
|
|---|
| 1013 | pugi::xml_document doc;
|
|---|
| 1014 | pugi::xml_node node_data = doc.append_child("Global").append_child("NumProcesses").append_child(pugi::node_pcdata);
|
|---|
| 1015 | node_data.set_value(Helper::ToString(GlobalProcs().Product()).c_str());
|
|---|
| 1016 | doc.save(xml_header);
|
|---|
| 1017 |
|
|---|
| 1018 | char* xml_header_array = Helper::GetCharArray(xml_header.str());
|
|---|
| 1019 |
|
|---|
| 1020 | MPI_File_write_shared(file, xml_header_array, xml_header.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
|
|---|
| 1021 |
|
|---|
| 1022 | delete [] xml_header_array;
|
|---|
| 1023 | }
|
|---|
| 1024 |
|
|---|
| 1025 | MPI_File_write_ordered(file, str_array, xml_data.size(), MPI_CHAR, MPI_STATUS_IGNORE);
|
|---|
| 1026 | MPI_File_close(&file);
|
|---|
| 1027 |
|
|---|
| 1028 | delete [] filename_array;
|
|---|
| 1029 | delete [] str_array;
|
|---|
| 1030 | }
|
|---|
| 1031 |
|
|---|
| 1032 | void CommMPI::PrintAllSettings()
|
|---|
| 1033 | {
|
|---|
| 1034 | std::stringstream buf;
|
|---|
| 1035 | MPI_File file;
|
|---|
| 1036 | MPI_Status status;
|
|---|
| 1037 |
|
|---|
| 1038 | Multigrid* mg = MG::GetFactory().Get("SOL")->Cast<Multigrid>();
|
|---|
| 1039 |
|
|---|
| 1040 | buf << OutputPath() << "settings.txt";
|
|---|
| 1041 | char *filename = Helper::GetCharArray(buf.str());
|
|---|
| 1042 |
|
|---|
| 1043 | MPI_File_open(comm_global, filename, MPI_MODE_WRONLY|MPI_MODE_CREATE, MPI_INFO_NULL, &file);
|
|---|
| 1044 | MPI_File_set_size(file, 0);
|
|---|
| 1045 | MPI_File_close(&file);
|
|---|
| 1046 |
|
|---|
| 1047 | for (int i=mg->MinLevel(); i<=mg->MaxLevel(); ++i) {
|
|---|
| 1048 |
|
|---|
| 1049 | MPI_Comm comm = comm_info.GetCommunicator((*mg)(i));
|
|---|
| 1050 |
|
|---|
| 1051 | if (comm != MPI_COMM_NULL) {
|
|---|
| 1052 |
|
|---|
| 1053 | MPI_File_open(comm, filename, MPI_MODE_WRONLY|MPI_MODE_CREATE|MPI_MODE_APPEND, MPI_INFO_NULL, &file);
|
|---|
| 1054 |
|
|---|
| 1055 | int comm_rank;
|
|---|
| 1056 | MPI_Comm_rank(comm, &comm_rank);
|
|---|
| 1057 |
|
|---|
| 1058 | if (comm_rank == 0) {
|
|---|
| 1059 |
|
|---|
| 1060 | buf.str("");
|
|---|
| 1061 | buf << "PROCESS INDEPENDENT INFORMATION:" << std::endl
|
|---|
| 1062 | << "GRID LEVEL " << i << std::endl
|
|---|
| 1063 | << "BEGINFINEST: " << (*mg)(i).Global().BeginFinest() << std::endl
|
|---|
| 1064 | << "ENDFINEST: " << (*mg)(i).Global().EndFinest() << std::endl
|
|---|
| 1065 | << "SIZEFINEST: " << (*mg)(i).Global().SizeFinest() << std::endl
|
|---|
| 1066 | << "SIZEGLOBAL: " << (*mg)(i).Global().SizeGlobal() << std::endl
|
|---|
| 1067 | << "BOUNDARYTYPE: " << (*mg)(i).Global().BoundaryType() << std::endl
|
|---|
| 1068 | << "EXTENT" << std::endl
|
|---|
| 1069 | << "SIZE: " << (*mg)(i).Extent().Size() << std::endl
|
|---|
| 1070 | << "SIZEFACTOR: " << (*mg)(i).Extent().SizeFactor() << std::endl
|
|---|
| 1071 | << "BEGIN: " << (*mg)(i).Extent().Begin() << std::endl
|
|---|
| 1072 | << "END: " << (*mg)(i).Extent().End() << std::endl
|
|---|
| 1073 | << "MESHWIDTH: " << (*mg)(i).Extent().MeshWidth() << std::endl
|
|---|
| 1074 | << std::endl
|
|---|
| 1075 | << "PROCESS DEPENDENT INFORMATION:" << std::endl;
|
|---|
| 1076 |
|
|---|
| 1077 | char* char_buf = Helper::GetCharArray(buf.str());
|
|---|
| 1078 |
|
|---|
| 1079 | MPI_File_write_shared(file, char_buf, buf.str().size(), MPI_CHAR, &status);
|
|---|
| 1080 |
|
|---|
| 1081 | delete [] char_buf;
|
|---|
| 1082 |
|
|---|
| 1083 | }
|
|---|
| 1084 |
|
|---|
| 1085 | buf.str("");
|
|---|
| 1086 | buf << "PROCESS " << comm_rank << ":" << std::endl
|
|---|
| 1087 | << "GLOBAL" << std::endl
|
|---|
| 1088 | << "BEGINLOCAL: " << (*mg)(i).Global().BeginLocal() << std::endl
|
|---|
| 1089 | << "ENDLOCAL: " << (*mg)(i).Global().EndLocal() << std::endl
|
|---|
| 1090 | << "SIZELOCAL: " << (*mg)(i).Global().SizeLocal() << std::endl
|
|---|
| 1091 | << "LOCAL" << std::endl
|
|---|
| 1092 | << "SIZE: " << (*mg)(i).Local().Size() << std::endl
|
|---|
| 1093 | << "SIZETOTAL: " << (*mg)(i).Local().SizeTotal() << std::endl
|
|---|
| 1094 | << "BEGIN: " << (*mg)(i).Local().Begin() << std::endl
|
|---|
| 1095 | << "END: " << (*mg)(i).Local().End() << std::endl
|
|---|
| 1096 | << "HALOBEGIN1: " << (*mg)(i).Local().HaloBegin1() << std::endl
|
|---|
| 1097 | << "HALOEND1: " << (*mg)(i).Local().HaloEnd1() << std::endl
|
|---|
| 1098 | << "HALOBEGIN2: " << (*mg)(i).Local().HaloBegin2() << std::endl
|
|---|
| 1099 | << "HALOEND2: " << (*mg)(i).Local().HaloEnd2() << std::endl
|
|---|
| 1100 | << "BOUNDARYBEGIN1: " << (*mg)(i).Local().BoundaryBegin1() << std::endl
|
|---|
| 1101 | << "BOUNDARYEND1: " << (*mg)(i).Local().BoundaryEnd1() << std::endl
|
|---|
| 1102 | << "BOUNDARYBEGIN2: " << (*mg)(i).Local().BoundaryBegin2() << std::endl
|
|---|
| 1103 | << "BOUNDARYEND2: " << (*mg)(i).Local().BoundaryEnd2() << std::endl
|
|---|
| 1104 | << "ALIGNMENTBEGIN: " << (*mg)(i).Local().AlignmentBegin() << std::endl
|
|---|
| 1105 | << "ALIGNMENTEND: " << (*mg)(i).Local().AlignmentEnd() << std::endl
|
|---|
| 1106 | << std::endl;
|
|---|
| 1107 |
|
|---|
| 1108 | char* char_buf = Helper::GetCharArray(buf.str());
|
|---|
| 1109 |
|
|---|
| 1110 | MPI_File_write_ordered(file, char_buf, buf.str().size(), MPI_CHAR, &status);
|
|---|
| 1111 |
|
|---|
| 1112 | delete [] char_buf;
|
|---|
| 1113 |
|
|---|
| 1114 | MPI_File_close(&file);
|
|---|
| 1115 |
|
|---|
| 1116 | }
|
|---|
| 1117 | }
|
|---|
| 1118 |
|
|---|
| 1119 | delete [] filename;
|
|---|
| 1120 |
|
|---|
| 1121 | }
|
|---|
| 1122 |
|
|---|
| 1123 | void CommMPI::PrintDefect(Grid& sol, Grid& rhs, const char* information)
|
|---|
| 1124 | {
|
|---|
| 1125 | TempGrid *temp = MG::GetTempGrid();
|
|---|
| 1126 | temp->SetProperties(sol);
|
|---|
| 1127 | temp->ImportFromResidual(sol, rhs);
|
|---|
| 1128 |
|
|---|
| 1129 | PrintInnerGrid(*temp, information);
|
|---|
| 1130 | }
|
|---|
| 1131 |
|
|---|
| 1132 | void CommMPI::PrintGrid(Grid& grid, const char* information)
|
|---|
| 1133 | {
|
|---|
| 1134 | Index i;
|
|---|
| 1135 | std::stringstream buf;
|
|---|
| 1136 |
|
|---|
| 1137 | Index begin = 0;
|
|---|
| 1138 | Index end = grid.Local().End();
|
|---|
| 1139 | Index begin_local, end_local, begin_global, end_global;
|
|---|
| 1140 |
|
|---|
| 1141 | CommToGhosts(grid);
|
|---|
| 1142 |
|
|---|
| 1143 | for (int j=0; j<3; ++j)
|
|---|
| 1144 | if (grid.Global().EndLocal()[j] == grid.Global().SizeGlobal()[j])
|
|---|
| 1145 | end[j] = grid.Local().SizeTotal()[j];
|
|---|
| 1146 |
|
|---|
| 1147 | for (i.Z()=begin.Z(); i.Z()<end.Z(); ++i.Z())
|
|---|
| 1148 | for (i.Y()=begin.Y(); i.Y()<end.Y(); ++i.Y())
|
|---|
| 1149 | for (i.X()=begin.X(); i.X()<end.X(); ++i.X())
|
|---|
| 1150 | buf << std::scientific << grid.GetVal(i) << " ";
|
|---|
| 1151 |
|
|---|
| 1152 | begin_local = grid.Global().BeginLocal();
|
|---|
| 1153 | end_local = grid.Global().EndLocal();
|
|---|
| 1154 | begin_global = 0;
|
|---|
| 1155 |
|
|---|
| 1156 | for (int i=0; i<3; ++i) {
|
|---|
| 1157 |
|
|---|
| 1158 | if (BoundaryConditions()[i] == Dirichlet ||
|
|---|
| 1159 | BoundaryConditions()[i] == Open) {
|
|---|
| 1160 |
|
|---|
| 1161 | if (begin_local[i] != 0)
|
|---|
| 1162 | begin_local[i] -= grid.Local().HaloEnd1()[i] - grid.Local().HaloBegin1()[i];
|
|---|
| 1163 |
|
|---|
| 1164 | end_global[i] = grid.Global().SizeGlobal()[i];
|
|---|
| 1165 |
|
|---|
| 1166 | }else if (BoundaryConditions()[i] == Periodic) {
|
|---|
| 1167 |
|
|---|
| 1168 | if (end_local[i] == grid.Global().SizeGlobal()[i])
|
|---|
| 1169 | end_local[i] += grid.Local().HaloEnd1()[i] - grid.Local().HaloBegin1()[i] + grid.Local().HaloEnd2()[i] - grid.Local().HaloBegin2()[i];
|
|---|
| 1170 | else
|
|---|
| 1171 | end_local[i] += grid.Local().HaloEnd1()[i] - grid.Local().HaloBegin1()[i];
|
|---|
| 1172 |
|
|---|
| 1173 | end_global[i] = grid.Global().SizeGlobal()[i] +
|
|---|
| 1174 | grid.Local().HaloEnd1()[i] - grid.Local().HaloBegin1()[i] +
|
|---|
| 1175 | grid.Local().HaloEnd2()[i] - grid.Local().HaloBegin2()[i];
|
|---|
| 1176 |
|
|---|
| 1177 | }
|
|---|
| 1178 |
|
|---|
| 1179 | }
|
|---|
| 1180 |
|
|---|
| 1181 | CreateOutputFiles(grid, buf, information,
|
|---|
| 1182 | begin_global, end_global,
|
|---|
| 1183 | begin_local, end_local);
|
|---|
| 1184 | }
|
|---|
| 1185 |
|
|---|
| 1186 | void CommMPI::PrintCorrectedGrid(Grid& grid, const char* information)
|
|---|
| 1187 | {
|
|---|
| 1188 | Index i;
|
|---|
| 1189 |
|
|---|
| 1190 | std::stringstream buf;
|
|---|
| 1191 |
|
|---|
| 1192 | Index begin = 0;
|
|---|
| 1193 | Index end = grid.Local().End();
|
|---|
| 1194 | Index begin_local, end_local, begin_global, end_global;
|
|---|
| 1195 |
|
|---|
| 1196 | CommToGhosts(grid);
|
|---|
| 1197 |
|
|---|
| 1198 | for (int j=0; j<3; ++j)
|
|---|
| 1199 | if (grid.Global().EndLocal()[j] == grid.Global().SizeGlobal()[j])
|
|---|
| 1200 | end[j] = grid.Local().SizeTotal()[j];
|
|---|
| 1201 |
|
|---|
| 1202 | for (i.Z()=begin.Z(); i.Z()<end.Z(); ++i.Z())
|
|---|
| 1203 | for (i.Y()=begin.Y(); i.Y()<end.Y(); ++i.Y())
|
|---|
| 1204 | for (i.X()=begin.X(); i.X()<end.X(); ++i.X())
|
|---|
| 1205 | buf << std::scientific << grid.GetVal(i) - Grid::Correction() << " ";
|
|---|
| 1206 |
|
|---|
| 1207 | begin_local = grid.Global().BeginLocal();
|
|---|
| 1208 | end_local = grid.Global().EndLocal();
|
|---|
| 1209 | begin_global = 0;
|
|---|
| 1210 |
|
|---|
| 1211 | for (int i=0; i<3; ++i) {
|
|---|
| 1212 |
|
|---|
| 1213 | if (BoundaryConditions()[i] == Dirichlet ||
|
|---|
| 1214 | BoundaryConditions()[i] == Open) {
|
|---|
| 1215 |
|
|---|
| 1216 | if (begin_local[i] != 0)
|
|---|
| 1217 | --begin_local[i];
|
|---|
| 1218 |
|
|---|
| 1219 | end_global[i] = grid.Global().SizeGlobal()[i];
|
|---|
| 1220 |
|
|---|
| 1221 | }else if (BoundaryConditions()[i] == Periodic) {
|
|---|
| 1222 |
|
|---|
| 1223 | if (begin_local[i] != 0)
|
|---|
| 1224 | --begin_local[i];
|
|---|
| 1225 |
|
|---|
| 1226 | if (end_local[i] == grid.Global().SizeGlobal()[i])
|
|---|
| 1227 | end_local[i] += 2;
|
|---|
| 1228 |
|
|---|
| 1229 | end_global[i] = grid.Global().SizeGlobal()[i] + 2;
|
|---|
| 1230 |
|
|---|
| 1231 | }
|
|---|
| 1232 |
|
|---|
| 1233 | }
|
|---|
| 1234 |
|
|---|
| 1235 | CreateOutputFiles(grid, buf, information,
|
|---|
| 1236 | 0, grid.Global().SizeGlobal()+2,
|
|---|
| 1237 | begin, end);
|
|---|
| 1238 | }
|
|---|
| 1239 |
|
|---|
| 1240 | void CommMPI::PrintInnerGrid(Grid& grid, const char* information)
|
|---|
| 1241 | {
|
|---|
| 1242 | Index i;
|
|---|
| 1243 | Index begin_local, end_local, begin_global, end_global;
|
|---|
| 1244 | std::stringstream buf;
|
|---|
| 1245 |
|
|---|
| 1246 | CommToGhosts(grid);
|
|---|
| 1247 |
|
|---|
| 1248 | Index begin = 0;
|
|---|
| 1249 | Index end = grid.Local().End();
|
|---|
| 1250 |
|
|---|
| 1251 | /* VTK parallel ImageData files seem to require some overlap. */
|
|---|
| 1252 | for (int j=0; j<3; ++j)
|
|---|
| 1253 | if (grid.Global().BeginLocal()[j] == 0)
|
|---|
| 1254 | ++begin[j];
|
|---|
| 1255 |
|
|---|
| 1256 | /* Write grid values into buffer */
|
|---|
| 1257 | for (i.Z()=begin.Z(); i.Z()<end.Z(); ++i.Z())
|
|---|
| 1258 | for (i.Y()=begin.Y(); i.Y()<end.Y(); ++i.Y())
|
|---|
| 1259 | for (i.X()=begin.X(); i.X()<end.X(); ++i.X())
|
|---|
| 1260 | buf << std::scientific << grid.GetVal(i) << " ";
|
|---|
| 1261 |
|
|---|
| 1262 | for (int j=0; j<3; ++j) {
|
|---|
| 1263 |
|
|---|
| 1264 | if (BoundaryConditions()[j] == Dirichlet ||
|
|---|
| 1265 | BoundaryConditions()[j] == Open) {
|
|---|
| 1266 |
|
|---|
| 1267 | begin_global[j] = 1;
|
|---|
| 1268 | end_global[j] = grid.Global().SizeGlobal()[j] - 1;
|
|---|
| 1269 |
|
|---|
| 1270 | if (grid.Global().BeginLocal()[j] == 0)
|
|---|
| 1271 | begin_local[j] = 1;
|
|---|
| 1272 | else
|
|---|
| 1273 | begin_local[j] = grid.Global().BeginLocal()[j] - 1;
|
|---|
| 1274 |
|
|---|
| 1275 | if (grid.Global().EndLocal()[j] == grid.Global().SizeGlobal()[j])
|
|---|
| 1276 | end_local[j] = grid.Global().EndLocal()[j] - 1;
|
|---|
| 1277 | else
|
|---|
| 1278 | end_local[j] = grid.Global().EndLocal()[j];
|
|---|
| 1279 |
|
|---|
| 1280 | }else if (BoundaryConditions()[j] == Periodic) {
|
|---|
| 1281 |
|
|---|
| 1282 | begin_global[j] = 0;
|
|---|
| 1283 | end_global[j] = grid.Global().SizeGlobal()[j];
|
|---|
| 1284 |
|
|---|
| 1285 | if (grid.Global().BeginLocal()[j] == 0)
|
|---|
| 1286 | begin_local[j] = grid.Global().BeginLocal()[j];
|
|---|
| 1287 | else
|
|---|
| 1288 | begin_local[j] = grid.Global().BeginLocal()[j] - 1;
|
|---|
| 1289 |
|
|---|
| 1290 | end_local[j] = grid.Global().EndLocal()[j];
|
|---|
| 1291 |
|
|---|
| 1292 | }else {
|
|---|
| 1293 |
|
|---|
| 1294 | assert(0 == "Boundary condition not handled properly.");
|
|---|
| 1295 |
|
|---|
| 1296 | }
|
|---|
| 1297 |
|
|---|
| 1298 | }
|
|---|
| 1299 |
|
|---|
| 1300 | CreateOutputFiles(grid, buf, information,
|
|---|
| 1301 | begin_global, end_global,
|
|---|
| 1302 | begin_local, end_local);
|
|---|
| 1303 | }
|
|---|
| 1304 |
|
|---|
| 1305 | void CommMPI::PrintInnerCorrectedGrid(Grid& grid, const char* information)
|
|---|
| 1306 | {
|
|---|
| 1307 | Index i;
|
|---|
| 1308 | const Index begin = grid.Local().Begin();
|
|---|
| 1309 | const Index end = grid.Local().End();
|
|---|
| 1310 | std::stringstream buf;
|
|---|
| 1311 |
|
|---|
| 1312 | for (i.Z()=begin.Z(); i.Z()<end.Z(); ++i.Z())
|
|---|
| 1313 | for (i.Y()=begin.Y(); i.Y()<end.Y(); ++i.Y())
|
|---|
| 1314 | for (i.X()=begin.X(); i.X()<end.X(); ++i.X())
|
|---|
| 1315 | buf << std::scientific << grid.GetVal(i) - Grid::Correction() << " ";
|
|---|
| 1316 |
|
|---|
| 1317 | CreateOutputFiles(grid, buf, information,
|
|---|
| 1318 | 0, grid.Global().SizeGlobal(),
|
|---|
| 1319 | grid.Global().BeginLocal(), grid.Global().EndLocal());
|
|---|
| 1320 | }
|
|---|
| 1321 |
|
|---|
| 1322 | void CommMPI::CreateOutputFiles(const Grid& grid, const std::stringstream& serial_data, const char* information,
|
|---|
| 1323 | const Index& begin_global, const Index& end_global,
|
|---|
| 1324 | const Index& begin_local, const Index& end_local)
|
|---|
| 1325 | {
|
|---|
| 1326 | int output_count = OutputCount();
|
|---|
| 1327 |
|
|---|
| 1328 | MPI_Comm comm = comm_info.GetCommunicator(grid);
|
|---|
| 1329 |
|
|---|
| 1330 | if (comm != MPI_COMM_NULL) {
|
|---|
| 1331 |
|
|---|
| 1332 | MPI_File file;
|
|---|
| 1333 | std::string conv_information = Helper::ReplaceWhitespaces(information, "_");
|
|---|
| 1334 |
|
|---|
| 1335 | CreateParallelOutputFile(grid, comm, output_count, conv_information.c_str(),
|
|---|
| 1336 | begin_global, end_global, begin_local, end_local);
|
|---|
| 1337 |
|
|---|
| 1338 | file = CreateSerialOutputFile(grid, comm, output_count, conv_information.c_str(),
|
|---|
| 1339 | begin_global, end_global, begin_local, end_local);
|
|---|
| 1340 |
|
|---|
| 1341 | char *char_buf = Helper::GetCharArray(serial_data.str());
|
|---|
| 1342 | MPI_File_write(file, char_buf, serial_data.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
|
|---|
| 1343 | delete [] char_buf;
|
|---|
| 1344 |
|
|---|
| 1345 | FinalizeSerialOutputFile(file);
|
|---|
| 1346 |
|
|---|
| 1347 | }
|
|---|
| 1348 | }
|
|---|
| 1349 |
|
|---|
| 1350 | void CommMPI::CreateParallelOutputFile(const Grid& grid, MPI_Comm& comm,
|
|---|
| 1351 | const int& output_count, const char* information,
|
|---|
| 1352 | const Index& begin_global, const Index& end_global,
|
|---|
| 1353 | const Index& begin_local, const Index& end_local)
|
|---|
| 1354 | {
|
|---|
| 1355 | int rank;
|
|---|
| 1356 | MPI_File file;
|
|---|
| 1357 | char parallel_filename[513], serial_filename[513];
|
|---|
| 1358 | std::stringstream buf;
|
|---|
| 1359 |
|
|---|
| 1360 | MPI_Comm_rank(comm, &rank);
|
|---|
| 1361 |
|
|---|
| 1362 | sprintf(parallel_filename, "%s%04d.pvti", OutputPath().c_str(), output_count);
|
|---|
| 1363 | sprintf(serial_filename, "%04d_%d.vti", output_count, rank);
|
|---|
| 1364 |
|
|---|
| 1365 | MPI_File_open(comm, parallel_filename, MPI_MODE_WRONLY|MPI_MODE_CREATE, MPI_INFO_NULL, &file);
|
|---|
| 1366 | MPI_File_set_size(file, 0);
|
|---|
| 1367 |
|
|---|
| 1368 | if (rank == 0) {
|
|---|
| 1369 |
|
|---|
| 1370 | buf << "<?xml version=\"1.0\"?>" << std::endl
|
|---|
| 1371 | << "<VTKFile type=\"PImageData\" version=\"0.1\" byte_order=\"LittleEndian\">" << std::endl
|
|---|
| 1372 | << " <PImageData WholeExtent=\"";
|
|---|
| 1373 |
|
|---|
| 1374 | for (int i=0; i<3; ++i)
|
|---|
| 1375 | buf << begin_global[i] << " " << end_global[i] - 1 << " ";
|
|---|
| 1376 |
|
|---|
| 1377 | buf << "\"" << std::endl
|
|---|
| 1378 | << " GhostLevel=\"0\" Origin=\"0 0 0\" Spacing=\"";
|
|---|
| 1379 |
|
|---|
| 1380 | for (int i=0; i<3; ++i)
|
|---|
| 1381 | buf << grid.Extent().MeshWidth()[i] << " ";
|
|---|
| 1382 |
|
|---|
| 1383 | buf << "\">" << std::endl
|
|---|
| 1384 | << " <PPointData Scalars=\"" << information << "\">" << std::endl
|
|---|
| 1385 | << " <PDataArray type=\"Float32\" Name=\"" << information << "\"/>" << std::endl
|
|---|
| 1386 | << " </PPointData>" << std::endl;
|
|---|
| 1387 |
|
|---|
| 1388 | char* char_buf = Helper::GetCharArray(buf.str());
|
|---|
| 1389 |
|
|---|
| 1390 | MPI_File_write_shared(file, char_buf, buf.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
|
|---|
| 1391 |
|
|---|
| 1392 | delete [] char_buf;
|
|---|
| 1393 | }
|
|---|
| 1394 |
|
|---|
| 1395 | buf.str("");
|
|---|
| 1396 | buf << " <Piece Extent=\"";
|
|---|
| 1397 |
|
|---|
| 1398 | for (int i=0; i<3; ++i)
|
|---|
| 1399 | buf << begin_local[i] << " " << end_local[i] - 1 << " ";
|
|---|
| 1400 |
|
|---|
| 1401 | buf << "\" Source=\"" << serial_filename << "\"/>" << std::endl;
|
|---|
| 1402 |
|
|---|
| 1403 | char* char_buf = Helper::GetCharArray(buf.str());
|
|---|
| 1404 |
|
|---|
| 1405 | MPI_File_write_ordered(file, char_buf, buf.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
|
|---|
| 1406 |
|
|---|
| 1407 | delete [] char_buf;
|
|---|
| 1408 |
|
|---|
| 1409 | if (rank == 0) {
|
|---|
| 1410 |
|
|---|
| 1411 | buf.str("");
|
|---|
| 1412 |
|
|---|
| 1413 | buf << " </PImageData>" << std::endl
|
|---|
| 1414 | << "</VTKFile>" << std::endl;
|
|---|
| 1415 |
|
|---|
| 1416 | char* char_buf = Helper::GetCharArray(buf.str());
|
|---|
| 1417 |
|
|---|
| 1418 | MPI_File_write_shared(file, char_buf, buf.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
|
|---|
| 1419 |
|
|---|
| 1420 | delete [] char_buf;
|
|---|
| 1421 | }
|
|---|
| 1422 |
|
|---|
| 1423 | MPI_File_close(&file);
|
|---|
| 1424 | }
|
|---|
| 1425 |
|
|---|
| 1426 | MPI_File CommMPI::CreateSerialOutputFile(const Grid& grid, MPI_Comm& comm,
|
|---|
| 1427 | const int& output_count, const char* information,
|
|---|
| 1428 | const Index& begin_global, const Index& end_global,
|
|---|
| 1429 | const Index& begin_local, const Index& end_local)
|
|---|
| 1430 | {
|
|---|
| 1431 | char serial_filename[513];
|
|---|
| 1432 | int rank;
|
|---|
| 1433 | MPI_File file;
|
|---|
| 1434 | std::stringstream buf;
|
|---|
| 1435 |
|
|---|
| 1436 | MPI_Comm_rank(comm_global, &rank);
|
|---|
| 1437 |
|
|---|
| 1438 | sprintf(serial_filename, "%s%04d_%d.vti", OutputPath().c_str(), output_count, rank);
|
|---|
| 1439 |
|
|---|
| 1440 | MPI_File_open(MPI_COMM_SELF, serial_filename, MPI_MODE_WRONLY|MPI_MODE_CREATE, MPI_INFO_NULL, &file);
|
|---|
| 1441 |
|
|---|
| 1442 | buf << "<?xml version=\"1.0\"?>" << std::endl
|
|---|
| 1443 | << "<VTKFile type=\"ImageData\" version=\"0.1\" byte_order=\"LittleEndian\">" << std::endl
|
|---|
| 1444 | << " <ImageData WholeExtent=\"";
|
|---|
| 1445 |
|
|---|
| 1446 | for (int i=0; i<3; ++i)
|
|---|
| 1447 | buf << begin_global[i] << " " << end_global[i] - 1 << " ";
|
|---|
| 1448 |
|
|---|
| 1449 | buf << "\"" << std::endl
|
|---|
| 1450 | << " Origin=\"0 0 0\" Spacing=\"";
|
|---|
| 1451 |
|
|---|
| 1452 | for (int i=0; i<3; ++i)
|
|---|
| 1453 | buf << grid.Extent().MeshWidth()[i] << " ";
|
|---|
| 1454 |
|
|---|
| 1455 | buf << "\">" << std::endl
|
|---|
| 1456 | << " <Piece Extent=\"";
|
|---|
| 1457 |
|
|---|
| 1458 | for (int i=0; i<3; ++i)
|
|---|
| 1459 | buf << begin_local[i] << " " << end_local[i] - 1 << " ";
|
|---|
| 1460 |
|
|---|
| 1461 | buf << "\">" << std::endl
|
|---|
| 1462 | << " <PointData Scalars=\"" << information << "\">" << std::endl
|
|---|
| 1463 | << " <DataArray type=\"Float32\" Name=\"" << information << "\" format=\"ascii\">" << std::endl
|
|---|
| 1464 | << " ";
|
|---|
| 1465 |
|
|---|
| 1466 | char* char_buf = Helper::GetCharArray(buf.str());
|
|---|
| 1467 | MPI_File_write(file, char_buf, buf.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
|
|---|
| 1468 | delete [] char_buf;
|
|---|
| 1469 |
|
|---|
| 1470 | return file;
|
|---|
| 1471 | }
|
|---|
| 1472 |
|
|---|
| 1473 | void CommMPI::FinalizeSerialOutputFile(MPI_File& file)
|
|---|
| 1474 | {
|
|---|
| 1475 | std::stringstream buf;
|
|---|
| 1476 |
|
|---|
| 1477 | buf << std::endl
|
|---|
| 1478 | << " </DataArray>" << std::endl
|
|---|
| 1479 | << " </PointData>" << std::endl
|
|---|
| 1480 | << " </Piece>" << std::endl
|
|---|
| 1481 | << " </ImageData>" << std::endl
|
|---|
| 1482 | << "</VTKFile>" << std::endl;
|
|---|
| 1483 |
|
|---|
| 1484 | char* char_buf = Helper::GetCharArray(buf.str());
|
|---|
| 1485 | MPI_File_write(file, char_buf, buf.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
|
|---|
| 1486 | delete [] char_buf;
|
|---|
| 1487 |
|
|---|
| 1488 | MPI_File_close(&file);
|
|---|
| 1489 | }
|
|---|
| 1490 |
|
|---|
| 1491 | int CommMPI::GlobalRank() const
|
|---|
| 1492 | {
|
|---|
| 1493 | int rank;
|
|---|
| 1494 |
|
|---|
| 1495 | MPI_Comm_rank(comm_global, &rank);
|
|---|
| 1496 |
|
|---|
| 1497 | return rank;
|
|---|
| 1498 | }
|
|---|
| 1499 |
|
|---|
| 1500 | Index CommMPI::GlobalPos() const
|
|---|
| 1501 | {
|
|---|
| 1502 | Index dims, periods, coords;
|
|---|
| 1503 |
|
|---|
| 1504 | MPI_Cart_get(comm_global, 3, dims.vec(), periods.vec(), coords.vec());
|
|---|
| 1505 |
|
|---|
| 1506 | return coords;
|
|---|
| 1507 | }
|
|---|
| 1508 |
|
|---|
| 1509 | Index CommMPI::GlobalProcs() const
|
|---|
| 1510 | {
|
|---|
| 1511 | Index dims, periods, coords;
|
|---|
| 1512 |
|
|---|
| 1513 | MPI_Cart_get(comm_global, 3, dims.vec(), periods.vec(), coords.vec());
|
|---|
| 1514 |
|
|---|
| 1515 | return dims;
|
|---|
| 1516 | }
|
|---|
| 1517 |
|
|---|
| 1518 | void CommMPI::InitCommMPI(const MPI_Comm& comm)
|
|---|
| 1519 | {
|
|---|
| 1520 | int status, size, rank;
|
|---|
| 1521 | int dims[3] = {0, 0, 0};
|
|---|
| 1522 | int periods[3];
|
|---|
| 1523 |
|
|---|
| 1524 | global_coarse_grid = NULL;
|
|---|
| 1525 |
|
|---|
| 1526 | for (int i=0; i<3; ++i)
|
|---|
| 1527 | periods[i] = (BoundaryConditions()[i] == Periodic ? 1 : 0);
|
|---|
| 1528 |
|
|---|
| 1529 | MPI_Comm_size(comm, &size);
|
|---|
| 1530 | MPI_Comm_rank(comm, &rank);
|
|---|
| 1531 |
|
|---|
| 1532 | MPI_Topo_test(comm, &status);
|
|---|
| 1533 |
|
|---|
| 1534 | if (status == MPI_CART) {
|
|---|
| 1535 |
|
|---|
| 1536 | comm_global = comm;
|
|---|
| 1537 |
|
|---|
| 1538 | }else {
|
|---|
| 1539 |
|
|---|
| 1540 | MPI_Dims_create(size, 3, dims);
|
|---|
| 1541 |
|
|---|
| 1542 | MPI_Cart_create(comm, 3, dims, periods, 1, &comm_global);
|
|---|
| 1543 |
|
|---|
| 1544 | }
|
|---|
| 1545 |
|
|---|
| 1546 | comm_info.SetCommGlobal(comm_global);
|
|---|
| 1547 |
|
|---|
| 1548 | MPI_Info_create(&info);
|
|---|
| 1549 | char key[] = "no_locks";
|
|---|
| 1550 | char val[] = "true";
|
|---|
| 1551 | MPI_Info_set(info, key, val);
|
|---|
| 1552 |
|
|---|
| 1553 | }
|
|---|
| 1554 |
|
|---|
| 1555 | CommMPI::~CommMPI()
|
|---|
| 1556 | {
|
|---|
| 1557 | std::map<const Grid*, TempGrid*>::iterator iter;
|
|---|
| 1558 |
|
|---|
| 1559 | for (iter=finer_grids.begin(); iter!=finer_grids.end(); ++iter)
|
|---|
| 1560 | delete iter->second;
|
|---|
| 1561 |
|
|---|
| 1562 | for (iter=coarser_grids.begin(); iter!=coarser_grids.end(); ++iter)
|
|---|
| 1563 | delete iter->second;
|
|---|
| 1564 |
|
|---|
| 1565 | MPI_Comm_free(&comm_global);
|
|---|
| 1566 |
|
|---|
| 1567 | if (win_created)
|
|---|
| 1568 | MPI_Win_free(&win);
|
|---|
| 1569 |
|
|---|
| 1570 | MPI_Info_free(&info);
|
|---|
| 1571 |
|
|---|
| 1572 | delete global_coarse_grid;
|
|---|
| 1573 | }
|
|---|
| 1574 |
|
|---|
| 1575 | std::string CommMPI::CreateOutputDirectory()
|
|---|
| 1576 | {
|
|---|
| 1577 | #ifdef HAVE_BOOST_FILESYSTEM
|
|---|
| 1578 | std::string path, unique_path;
|
|---|
| 1579 | std::stringstream unique_suffix;
|
|---|
| 1580 | int suffix_counter = 0;
|
|---|
| 1581 | char buffer[129];
|
|---|
| 1582 | time_t rawtime;
|
|---|
| 1583 | struct tm *timeinfo;
|
|---|
| 1584 | int path_size;
|
|---|
| 1585 | char* path_array;
|
|---|
| 1586 |
|
|---|
| 1587 | if (GlobalRank() == 0) {
|
|---|
| 1588 |
|
|---|
| 1589 | time(&rawtime);
|
|---|
| 1590 | timeinfo = localtime(&rawtime);
|
|---|
| 1591 | strftime(buffer, 128, "./output/%Y_%m_%d_%H_%M_%S/", timeinfo);
|
|---|
| 1592 | path = buffer;
|
|---|
| 1593 |
|
|---|
| 1594 | if (!fs::exists("output"))
|
|---|
| 1595 | fs::create_directory("output");
|
|---|
| 1596 |
|
|---|
| 1597 | unique_path = path;
|
|---|
| 1598 |
|
|---|
| 1599 | while (fs::exists(unique_path.c_str())) {
|
|---|
| 1600 |
|
|---|
| 1601 | unique_suffix.str("");
|
|---|
| 1602 | unique_suffix << "_" << suffix_counter++ << "/";
|
|---|
| 1603 |
|
|---|
| 1604 | unique_path = path;
|
|---|
| 1605 | unique_path.replace(unique_path.size()-1, 1, unique_suffix.str());
|
|---|
| 1606 |
|
|---|
| 1607 | }
|
|---|
| 1608 |
|
|---|
| 1609 | fs::create_directory(unique_path.c_str());
|
|---|
| 1610 |
|
|---|
| 1611 | path_size = unique_path.size() + 1;
|
|---|
| 1612 | path_array = Helper::GetCharArray(unique_path);
|
|---|
| 1613 |
|
|---|
| 1614 | MPI_Bcast(&path_size, 1, MPI_INT, 0, comm_global);
|
|---|
| 1615 |
|
|---|
| 1616 | }else {
|
|---|
| 1617 |
|
|---|
| 1618 | MPI_Bcast(&path_size, 1, MPI_INT, 0, comm_global);
|
|---|
| 1619 | path_array = new char[path_size];
|
|---|
| 1620 |
|
|---|
| 1621 | }
|
|---|
| 1622 |
|
|---|
| 1623 | MPI_Bcast(path_array, path_size, MPI_CHAR, 0, comm_global);
|
|---|
| 1624 |
|
|---|
| 1625 | unique_path = path_array;
|
|---|
| 1626 |
|
|---|
| 1627 | delete [] path_array;
|
|---|
| 1628 |
|
|---|
| 1629 | return unique_path;
|
|---|
| 1630 |
|
|---|
| 1631 | #else
|
|---|
| 1632 |
|
|---|
| 1633 | return "./";
|
|---|
| 1634 |
|
|---|
| 1635 | #endif
|
|---|
| 1636 | }
|
|---|
| 1637 |
|
|---|
| 1638 |
|
|---|
| 1639 | #endif /* HAVE_MPI */
|
|---|