ug4
Loading...
Searching...
No Matches
space_time_communicator.hpp
Go to the documentation of this file.
1/*
2 * Copyright (c) 2025: G-CSC, Goethe University Frankfurt
3 * Author: Arne Nägel, Martin Parnet
4
5 * Copyright (c) 2026: KAUST, King Abdullah University of Science and Technology
6 * Author: Daniel Gonzalez
7 *
8 * This file is part of UG4.
9 *
10 * UG4 is free software: you can redistribute it and/or modify it under the
11 * terms of the GNU Lesser General Public License version 3 (as published by the
12 * Free Software Foundation) with the following additional attribution
13 * requirements (according to LGPL/GPL v3 §7):
14 *
15 * (1) The following notice must be displayed in the Appropriate Legal Notices
16 * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
17 *
18 * (2) The following notice must be displayed at a prominent place in the
19 * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
20 *
21 * (3) The following bibliography is recommended for citation and must be
22 * preserved in all covered files:
23 * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
24 * parallel geometric multigrid solver on hierarchically distributed grids.
25 * Computing and visualization in science 16, 4 (2013), 151-164"
26 * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
27 * flexible software system for simulating pde based models on high performance
28 * computers. Computing and visualization in science 16, 4 (2013), 165-179"
29 *
30 * This program is distributed in the hope that it will be useful,
31 * but WITHOUT ANY WARRANTY; without even the implied warranty of
32 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
33 * GNU Lesser General Public License for more details.
34 */
35
36#ifndef __H__PCL__PCL_SPACE_TIME_COMMUNICATOR_H
37#define __H__PCL__PCL_SPACE_TIME_COMMUNICATOR_H
38
39
40#include "pcl/pcl_comm_world.h"
41
42namespace pcl{
43
45 public:
46
47 //--------------------------------------------------------------------------------------------------------------
48
50 virtual ~SpaceTimeCommunicator() = default;
51
52 //--------------------------------------------------------------------------------------------------------------
53
54 void split(int numTemporalProcesses)
55 {
56
57 int world_size, myid;
58 MPI_Comm_size(PCL_COMM_WORLD, &world_size);
59
60 globalsize_ = world_size;
61 spatialsize_ = world_size / numTemporalProcesses;
62 temporalsize_ = numTemporalProcesses;
63
64 if(world_size % numTemporalProcesses != 0 )
65 UG_THROW("SpaceTimeCommunicator: world_size is not a multiple of # temporal processes:"
66 " Change the number of processors to "
67 << numTemporalProcesses * spatialsize_ <<" or "
68 << (numTemporalProcesses) * (spatialsize_ + 1 )<<".\n");
69
70
72
73
74
75 MPI_Comm_rank(GLOBAL, &myid);
76 const int xcolor = myid / spatialsize_;
77 const int tcolor = myid % spatialsize_;
78
79 MPI_Comm_split(GLOBAL, xcolor, myid, &SPATIAL);
80 MPI_Comm_split(GLOBAL, tcolor, myid, &TEMPORAL);
81
82 if (verbose_)
83 {
84 UG_LOG("World size after splitting is:\t" << world_size );
85 UG_LOG( "... with temporal world size:\t" << temporalsize_ );
86 UG_LOG( "... and spatial world size:\t" << spatialsize_ );
87 }
88
89
90 PCL_COMM_WORLD = SPATIAL; // replaces ugs world communicator with the communicator for spatial
91 }
92
93 void unsplit() {
94 PCL_COMM_WORLD = GLOBAL; // reset the world communicator
95
96 MPI_Comm_free(&SPATIAL); // free the spatial communicator
98
99 MPI_Comm_free(&TEMPORAL);// free the temporal communicator
101 }
102
103 int get_global_size() const {
104 return globalsize_;
105 }
106
107 int get_temporal_size() const {
108 return temporalsize_;
109 }
110
111 int get_spatial_size() const {
112 return spatialsize_;
113 }
114
115 int get_temporal_rank() const {
116 int rank = 0;
117 MPI_Comm_rank(TEMPORAL, &rank);
118 return rank;
119 }
120
121 int get_spatial_rank() const {
122 int rank = 0;
123 MPI_Comm_rank(SPATIAL, &rank);
124 return rank;
125 }
126
127 int get_global_rank() const {
128 int rank = 0;
129 MPI_Comm_rank(GLOBAL, &rank);
130 return rank;
131 }
132
133 void sleep(int microseconds) {
134 usleep(microseconds);
135 }
136
137 //--------------------------------------------------------------------------------------------------------------
141
142 bool verbose_ = true;
143 int globalsize_ = 1;
146
147 //--------------------------------------------------------------------------------------------------------------
148 };
149}
150#endif
Definition space_time_communicator.hpp:44
MPI_Comm TEMPORAL
Definition space_time_communicator.hpp:139
int temporalsize_
Definition space_time_communicator.hpp:144
int get_global_size() const
Definition space_time_communicator.hpp:103
int globalsize_
Definition space_time_communicator.hpp:143
void split(int numTemporalProcesses)
Definition space_time_communicator.hpp:54
int get_temporal_size() const
Definition space_time_communicator.hpp:107
void sleep(int microseconds)
Definition space_time_communicator.hpp:133
void unsplit()
Definition space_time_communicator.hpp:93
int spatialsize_
Definition space_time_communicator.hpp:145
int get_spatial_size() const
Definition space_time_communicator.hpp:111
int get_spatial_rank() const
Definition space_time_communicator.hpp:121
int get_temporal_rank() const
Definition space_time_communicator.hpp:115
MPI_Comm GLOBAL
Definition space_time_communicator.hpp:138
MPI_Comm SPATIAL
Definition space_time_communicator.hpp:140
bool verbose_
Definition space_time_communicator.hpp:142
virtual ~SpaceTimeCommunicator()=default
int get_global_rank() const
Definition space_time_communicator.hpp:127
#define UG_THROW(msg)
Definition error.h:57
#define UG_LOG(msg)
Definition log.h:367
Definition parallel_grid_layout.h:46
MPI_Comm PCL_COMM_WORLD
Definition pcl_comm_world.cpp:34