7.2
general documentation
cs_parall.h
Go to the documentation of this file.
1 #ifndef __CS_PARALL_H__
2 #define __CS_PARALL_H__
3 
4 /*============================================================================
5  * Functions dealing with parallelism
6  *============================================================================*/
7 
8 /*
9  This file is part of code_saturne, a general-purpose CFD tool.
10 
11  Copyright (C) 1998-2022 EDF S.A.
12 
13  This program is free software; you can redistribute it and/or modify it under
14  the terms of the GNU General Public License as published by the Free Software
15  Foundation; either version 2 of the License, or (at your option) any later
16  version.
17 
18  This program is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
20  FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
21  details.
22 
23  You should have received a copy of the GNU General Public License along with
24  this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
25  Street, Fifth Floor, Boston, MA 02110-1301, USA.
26 */
27 
28 /*----------------------------------------------------------------------------*/
29 
30 /*----------------------------------------------------------------------------
31  * Local headers
32  *----------------------------------------------------------------------------*/
33 
34 #include "cs_defs.h"
35 
36 /*----------------------------------------------------------------------------*/
37 
39 
40 /*=============================================================================
41  * Public function prototypes
42  *============================================================================*/
43 
44 /*----------------------------------------------------------------------------*/
51 /*----------------------------------------------------------------------------*/
52 
53 #if defined(HAVE_MPI_IN_PLACE)
54 
55 inline static void
57  const int n)
58 {
59  if (cs_glob_n_ranks > 1) {
60  MPI_Allreduce(MPI_IN_PLACE, cpt, n, CS_MPI_GNUM, MPI_SUM,
62  }
63 }
64 
65 #elif defined(HAVE_MPI)
66 
67 void
69  const int n);
70 
71 #else
72 
73 #define cs_parall_counter(_cpt, _n)
74 
75 #endif
76 
77 /*----------------------------------------------------------------------------*/
84 /*----------------------------------------------------------------------------*/
85 
86 #if defined(HAVE_MPI_IN_PLACE)
87 
88 inline static void
90  const int n)
91 {
92  if (cs_glob_n_ranks > 1) {
93  MPI_Allreduce(MPI_IN_PLACE, cpt, n, CS_MPI_LNUM, MPI_MAX,
95  }
96 }
97 
98 #elif defined(HAVE_MPI)
99 
100 void
102  const int n);
103 
104 #else
105 
106 #define cs_parall_counter_max(_cpt, _n)
107 
108 #endif
109 
110 /*----------------------------------------------------------------------------*/
118 /*----------------------------------------------------------------------------*/
119 
120 #if defined(HAVE_MPI_IN_PLACE)
121 
122 inline static void
124  cs_datatype_t datatype,
125  void *val)
126 {
127  if (cs_glob_n_ranks > 1) {
128  MPI_Allreduce(MPI_IN_PLACE, val, n, cs_datatype_to_mpi[datatype], MPI_SUM,
130  }
131 }
132 
133 #elif defined(HAVE_MPI)
134 
135 void
136 cs_parall_sum(int n,
137  cs_datatype_t datatype,
138  void *val);
139 
140 #else
141 
142 #define cs_parall_sum(_n, _datatype, _val) { };
143 
144 #endif
145 
146 /*----------------------------------------------------------------------------*/
155 /*----------------------------------------------------------------------------*/
156 
157 #if defined(HAVE_MPI_IN_PLACE)
158 
159 inline static void
161  cs_datatype_t datatype,
162  void *val)
163 {
164  if (cs_glob_n_ranks > 1) {
165  MPI_Allreduce(MPI_IN_PLACE, val, n, cs_datatype_to_mpi[datatype], MPI_MAX,
167  }
168 }
169 
170 #elif defined(HAVE_MPI)
171 
172 void
173 cs_parall_max(int n,
174  cs_datatype_t datatype,
175  void *val);
176 
177 #else
178 
179 #define cs_parall_max(_n, _datatype, _val);
180 
181 #endif
182 
183 /*----------------------------------------------------------------------------*/
192 /*----------------------------------------------------------------------------*/
193 
194 #if defined(HAVE_MPI_IN_PLACE)
195 
196 inline static void
198  cs_datatype_t datatype,
199  void *val)
200 {
201  if (cs_glob_n_ranks > 1) {
202  MPI_Allreduce(MPI_IN_PLACE, val, n, cs_datatype_to_mpi[datatype], MPI_MIN,
204  }
205 }
206 
207 #elif defined(HAVE_MPI)
208 
209 void
210 cs_parall_min(int n,
211  cs_datatype_t datatype,
212  void *val);
213 
214 #else
215 
216 #define cs_parall_min(_n, _datatype, _val);
217 
218 #endif
219 
220 /*----------------------------------------------------------------------------*/
231 /*----------------------------------------------------------------------------*/
232 
233 #if defined(HAVE_MPI)
234 
235 inline static void
236 cs_parall_bcast(int root_rank,
237  int n,
238  cs_datatype_t datatype,
239  void *val)
240 {
241  if (cs_glob_n_ranks > 1)
242  MPI_Bcast(val, n, cs_datatype_to_mpi[datatype], root_rank,
244 }
245 
246 #else
247 
248 #define cs_parall_bcast(_root_rank, _n, _datatype, _val);
249 
250 #endif
251 
252 /*----------------------------------------------------------------------------*/
268 /*----------------------------------------------------------------------------*/
269 
270 void
271 cs_parall_allgather_r(int n_elts,
272  int n_g_elts,
273  cs_real_t array[],
274  cs_real_t g_array[]);
275 
276 /*----------------------------------------------------------------------------*/
296 /*----------------------------------------------------------------------------*/
297 
298 void
300  int n_g_elts,
301  int stride,
302  cs_real_t o_key[],
303  cs_real_t array[],
304  cs_real_t g_array[]);
305 
306 /*----------------------------------------------------------------------------*/
324 /*----------------------------------------------------------------------------*/
325 
326 void
327 cs_parall_gather_r(int root_rank,
328  int n_elts,
329  int n_g_elts,
330  const cs_real_t array[],
331  cs_real_t g_array[]);
332 
333 /*----------------------------------------------------------------------------*/
355 /*----------------------------------------------------------------------------*/
356 
357 void
358 cs_parall_gather_ordered_r(int root_rank,
359  int n_elts,
360  int n_g_elts,
361  int stride,
362  cs_real_t o_key[],
363  cs_real_t array[],
364  cs_real_t g_array[]);
365 
366 /*----------------------------------------------------------------------------*/
384 /*----------------------------------------------------------------------------*/
385 
386 void
387 cs_parall_scatter_r(int root_rank,
388  int n_elts,
389  int n_g_elts,
390  const cs_real_t g_array[],
391  cs_real_t array[]);
392 
393 /*----------------------------------------------------------------------------*/
412 /*----------------------------------------------------------------------------*/
413 
414 void
415 cs_parall_gather_f(int root_rank,
416  int n_elts,
417  int n_g_elts,
418  const float array[],
419  float g_array[]);
420 
421 /*----------------------------------------------------------------------------*/
440 /*----------------------------------------------------------------------------*/
441 
442 void
443 cs_parall_scatter_f(int root_rank,
444  int n_elts,
445  int n_g_elts,
446  const float g_array[],
447  float array[]);
448 
449 /*----------------------------------------------------------------------------*/
459 /*----------------------------------------------------------------------------*/
460 
461 void
463  cs_real_t *max,
464  cs_real_t max_loc_vals[]);
465 
466 /*----------------------------------------------------------------------------*/
476 /*----------------------------------------------------------------------------*/
477 
478 void
480  cs_real_t *min,
481  cs_real_t min_loc_vals[]);
482 
483 /*----------------------------------------------------------------------------*/
494 /*----------------------------------------------------------------------------*/
495 
496 void
498  int *rank_id,
499  cs_real_t val);
500 
501 /*----------------------------------------------------------------------------*/
510 /*----------------------------------------------------------------------------*/
511 
512 size_t
514 
515 /*----------------------------------------------------------------------------*/
525 /*----------------------------------------------------------------------------*/
526 
527 void
528 cs_parall_set_min_coll_buf_size(size_t buffer_size);
529 
530 /*----------------------------------------------------------------------------*/
543 /*----------------------------------------------------------------------------*/
544 
545 inline static void
547  size_t type_size,
548  cs_lnum_t *s_id,
549  cs_lnum_t *e_id)
550 {
551 #if defined(HAVE_OPENMP)
552  const int t_id = omp_get_thread_num();
553  const int n_t = omp_get_num_threads();
554  const cs_lnum_t t_n = (n + n_t - 1) / n_t;
555  const cs_lnum_t cl_m = CS_CL_SIZE / type_size; /* Cache line multiple */
556 
557  *s_id = t_id * t_n;
558  *e_id = (t_id+1) * t_n;
559  *s_id = cs_align(*s_id, cl_m);
560  *e_id = cs_align(*e_id, cl_m);
561  if (*e_id > n) *e_id = n;
562 #else
563  CS_UNUSED(type_size); /* avoid compiler warning */
564  *s_id = 0;
565  *e_id = n;
566 #endif
567 }
568 
569 /*----------------------------------------------------------------------------*/
570 
572 
573 #endif /* __CS_PARALL_H__ */
static void cs_parall_counter(cs_gnum_t cpt[], const int n)
Sum values of a counter on all default communicator processes.
Definition: cs_parall.h:56
cs_datatype_t
Definition: cs_defs.h:275
unsigned long cs_gnum_t
global mesh entity number
Definition: cs_defs.h:301
void cs_parall_gather_r(int root_rank, int n_elts, int n_g_elts, const cs_real_t array[], cs_real_t g_array[])
Build a global array on the given root rank from all local arrays.
Definition: cs_parall.c:1008
#define CS_CL_SIZE
Definition: cs_defs.h:466
void cs_parall_allgather_ordered_r(int n_elts, int n_g_elts, int stride, cs_real_t o_key[], cs_real_t array[], cs_real_t g_array[])
Build an ordered global array from each local array in each domain.
Definition: cs_parall.c:964
#define BEGIN_C_DECLS
Definition: cs_defs.h:510
void cs_parall_scatter_r(int root_rank, int n_elts, int n_g_elts, const cs_real_t g_array[], cs_real_t array[])
Distribute a global array from a given root rank over all ranks. Each rank receive the part related t...
Definition: cs_parall.c:1124
static void cs_parall_thread_range(cs_lnum_t n, size_t type_size, cs_lnum_t *s_id, cs_lnum_t *e_id)
Compute array index bounds for a local thread. When called inside an OpenMP parallel section...
Definition: cs_parall.h:546
#define CS_UNUSED(x)
Definition: cs_defs.h:496
void cs_parall_gather_f(int root_rank, int n_elts, int n_g_elts, const float array[], float g_array[])
Build a global array on the given root rank from all local arrays. Function dealing with single-preci...
Definition: cs_parall.c:1189
static void cs_parall_counter_max(cs_lnum_t cpt[], const int n)
Maximum values of a counter on all default communicator processes.
Definition: cs_parall.h:89
int cs_glob_n_ranks
Definition: cs_defs.c:175
static void cs_parall_bcast(int root_rank, int n, cs_datatype_t datatype, void *val)
Broadcast values of a given datatype to all default communicator processes.
Definition: cs_parall.h:236
double cs_real_t
Floating-point value.
Definition: cs_defs.h:322
void cs_parall_min_loc_vals(int n, cs_real_t *min, cs_real_t min_loc_vals[])
Minimum value of a real and the value of related array on all default communicator processes...
Definition: cs_parall.c:794
void cs_parall_min_id_rank_r(cs_lnum_t *elt_id, int *rank_id, cs_real_t val)
Given an (id, rank, value) tuple, return the local id and rank corresponding to the global minimum va...
Definition: cs_parall.c:833
MPI_Datatype cs_datatype_to_mpi[]
Definition: cs_defs.c:157
size_t cs_parall_get_min_coll_buf_size(void)
Return minimum recommended scatter or gather buffer size.
Definition: cs_parall.c:1309
void cs_parall_scatter_f(int root_rank, int n_elts, int n_g_elts, const float g_array[], float array[])
Distribute a global array from a given root rank over all ranks. Each rank receive the part related t...
Definition: cs_parall.c:1254
#define CS_MPI_GNUM
Definition: cs_defs.h:386
void cs_parall_gather_ordered_r(int root_rank, int n_elts, int n_g_elts, int stride, cs_real_t o_key[], cs_real_t array[], cs_real_t g_array[])
Build an ordered global array on the given root rank from all local arrays.
Definition: cs_parall.c:1075
void cs_parall_set_min_coll_buf_size(size_t buffer_size)
Define minimum recommended scatter or gather buffer size.
Definition: cs_parall.c:1331
static void cs_parall_max(int n, cs_datatype_t datatype, void *val)
Maximum values of a given datatype on all default communicator processes.
Definition: cs_parall.h:160
MPI_Comm cs_glob_mpi_comm
Definition: cs_defs.c:183
int cs_lnum_t
local mesh entity id
Definition: cs_defs.h:316
void cs_parall_max_loc_vals(int n, cs_real_t *max, cs_real_t max_loc_vals[])
Maximum value of a real and the value of related array on all default communicator processes...
Definition: cs_parall.c:756
#define CS_MPI_LNUM
Definition: cs_defs.h:406
#define END_C_DECLS
Definition: cs_defs.h:511
static void cs_parall_sum(int n, cs_datatype_t datatype, void *val)
Sum values of a given datatype on all default communicator processes.
Definition: cs_parall.h:123
void cs_parall_allgather_r(int n_elts, int n_g_elts, cs_real_t array[], cs_real_t g_array[])
Build a global array from each local array in each domain.
Definition: cs_parall.c:887
static void cs_parall_min(int n, cs_datatype_t datatype, void *val)
Minimum values of a given datatype on all default communicator processes.
Definition: cs_parall.h:197
static cs_lnum_t cs_align(cs_lnum_t i, cs_lnum_t m)
Given a base index i, return the next index aligned with a size m.
Definition: cs_defs.h:592