source: src/comm/comm_mpi.cpp@ f57182

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

Open boundary conditions.

Conflicts:

lib/vmg/src/Makefile.am
lib/vmg/src/base/factory.cpp
lib/vmg/test/unit_test/library/dirichlet_fas_lr_mpi.cpp

  • Property mode set to 100644
File size: 29.7 KB
Line 
1/*
2 * vmg - a versatile multigrid solver
3 * Copyright (C) 2012 Institute for Numerical Simulation, University of Bonn
4 *
5 * vmg is free software: you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation, either version 3 of the License, or
8 * (at your option) any later version.
9 *
10 * vmg is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program. If not, see <http://www.gnu.org/licenses/>.
17 */
18
19/**
20 * @file comm_mpi.cpp
21 * @author Julian Iseringhausen <isering@ins.uni-bonn.de>
22 * @date Wed Jun 16 13:21:06 2010
23 *
24 * @brief Class for MPI-based communication.
25 *
26 */
27
28#ifdef HAVE_CONFIG_H
29#include <config.h>
30#endif
31
32#ifdef HAVE_MPI
33
34#include <mpi.h>
35#ifdef HAVE_MARMOT
36#include <enhancempicalls.h>
37#include <sourceinfompicalls.h>
38#endif
39
40#ifdef HAVE_BOOST_FILESYSTEM
41#include <boost/filesystem.hpp>
42namespace fs = boost::filesystem;
43#endif
44
45#ifdef HAVE_VTK
46#include <vtkAbstractArray.h>
47#include <vtkImageData.h>
48#include <vtkPointData.h>
49#include <vtkSmartPointer.h>
50#include <vtkXMLImageDataWriter.h>
51#endif
52
53#include <cstring>
54#include <sstream>
55
56#include "base/helper.hpp"
57#include "base/tuple.hpp"
58#include "comm/comm_mpi.hpp"
59#include "comm/mpi/datatypes_local.hpp"
60#include "grid/grid.hpp"
61#include "grid/multigrid.hpp"
62#include "grid/tempgrid.hpp"
63#include "mg.hpp"
64#include "base/timer.hpp"
65
66static char print_buffer[512];
67
68using namespace VMG;
69
70void CommMPI::IsendAll(Grid& grid, std::vector<VMG::MPI::Datatype>& types, const MPI_Comm& comm, const int& tag_start)
71{
72 for (std::vector<VMG::MPI::Datatype>::const_iterator iter=types.begin(); iter!=types.end(); ++iter)
73 iter->Isend(grid, tag_start, comm, Request());
74}
75
76void CommMPI::IrecvAll(Grid& grid, std::vector<VMG::MPI::Datatype>& types, const MPI_Comm& comm, const int& tag_start)
77{
78 for (std::vector<VMG::MPI::Datatype>::const_iterator iter=types.begin(); iter!=types.end(); ++iter)
79 iter->Irecv(grid, tag_start, comm, Request());
80}
81
82void CommMPI::IsendAllBuffered(const Grid& grid, std::vector<VMG::MPI::Datatype>& types, const MPI_Comm& comm, const int& tag_start)
83{
84 for (std::vector<VMG::MPI::Datatype>::iterator iter=types.begin(); iter!=types.end(); ++iter)
85 iter->IsendBuffered(grid, tag_start, comm, Request());
86}
87
88void CommMPI::IrecvAllBuffered(std::vector<VMG::MPI::Datatype>& types, const MPI_Comm& comm, const int& tag_start)
89{
90 for (std::vector<VMG::MPI::Datatype>::iterator iter=types.begin(); iter!=types.end(); ++iter)
91 iter->IrecvBuffered(tag_start, comm, Request());
92}
93
94void CommMPI::ReplaceBufferAll(Grid& grid, const std::vector<VMG::MPI::Datatype>& types)
95{
96 for (std::vector<VMG::MPI::Datatype>::const_iterator iter=types.begin(); iter!= types.end(); ++iter)
97 iter->GridReplace(grid);
98}
99
100void CommMPI::AddBufferAll(Grid& grid, const std::vector<VMG::MPI::Datatype>& types)
101{
102 for (std::vector<VMG::MPI::Datatype>::const_iterator iter=types.begin(); iter!= types.end(); ++iter)
103 iter->GridSum(grid);
104}
105
106void CommMPI::CommSubgrid(Grid& grid_old, Grid& grid_new, const int& direction)
107{
108 MPI_Comm comm = settings.CommunicatorGlobal(grid_old);
109 if (comm != MPI_COMM_NULL) {
110 VMG::MPI::DatatypesGlobal& datatypes = settings.DatatypesGlobal(grid_old, grid_new, direction);
111 IrecvAllBuffered(datatypes.Receive(), comm, 0411);
112 IsendAllBuffered(grid_old, datatypes.Send(), comm, 0411);
113 WaitAll();
114 ReplaceBufferAll(grid_new, datatypes.Receive());
115 }
116}
117
118void CommMPI::CommAddSubgrid(Grid& grid_old, Grid& grid_new, const int& direction)
119{
120 MPI_Comm comm = settings.CommunicatorGlobal(grid_old);
121 if (comm != MPI_COMM_NULL) {
122 VMG::MPI::DatatypesGlobal& datatypes = settings.DatatypesGlobal(grid_old, grid_new, direction);
123 IrecvAllBuffered(datatypes.Receive(), comm, 1806);
124 IsendAllBuffered(grid_old, datatypes.Send(), comm, 1806);
125 WaitAll();
126 AddBufferAll(grid_new, datatypes.Receive());
127 }
128}
129
130void CommMPI::CommToGhosts(Grid& grid)
131{
132 MPI_Comm comm = settings.CommunicatorLocal(grid);
133 if (comm != MPI_COMM_NULL) {
134 VMG::MPI::DatatypesLocal& types = settings.DatatypesLocal(grid);
135 IrecvAllBuffered(types.Halo(), comm, 2310);
136 IsendAllBuffered(grid, types.NB(), comm, 2310);
137 WaitAll();
138 ReplaceBufferAll(grid, types.Halo());
139 }
140}
141
142void CommMPI::CommToGhostsAsyncStart(Grid& grid)
143{
144 MPI_Comm comm = settings.CommunicatorLocal(grid);
145 if (comm != MPI_COMM_NULL) {
146 VMG::MPI::DatatypesLocal& types = settings.DatatypesLocal(grid);
147 IrecvAllBuffered(types.Halo(), comm, 2412);
148 IsendAllBuffered(grid, types.NB(), comm, 2412);
149 TestAll();
150 }
151}
152
153void CommMPI::CommToGhostsAsyncFinish(Grid& grid)
154{
155 WaitAll();
156 ReplaceBufferAll(grid, settings.DatatypesLocal(grid).Halo());
157}
158
159void CommMPI::CommFromGhosts(Grid& grid)
160{
161 MPI_Comm comm = settings.CommunicatorLocal(grid);
162 if (comm != MPI_COMM_NULL) {
163 VMG::MPI::DatatypesLocal& types = settings.DatatypesLocal(grid);
164 IrecvAllBuffered(types.NB(), comm, 1337);
165 IsendAllBuffered(grid, types.Halo(), comm, 1337);
166 WaitAll();
167 AddBufferAll(grid, types.NB());
168 }
169}
170
171void CommMPI::CommFromGhostsAsyncStart(Grid& grid)
172{
173 MPI_Comm comm = settings.CommunicatorLocal(grid);
174 if (comm != MPI_COMM_NULL) {
175 VMG::MPI::DatatypesLocal& types = settings.DatatypesLocal(grid);
176 IrecvAllBuffered(types.NB(), comm, 0xc0ffee);
177 IsendAllBuffered(grid, types.Halo(), comm, 0xc0ffee);
178 TestAll();
179 }
180}
181
182void CommMPI::CommFromGhostsAsyncFinish(Grid& grid)
183{
184 WaitAll();
185 AddBufferAll(grid, settings.DatatypesLocal(grid).NB());
186}
187
188void CommMPI::Barrier()
189{
190 MPI_Barrier(comm_global);
191}
192
193vmg_float CommMPI::GlobalSum(vmg_float value)
194{
195 vmg_float result = 0;
196 MPI_Allreduce(&value, &result, 1, MPI_DOUBLE, MPI_SUM, comm_global);
197 return result;
198}
199
200vmg_float CommMPI::GlobalSumRoot(vmg_float value)
201{
202 vmg_float result = 0;
203 MPI_Reduce(&value, &result, 1, MPI_DOUBLE, MPI_SUM, 0, comm_global);
204 return result;
205}
206
207void CommMPI::GlobalSumArray(vmg_float* array, const vmg_int& size)
208{
209 MPI_Allreduce(MPI_IN_PLACE, array, size, MPI_DOUBLE, MPI_SUM, comm_global);
210}
211
212void CommMPI::GlobalBroadcast(vmg_float& value)
213{
214 MPI_Bcast(&value, 1, MPI_DOUBLE, 0, comm_global);
215}
216
217void CommMPI::GlobalGather(vmg_float& value, vmg_float* array)
218{
219 MPI_Gather(&value, 1, MPI_DOUBLE, array, 1, MPI_DOUBLE, 0, comm_global);
220}
221
222vmg_int CommMPI::GlobalSum(vmg_int value)
223{
224 vmg_int result = 0;
225 MPI_Allreduce(&value, &result, 1, MPI_INT, MPI_SUM, comm_global);
226 return result;
227}
228
229vmg_int CommMPI::GlobalSumRoot(vmg_int value)
230{
231 vmg_int result = 0;
232 MPI_Reduce(&value, &result, 1, MPI_INT, MPI_SUM, 0, comm_global);
233 return result;
234}
235
236void CommMPI::GlobalSumArray(vmg_int* array, const vmg_int& size)
237{
238 MPI_Allreduce(MPI_IN_PLACE, array, size, MPI_INT, MPI_SUM, comm_global);
239}
240
241vmg_float CommMPI::GlobalMax(vmg_float value)
242{
243 vmg_float result = 0.0;
244 MPI_Allreduce(&value, &result, 1, MPI_DOUBLE, MPI_MAX, comm_global);
245 return result;
246}
247
248vmg_float CommMPI::GlobalMaxRoot(vmg_float value)
249{
250 vmg_float result = 0.0;
251 MPI_Reduce(&value, &result, 1, MPI_DOUBLE, MPI_MAX, 0, comm_global);
252 return result;
253}
254
255void CommMPI::GlobalMaxArray(vmg_float* array, const vmg_int& size)
256{
257 MPI_Allreduce(MPI_IN_PLACE, array, size, MPI_DOUBLE, MPI_MAX, comm_global);
258}
259
260vmg_int CommMPI::GlobalMax(vmg_int value)
261{
262 vmg_int result = 0;
263 MPI_Allreduce(&value, &result, 1, MPI_INT, MPI_MAX, comm_global);
264 return result;
265}
266
267vmg_int CommMPI::GlobalMaxRoot(vmg_int value)
268{
269 vmg_int result = 0;
270 MPI_Reduce(&value, &result, 1, MPI_INT, MPI_MAX, 0, comm_global);
271 return result;
272}
273
274void CommMPI::GlobalMaxArray(vmg_int* array, const vmg_int& size)
275{
276 MPI_Allreduce(MPI_IN_PLACE, array, size, MPI_INT, MPI_MAX, comm_global);
277}
278
279void CommMPI::GlobalBroadcast(vmg_int& value)
280{
281 MPI_Bcast(&value, 1, MPI_INT, 0, comm_global);
282}
283
284void CommMPI::GlobalGather(vmg_int& value, vmg_int* array)
285{
286 MPI_Gather(&value, 1, MPI_INT, array, 1, MPI_INT, 0, comm_global);
287}
288
289void CommMPI::GlobalBroadcast(char* str)
290{
291 int size = std::strlen(str) + 1;
292 MPI_Bcast(&size, 1, MPI_INT, 0, comm_global);
293 MPI_Bcast(str, size, MPI_CHAR, 0, comm_global);
294}
295
296vmg_float CommMPI::LevelSum(const Grid& grid, vmg_float value)
297{
298 vmg_float result = 0.0;
299 MPI_Comm comm = settings.CommunicatorLocal(grid);
300 assert(comm != MPI_COMM_NULL);
301 MPI_Allreduce(&value, &result, 1, MPI_DOUBLE, MPI_SUM, comm);
302 return result;
303}
304
305vmg_float CommMPI::LevelSumRoot(const Grid& grid, vmg_float value)
306{
307 vmg_float result = 0.0;
308 MPI_Comm comm = settings.CommunicatorLocal(grid);
309 assert(comm != MPI_COMM_NULL);
310 MPI_Reduce(&value, &result, 1, MPI_DOUBLE, MPI_SUM, 0, comm);
311 return result;
312}
313
314void CommMPI::LevelSumArray(const Grid& grid, vmg_float* array, const vmg_int& size)
315{
316 MPI_Comm comm = settings.CommunicatorLocal(grid);
317 assert(comm != MPI_COMM_NULL);
318 MPI_Allreduce(MPI_IN_PLACE, array, size, MPI_DOUBLE, MPI_SUM, comm);
319}
320
321vmg_int CommMPI::LevelSum(const Grid& grid, vmg_int value)
322{
323 vmg_int result = 0;
324 MPI_Comm comm = settings.CommunicatorLocal(grid);
325 assert(comm != MPI_COMM_NULL);
326 MPI_Allreduce(&value, &result, 1, MPI_INT, MPI_SUM, comm);
327 return result;
328}
329
330vmg_int CommMPI::LevelSumRoot(const Grid& grid, vmg_int value)
331{
332 vmg_int result = 0;
333 MPI_Comm comm = settings.CommunicatorLocal(grid);
334 assert(comm != MPI_COMM_NULL);
335 MPI_Reduce(&value, &result, 1, MPI_INT, MPI_SUM, 0, comm);
336 return result;
337}
338
339void CommMPI::LevelSumArray(const Grid& grid, vmg_int* array, const vmg_int& size)
340{
341 MPI_Comm comm = settings.CommunicatorLocal(grid);
342 assert(comm != MPI_COMM_NULL);
343 MPI_Allreduce(MPI_IN_PLACE, array, size, MPI_INT, MPI_SUM, comm);
344}
345
346void CommMPI::Print(const OutputLevel level, const char* format, ...)
347{
348 bool print = (level == Output);
349
350#ifdef OUTPUT_INFO
351 print |= (level == Info);
352#endif
353#ifdef OUTPUT_DEBUG
354 print |= (level == Debug);
355#endif
356#ifdef OUTPUT_TIMING
357 print |= (level == Timing);
358#endif
359
360 if (print) {
361 va_list args;
362 va_start(args, format);
363 vsprintf(print_buffer, format, args);
364 printf("VMG: Rank %d: %s\n", GlobalRank(), print_buffer);
365 va_end(args);
366 }
367}
368
369void CommMPI::PrintOnce(const OutputLevel level, const char* format, ...)
370{
371 bool print = (level == Output);
372
373#ifdef OUTPUT_INFO
374 print |= (level == Info);
375#endif
376#ifdef OUTPUT_DEBUG
377 print |= (level == Debug);
378#endif
379#ifdef OUTPUT_TIMING
380 print |= (level == Timing);
381#endif
382
383 if (GlobalRank() == 0 && print) {
384 va_list args;
385 va_start(args, format);
386 vsprintf(print_buffer, format, args);
387 printf("VMG: Rank %d: %s\n", GlobalRank(), print_buffer);
388 va_end(args);
389 }
390}
391
392void CommMPI::PrintXML(const std::string& filename, const std::string& xml_data)
393{
394 MPI_File file;
395 std::stringstream path, xml_header;
396
397 pugi::xml_document doc;
398 pugi::xml_node node_data = doc.append_child("Global").append_child("NumProcesses").append_child(pugi::node_pcdata);
399 node_data.set_value(Helper::ToString(GlobalProcs().Product()).c_str());
400 doc.save(xml_header);
401
402 path << OutputPath() << filename;
403
404 char* filename_array = Helper::GetCharArray(path.str());
405 char* xml_header_array = Helper::GetCharArray(xml_header.str());
406 char* str_array = Helper::GetCharArray(xml_data);
407
408 MPI_File_open(MPI_COMM_SELF, filename_array, MPI_MODE_WRONLY|MPI_MODE_CREATE, MPI_INFO_NULL, &file);
409 MPI_File_set_size(file, 0);
410 MPI_File_write(file, xml_header_array, xml_header.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
411 MPI_File_write(file, str_array, xml_data.size(), MPI_CHAR, MPI_STATUS_IGNORE);
412 MPI_File_close(&file);
413
414 delete [] filename_array;
415 delete [] xml_header_array;
416 delete [] str_array;
417}
418
419void CommMPI::PrintXMLAll(const std::string& filename, const std::string& xml_data)
420{
421 MPI_File file;
422 std::stringstream path;
423
424 path << OutputPath() << filename;
425
426 char* filename_array = Helper::GetCharArray(path.str());
427 char* str_array = Helper::GetCharArray(xml_data);
428
429 MPI_File_open(comm_global, filename_array, MPI_MODE_WRONLY|MPI_MODE_CREATE, MPI_INFO_NULL, &file);
430 MPI_File_set_size(file, 0);
431
432 if (GlobalRank() == 0) {
433 std::stringstream xml_header;
434 pugi::xml_document doc;
435 pugi::xml_node node_data = doc.append_child("Global").append_child("NumProcesses").append_child(pugi::node_pcdata);
436 node_data.set_value(Helper::ToString(GlobalProcs().Product()).c_str());
437 doc.save(xml_header);
438
439 char* xml_header_array = Helper::GetCharArray(xml_header.str());
440
441 MPI_File_write_shared(file, xml_header_array, xml_header.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
442
443 delete [] xml_header_array;
444 }
445
446 MPI_File_write_ordered(file, str_array, xml_data.size(), MPI_CHAR, MPI_STATUS_IGNORE);
447 MPI_File_close(&file);
448
449 delete [] filename_array;
450 delete [] str_array;
451}
452
453void CommMPI::PrintGridInformation(const Grid& grid, char* filename, const std::string& name)
454{
455 std::stringstream buf;
456 MPI_File file;
457 int rank, size;
458 int size_local, size_local_max;
459
460 MPI_Comm comm = settings.CommunicatorGlobal(grid);
461 MPI_Comm comm_local = settings.CommunicatorLocal(grid);
462
463 if (comm_local != MPI_COMM_NULL)
464 MPI_Comm_size(comm_local, &size_local);
465 else
466 size_local = 0;
467
468 if (comm != MPI_COMM_NULL) {
469
470 MPI_Reduce(&size_local, &size_local_max, 1, MPI_INT, MPI_MAX, 0, comm);
471
472 MPI_File_open(comm, filename, MPI_MODE_WRONLY|MPI_MODE_CREATE|MPI_MODE_APPEND, MPI_INFO_NULL, &file);
473
474 MPI_Comm_rank(comm, &rank);
475 MPI_Comm_size(comm, &size);
476
477 if (rank == 0) {
478
479 buf << "###########################################################" << std::endl
480 << "GLOBAL INFORMATION:" << std::endl
481 << " NAME: " << name << std::endl
482 << " LEVEL: " << grid.Level() << std::endl
483 << " COMM_SIZE_GLOBAL: " << size << std::endl
484 << " COMM_SIZE_LOCAL: " << size_local_max << std::endl
485 << " GLOBAL:" << std::endl
486 << " GLOBAL_FINER_BEGIN: " << grid.Global().GlobalFinerBegin() << std::endl
487 << " GLOBAL_FINER_END: " << grid.Global().GlobalFinerEnd() << std::endl
488 << " GLOBAL_FINER_SIZE: " << grid.Global().GlobalFinerSize() << std::endl
489 << " GLOBAL_SIZE: " << grid.Global().GlobalSize() << std::endl
490 << " EXTENT:" << std::endl
491 << " BEGIN: " << grid.Extent().Begin() << std::endl
492 << " END: " << grid.Extent().End() << std::endl
493 << " SIZE: " << grid.Extent().Size() << std::endl
494 << " MESH_WIDTH: " << grid.Extent().MeshWidth() << std::endl
495 << std::endl
496 << "LOCAL INFORMATION:" << std::endl;
497 }
498
499 buf << "RANK " << rank << ":" << std::endl
500 << " GLOBAL:" << std::endl
501 << " LOCAL_BEGIN: " << grid.Global().LocalBegin() << std::endl
502 << " LOCAL_END: " << grid.Global().LocalEnd() << std::endl
503 << " LOCAL_SIZE: " << grid.Global().LocalSize() << std::endl
504 << " LOCAL_FINER_BEGIN: " << grid.Global().LocalFinerBegin() << std::endl
505 << " LOCAL_FINER_END: " << grid.Global().LocalFinerEnd() << std::endl
506 << " LOCAL_FINER_SIZE: " << grid.Global().LocalFinerSize() << std::endl
507 << " BOUNDARY_TYPE: " << grid.Global().BoundaryType() << std::endl
508 << " LOCAL:" << std::endl
509 << " BEGIN: " << grid.Local().Begin() << std::endl
510 << " END: " << grid.Local().End() << std::endl
511 << " SIZE: " << grid.Local().Size() << std::endl
512 << " SIZE_TOTAL: " << grid.Local().SizeTotal() << std::endl
513 << " HALO_BEGIN_1: " << grid.Local().HaloBegin1() << std::endl
514 << " HALO_END_1: " << grid.Local().HaloEnd1() << std::endl
515 << " HALO_SIZE_1: " << grid.Local().HaloSize1() << std::endl
516 << " HALO_BEGIN_2: " << grid.Local().HaloBegin2() << std::endl
517 << " HALO_END_2: " << grid.Local().HaloEnd2() << std::endl
518 << " HALO_SIZE_2: " << grid.Local().HaloSize2() << std::endl
519 << " BOUNDARY_BEGIN_1: " << grid.Local().BoundaryBegin1() << std::endl
520 << " BOUNDARY_END_1: " << grid.Local().BoundaryEnd1() << std::endl
521 << " BOUNDARY_SIZE_1: " << grid.Local().BoundarySize1() << std::endl
522 << " BOUNDARY_BEGIN_2: " << grid.Local().BoundaryBegin2() << std::endl
523 << " BOUNDARY_END_2: " << grid.Local().BoundaryEnd2() << std::endl
524 << " BOUNDARY_SIZE_2: " << grid.Local().BoundarySize2() << std::endl
525 << " FINER_BEGIN: " << grid.Local().FinerBegin() << std::endl
526 << " FINER_END: " << grid.Local().FinerEnd() << std::endl
527 << " FINER_SIZE: " << grid.Local().FinerSize() << std::endl;
528
529 if (rank == size-1)
530 buf << "###########################################################" << std::endl;
531
532 char* char_buf = Helper::GetCharArray(buf.str());
533 MPI_File_write_ordered(file, char_buf, buf.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
534 delete [] char_buf;
535
536 MPI_File_close(&file);
537
538 }
539}
540
541void CommMPI::PrintAllSettings()
542{
543 std::stringstream buf;
544 MPI_File file;
545
546 const Multigrid& mg = *MG::GetSol();
547
548 buf << OutputPath() << "settings.txt";
549 char *filename = Helper::GetCharArray(buf.str());
550
551 MPI_File_open(comm_global, filename, MPI_MODE_WRONLY|MPI_MODE_CREATE, MPI_INFO_NULL, &file);
552 MPI_File_set_size(file, 0);
553 MPI_File_close(&file);
554
555 for (int i=mg.MinLevel(); i<=mg.MaxLevel(); ++i)
556 PrintGridInformation(mg(i), filename, "MULTIGRID");
557
558 for (int i=mg.MinLevel()+1; i<=mg.MaxLevel(); ++i)
559 PrintGridInformation(settings.CoarserGrid(mg(i)), filename, "COARSER_GRID");
560
561 for (int i=mg.MinLevel(); i<mg.MaxLevel(); ++i)
562 PrintGridInformation(settings.FinerGrid(mg(i)), filename, "FINER_GRID");
563
564 delete [] filename;
565
566}
567
568void CommMPI::PrintGrid(Grid& grid, const char* information)
569{
570 int output_count = OutputCount();
571
572#ifdef HAVE_VTK
573
574 if (settings.CommunicatorLocal(grid) != MPI_COMM_NULL) {
575
576 vtkSmartPointer<vtkImageData> image = vtkSmartPointer<vtkImageData>::New();
577 image->SetExtent(grid.Global().LocalBegin().X(), grid.Global().LocalEnd().X()-1,
578 grid.Global().LocalBegin().Y(), grid.Global().LocalEnd().Y()-1,
579 grid.Global().LocalBegin().Z(), grid.Global().LocalEnd().Z()-1);
580 image->SetSpacing(grid.Extent().MeshWidth().vec());
581 image->SetOrigin(grid.Extent().Begin().vec());
582 image->SetScalarTypeToDouble();
583 image->SetNumberOfScalarComponents(1);
584 image->AllocateScalars();
585 image->GetPointData()->GetScalars()->SetName(information);
586
587 Index begin, end, i;
588 for (int j=0; j<3; ++j) {
589 begin[j] = (grid.Local().BoundarySize1()[j] > 0 ? grid.Local().BoundaryBegin1()[j] : grid.Local().Begin()[j]);
590 end[j] = (grid.Local().BoundarySize2()[j] > 0 ? grid.Local().BoundaryEnd2()[j] : grid.Local().End()[j]);
591 }
592
593 for (i.X() = begin.X(); i.X() < end.X(); ++i.X())
594 for (i.Y() = begin.Y(); i.Y() < end.Y(); ++i.Y())
595 for (i.Z() = begin.Z(); i.Z() < end.Z(); ++i.Z())
596 image->SetScalarComponentFromDouble(i.X() - begin.X() + grid.Global().LocalBegin().X(),
597 i.Y() - begin.Y() + grid.Global().LocalBegin().Y(),
598 i.Z() - begin.Z() + grid.Global().LocalBegin().Z(),
599 0, grid.GetVal(i));
600
601 image->Update();
602
603 int rank, size;
604 MPI_Comm_rank(comm_global, &rank);
605 MPI_Comm_size(comm_global, &size);
606
607 std::stringstream filename;
608 filename << OutputPath() << std::setw(4) << std::setfill('0') << output_count << "_" << rank << ".vti";
609
610 vtkSmartPointer<vtkXMLImageDataWriter> writer = vtkSmartPointer<vtkXMLImageDataWriter>::New();
611 writer->SetFileName(filename.str().c_str());
612 writer->SetCompressorTypeToNone();
613 writer->SetDataModeToAscii();
614 writer->SetInput(image);
615 writer->Update();
616 writer->Write();
617
618 }
619
620#else /* HAVE_VTK */
621 Index i;
622 std::stringstream buf;
623
624 Index begin, end;
625 Index begin_local, end_local, begin_global, end_global;
626
627 CommToGhosts(grid);
628
629 for (int i=0; i<3; ++i) {
630 end[i] = grid.Local().End()[i] + (grid.Global().LocalEnd()[i] == grid.Global().GlobalSize()[i] ? 0 : grid.Local().HaloSize1()[i]);
631 end_local[i] = grid.Global().LocalEnd()[i] - (grid.Global().LocalEnd()[i] == grid.Global().GlobalSize()[i] ? 1 : 0);
632 }
633
634 begin = grid.Local().Begin();
635 begin_local = grid.Global().LocalBegin();
636 begin_global = 0;
637 end_global = grid.Global().GlobalSize()-1;
638
639 for (i.Z()=begin.Z(); i.Z()<end.Z(); ++i.Z())
640 for (i.Y()=begin.Y(); i.Y()<end.Y(); ++i.Y())
641 for (i.X()=begin.X(); i.X()<end.X(); ++i.X())
642 buf << std::scientific << grid.GetVal(i) << " ";
643
644 CreateOutputFiles(grid, buf, information,
645 begin_global, end_global,
646 begin_local, end_local,
647 output_count);
648#endif /* HAVE_VTK */
649}
650
651void CommMPI::PrintDefect(Grid& sol, Grid& rhs, const char* information)
652{
653 TempGrid *temp = MG::GetTempGrid();
654 temp->SetProperties(sol);
655 temp->ImportFromResidual(sol, rhs);
656 PrintGrid(*temp, information);
657}
658
659void CommMPI::CreateOutputFiles(const Grid& grid, const std::stringstream& serial_data, const char* information,
660 const Index& begin_global, const Index& end_global,
661 const Index& begin_local, const Index& end_local,
662 const int& output_count)
663{
664 MPI_Comm comm = settings.CommunicatorGlobal(grid);
665
666 if (comm != MPI_COMM_NULL) {
667
668 MPI_File file;
669 std::string conv_information = Helper::ReplaceWhitespaces(information, "_");
670
671 CreateParallelOutputFile(grid, comm, output_count, conv_information.c_str(),
672 begin_global, end_global, begin_local, end_local);
673
674 file = CreateSerialOutputFile(grid, comm, output_count, conv_information.c_str(),
675 begin_global, end_global, begin_local, end_local);
676
677 char *char_buf = Helper::GetCharArray(serial_data.str());
678 MPI_File_write(file, char_buf, serial_data.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
679 delete [] char_buf;
680
681 FinalizeSerialOutputFile(file);
682
683 }
684}
685
686void CommMPI::CreateParallelOutputFile(const Grid& grid, MPI_Comm& comm,
687 const int& output_count, const char* information,
688 const Index& begin_global, const Index& end_global,
689 const Index& begin_local, const Index& end_local)
690{
691 int rank;
692 MPI_File file;
693 char parallel_filename[513], serial_filename[513];
694 std::stringstream buf;
695
696 MPI_Comm_rank(comm, &rank);
697
698 sprintf(parallel_filename, "%s%04d.pvti", OutputPath().c_str(), output_count);
699 sprintf(serial_filename, "%04d_%d.vti", output_count, rank);
700
701 MPI_File_open(comm, parallel_filename, MPI_MODE_WRONLY|MPI_MODE_CREATE, MPI_INFO_NULL, &file);
702 MPI_File_set_size(file, 0);
703
704 if (rank == 0) {
705
706 buf << "<?xml version=\"1.0\"?>" << std::endl
707 << "<VTKFile type=\"PImageData\" version=\"0.1\" byte_order=\"LittleEndian\">" << std::endl
708 << " <PImageData WholeExtent=\"";
709
710 for (int i=0; i<3; ++i)
711 buf << begin_global[i] << " " << end_global[i] << " ";
712
713 buf << "\"" << std::endl
714 << " GhostLevel=\"0\" Origin=\"0 0 0\" Spacing=\"";
715
716 for (int i=0; i<3; ++i)
717 buf << grid.Extent().MeshWidth()[i] << " ";
718
719 buf << "\">" << std::endl
720 << " <PPointData Scalars=\"" << information << "\">" << std::endl
721 << " <PDataArray type=\"Float32\" Name=\"" << information << "\"/>" << std::endl
722 << " </PPointData>" << std::endl;
723
724 char* char_buf = Helper::GetCharArray(buf.str());
725
726 MPI_File_write_shared(file, char_buf, buf.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
727
728 delete [] char_buf;
729 }
730
731 buf.str("");
732
733 if ((end_local-begin_local).Product() > 0) {
734 buf << " <Piece Extent=\"";
735
736 for (int i=0; i<3; ++i)
737 buf << begin_local[i] << " " << end_local[i] << " ";
738
739 buf << "\" Source=\"" << serial_filename << "\"/>" << std::endl;
740 }
741
742 char* char_buf = Helper::GetCharArray(buf.str());
743
744 MPI_File_write_ordered(file, char_buf, buf.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
745
746 delete [] char_buf;
747
748 if (rank == 0) {
749
750 buf.str("");
751
752 buf << " </PImageData>" << std::endl
753 << "</VTKFile>" << std::endl;
754
755 char* char_buf = Helper::GetCharArray(buf.str());
756
757 MPI_File_write_shared(file, char_buf, buf.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
758
759 delete [] char_buf;
760 }
761
762 MPI_File_close(&file);
763}
764
765MPI_File CommMPI::CreateSerialOutputFile(const Grid& grid, MPI_Comm& comm,
766 const int& output_count, const char* information,
767 const Index& begin_global, const Index& end_global,
768 const Index& begin_local, const Index& end_local)
769{
770 char serial_filename[513];
771 int rank;
772 MPI_File file;
773 std::stringstream buf;
774
775 MPI_Comm_rank(comm_global, &rank);
776
777 sprintf(serial_filename, "%s%04d_%d.vti", OutputPath().c_str(), output_count, rank);
778
779 MPI_File_open(MPI_COMM_SELF, serial_filename, MPI_MODE_WRONLY|MPI_MODE_CREATE, MPI_INFO_NULL, &file);
780
781 buf << "<?xml version=\"1.0\"?>" << std::endl
782 << "<VTKFile type=\"ImageData\" version=\"0.1\" byte_order=\"LittleEndian\">" << std::endl
783 << " <ImageData WholeExtent=\"";
784
785 for (int i=0; i<3; ++i)
786 buf << begin_global[i] << " " << end_global[i] << " ";
787
788 buf << "\"" << std::endl
789 << " Origin=\"0 0 0\" Spacing=\"";
790
791 for (int i=0; i<3; ++i)
792 buf << grid.Extent().MeshWidth()[i] << " ";
793
794 buf << "\">" << std::endl
795 << " <Piece Extent=\"";
796
797 for (int i=0; i<3; ++i)
798 buf << begin_local[i] << " " << end_local[i] << " ";
799
800 buf << "\">" << std::endl
801 << " <PointData Scalars=\"" << information << "\">" << std::endl
802 << " <DataArray type=\"Float32\" Name=\"" << information << "\" format=\"ascii\">" << std::endl
803 << " ";
804
805 char* char_buf = Helper::GetCharArray(buf.str());
806 MPI_File_write(file, char_buf, buf.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
807 delete [] char_buf;
808
809 return file;
810}
811
812void CommMPI::FinalizeSerialOutputFile(MPI_File& file)
813{
814 std::stringstream buf;
815
816 buf << std::endl
817 << " </DataArray>" << std::endl
818 << " </PointData>" << std::endl
819 << " </Piece>" << std::endl
820 << " </ImageData>" << std::endl
821 << "</VTKFile>" << std::endl;
822
823 char* char_buf = Helper::GetCharArray(buf.str());
824 MPI_File_write(file, char_buf, buf.str().size(), MPI_CHAR, MPI_STATUS_IGNORE);
825 delete [] char_buf;
826
827 MPI_File_close(&file);
828}
829
830int CommMPI::GlobalRank() const
831{
832 int rank;
833 MPI_Comm_rank(comm_global, &rank);
834 return rank;
835}
836
837int CommMPI::GlobalSize() const
838{
839 int size;
840 MPI_Comm_size(comm_global, &size);
841 return size;
842}
843
844Index CommMPI::GlobalPos() const
845{
846 Index dims, periods, coords;
847 MPI_Cart_get(comm_global, 3, dims.vec(), periods.vec(), coords.vec());
848 return coords;
849}
850
851Index CommMPI::GlobalProcs() const
852{
853 Index dims, periods, coords;
854 MPI_Cart_get(comm_global, 3, dims.vec(), periods.vec(), coords.vec());
855 return dims;
856}
857
858int CommMPI::Rank(const Grid& grid) const
859{
860 int rank;
861 MPI_Comm comm = settings.CommunicatorLocal(grid);
862 assert(comm != MPI_COMM_NULL);
863 MPI_Comm_rank(comm, &rank);
864 return rank;
865}
866
867int CommMPI::Size(const Grid& grid) const
868{
869 int size;
870 MPI_Comm comm = settings.CommunicatorLocal(grid);
871 assert(comm != MPI_COMM_NULL);
872 MPI_Comm_size(comm, &size);
873 return size;
874}
875
876Index CommMPI::Pos(const Grid& grid) const
877{
878 Index dims, periods, coords;
879 MPI_Comm comm = settings.CommunicatorLocal(grid);
880 assert(comm != MPI_COMM_NULL);
881 MPI_Cart_get(comm, 3, dims.vec(), periods.vec(), coords.vec());
882 return coords;
883}
884
885Index CommMPI::Procs(const Grid& grid) const
886{
887 Index dims, periods, coords;
888 MPI_Comm comm = settings.CommunicatorLocal(grid);
889 assert(comm != MPI_COMM_NULL);
890 MPI_Cart_get(comm, 3, dims.vec(), periods.vec(), coords.vec());
891 return dims;
892}
893
894void CommMPI::InitCommMPI(const MPI_Comm& comm)
895{
896 int status, size, rank;
897 int dims[3] = {0, 0, 0};
898 int periods[3];
899
900 for (int i=0; i<3; ++i)
901 periods[i] = (BoundaryConditions()[i] == Periodic ? 1 : 0);
902
903 MPI_Comm_size(comm, &size);
904 MPI_Comm_rank(comm, &rank);
905
906 MPI_Topo_test(comm, &status);
907
908 if (status == MPI_CART) {
909
910 comm_global = comm;
911
912 }else {
913
914 const int log2 = Helper::log_2(size);
915
916 if (Helper::intpow(2, log2) == size) {
917 for (int i=0; i<3; ++i)
918 dims[i] = Helper::intpow(2, log2 / 3 + (log2%3 > i ? 1 : 0));
919 }else {
920 MPI_Dims_create(size, 3, dims);
921 }
922
923#ifdef OUTPUT_DEBUG
924 if (rank == 0)
925 std::printf("Process grid: %d %d %d\n", dims[0], dims[1], dims[2]);
926#endif
927
928 MPI_Cart_create(comm, 3, dims, periods, 1, &comm_global);
929
930 }
931
932 MPI_Info_create(&info);
933 char key[] = "no_locks";
934 char val[] = "true";
935 MPI_Info_set(info, key, val);
936
937}
938
939CommMPI::~CommMPI()
940{
941 MPI_Comm_free(&comm_global);
942#ifdef VMG_ONE_SIDED
943 if (win_created)
944 MPI_Win_free(&win);
945#endif
946 MPI_Info_free(&info);
947}
948
949Grid& CommMPI::GetCoarserGrid(Multigrid& multigrid)
950{
951 return settings.CoarserGrid(multigrid(multigrid.Level()));
952}
953
954Grid& CommMPI::GetFinerGrid(Multigrid& multigrid)
955{
956 return settings.FinerGrid(multigrid(multigrid.Level()));
957}
958
959std::string CommMPI::CreateOutputDirectory()
960{
961#ifdef HAVE_BOOST_FILESYSTEM
962 std::string path, unique_path;
963 std::stringstream unique_suffix;
964 int suffix_counter = 0;
965 char buffer[129];
966 time_t rawtime;
967 struct tm *timeinfo;
968 int path_size;
969 char* path_array;
970
971 if (GlobalRank() == 0) {
972
973 time(&rawtime);
974 timeinfo = localtime(&rawtime);
975 strftime(buffer, 128, "./output/%Y_%m_%d_%H_%M_%S/", timeinfo);
976 path = buffer;
977
978 if (!fs::exists("output"))
979 fs::create_directory("output");
980
981 unique_path = path;
982
983 while (fs::exists(unique_path.c_str())) {
984
985 unique_suffix.str("");
986 unique_suffix << "_" << suffix_counter++ << "/";
987
988 unique_path = path;
989 unique_path.replace(unique_path.size()-1, 1, unique_suffix.str());
990
991 }
992
993 fs::create_directory(unique_path.c_str());
994
995 path_size = unique_path.size() + 1;
996 path_array = Helper::GetCharArray(unique_path);
997
998 MPI_Bcast(&path_size, 1, MPI_INT, 0, comm_global);
999
1000 }else {
1001
1002 MPI_Bcast(&path_size, 1, MPI_INT, 0, comm_global);
1003 path_array = new char[path_size];
1004
1005 }
1006
1007 MPI_Bcast(path_array, path_size, MPI_CHAR, 0, comm_global);
1008
1009 unique_path = path_array;
1010
1011 delete [] path_array;
1012
1013 return unique_path;
1014
1015#else
1016
1017 return "./";
1018
1019#endif
1020}
1021
1022
1023#endif /* HAVE_MPI */
Note: See TracBrowser for help on using the repository browser.