GADGET-4
setcomm.h
Go to the documentation of this file.
1/*******************************************************************************
2 * \copyright This file is part of the GADGET4 N-body/SPH code developed
3 * \copyright by Volker Springel. Copyright (C) 2014-2020 by Volker Springel
4 * \copyright (vspringel@mpa-garching.mpg.de) and all contributing authors.
5 *******************************************************************************/
6
12#ifndef SETCOMM_H
13#define SETCOMM_H
14
15#include <mpi.h>
16#include <stdarg.h>
17#include <stdio.h>
18#include <stdlib.h>
19#include <string.h>
20#include <algorithm>
21
23{
24 public:
25 setcomm(MPI_Comm Comm) { initcomm(Comm); }
26 setcomm(const char *str)
27 {
28 /* do nothing in this case, because we need to delay the initialization until MPI_Init has been executed */
29 }
30
31 MPI_Comm Communicator;
32 int NTask;
34 int PTask;
35
37 int NumNodes = 0;
42 long long MemoryOnNode;
44
45 void initcomm(MPI_Comm Comm)
46 {
47 Communicator = Comm;
48 MPI_Comm_rank(Communicator, &ThisTask);
49 MPI_Comm_size(Communicator, &NTask);
50
51 for(PTask = 0; NTask > (1 << PTask); PTask++)
52 ;
53 }
54
55 void mpi_printf(const char *fmt, ...)
56 {
57 if(ThisTask == 0)
58 {
59 va_list l;
60 va_start(l, fmt);
61 vprintf(fmt, l);
62 // myflush(stdout);
63 va_end(l);
64 }
65 }
66
67 private:
68 struct node_data
69 {
70 int task, this_node, first_task_in_this_node;
71 int first_index, rank_in_node, tasks_in_node;
72 char name[MPI_MAX_PROCESSOR_NAME];
73 };
74 node_data loc_node, *list_of_nodes;
75
76 static bool system_compare_hostname(const node_data &a, const node_data &b)
77 {
78 int res = strcmp(a.name, b.name);
79 if(res < 0)
80 return true;
81 if(res > 0)
82 return false;
83 return a.task < b.task;
84 }
85
86 static bool system_compare_first_task(const node_data &a, const node_data &b)
87 {
88 if(a.first_task_in_this_node < b.first_task_in_this_node)
89 return true;
90 if(a.first_task_in_this_node > b.first_task_in_this_node)
91 return false;
92 return a.task < b.task;
93 }
94
95 static bool system_compare_task(const node_data &a, const node_data &b) { return a.task < b.task; }
96
97 public:
99 {
100 int len, nodes, i, no, rank, first_index;
101
102 MPI_Get_processor_name(loc_node.name, &len);
103 loc_node.task = ThisTask;
104
105 list_of_nodes = (node_data *)malloc(
106 sizeof(node_data) * NTask); /* Note: Internal memory allocation routines are not yet available when this function is called */
107
108 MPI_Allgather(&loc_node, sizeof(node_data), MPI_BYTE, list_of_nodes, sizeof(node_data), MPI_BYTE, Communicator);
109
110 std::sort(list_of_nodes, list_of_nodes + NTask, system_compare_hostname);
111
112 list_of_nodes[0].first_task_in_this_node = list_of_nodes[0].task;
113
114 for(i = 1, nodes = 1; i < NTask; i++)
115 {
116 if(strcmp(list_of_nodes[i].name, list_of_nodes[i - 1].name) != 0)
117 {
118 list_of_nodes[i].first_task_in_this_node = list_of_nodes[i].task;
119 nodes++;
120 }
121 else
122 list_of_nodes[i].first_task_in_this_node = list_of_nodes[i - 1].first_task_in_this_node;
123 }
124
125 std::sort(list_of_nodes, list_of_nodes + NTask, system_compare_first_task);
126
127 for(i = 0; i < NTask; i++)
128 list_of_nodes[i].tasks_in_node = 0;
129
130 for(i = 0, no = 0, rank = 0, first_index = 0; i < NTask; i++)
131 {
132 if(i ? list_of_nodes[i].first_task_in_this_node != list_of_nodes[i - 1].first_task_in_this_node : 0)
133 {
134 no++;
135 rank = 0;
136 first_index = i;
137 }
138
139 list_of_nodes[i].first_index = first_index;
140 list_of_nodes[i].this_node = no;
141 list_of_nodes[i].rank_in_node = rank++;
142 list_of_nodes[first_index].tasks_in_node++;
143 }
144
145 int max_count = 0;
146 int min_count = (1 << 30);
147
148 for(i = 0; i < NTask; i++)
149 {
150 list_of_nodes[i].tasks_in_node = list_of_nodes[list_of_nodes[i].first_index].tasks_in_node;
151
152 if(list_of_nodes[i].tasks_in_node > max_count)
153 max_count = list_of_nodes[i].tasks_in_node;
154 if(list_of_nodes[i].tasks_in_node < min_count)
155 min_count = list_of_nodes[i].tasks_in_node;
156 }
157
158 std::sort(list_of_nodes, list_of_nodes + NTask, system_compare_task);
159
160 TasksInThisNode = list_of_nodes[ThisTask].tasks_in_node;
161 RankInThisNode = list_of_nodes[ThisTask].rank_in_node;
162
163 ThisNode = list_of_nodes[ThisTask].this_node;
164
165 NumNodes = nodes;
166 MinTasksPerNode = min_count;
167 MaxTasksPerNode = max_count;
168
169 free(list_of_nodes);
170 }
171};
172
173#endif
long long MemoryOnNode
Definition: setcomm.h:42
int ThisNode
Definition: setcomm.h:36
setcomm(const char *str)
Definition: setcomm.h:26
void mpi_printf(const char *fmt,...)
Definition: setcomm.h:55
int ThisTask
Definition: setcomm.h:33
int NTask
Definition: setcomm.h:32
int RankInThisNode
Definition: setcomm.h:39
void determine_compute_nodes(void)
Definition: setcomm.h:98
int NumNodes
Definition: setcomm.h:37
int PTask
Definition: setcomm.h:34
int MaxTasksPerNode
Definition: setcomm.h:41
MPI_Comm Communicator
Definition: setcomm.h:31
int TasksInThisNode
Definition: setcomm.h:38
long long SharedMemoryOnNode
Definition: setcomm.h:43
int MinTasksPerNode
Definition: setcomm.h:40
void initcomm(MPI_Comm Comm)
Definition: setcomm.h:45
setcomm(MPI_Comm Comm)
Definition: setcomm.h:25