9.1
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-2025 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 "base/cs_defs.h"
36
37/*----------------------------------------------------------------------------*/
38
40
41/*============================================================================
42 * General types and macros used throughout code_saturne
43 *============================================================================*/
44
45/*----------------------------------------------------------------------------
46 * Variable value type.
47 *----------------------------------------------------------------------------*/
48
52typedef enum {
53
63
64/*============================================================================
65 * Global variables
66 *============================================================================*/
67
68/* Preferred indexed sum option, adapted to shared-memory parallelism */
69
71
72/*=============================================================================
73 * Public function prototypes
74 *============================================================================*/
75
76/*----------------------------------------------------------------------------*/
83/*----------------------------------------------------------------------------*/
84
85inline static void
87 const int n)
88{
89#if defined(HAVE_MPI)
90
91 if (cs_glob_n_ranks > 1) {
92 MPI_Allreduce(MPI_IN_PLACE, cpt, n, CS_MPI_GNUM, MPI_SUM,
93 cs_glob_mpi_comm);
94 }
95
96#else
97
98 CS_UNUSED(cpt);
99 CS_UNUSED(n);
100
101#endif
102}
103
104/*----------------------------------------------------------------------------*/
111/*----------------------------------------------------------------------------*/
112
113inline static void
115 const int n)
116{
117#if defined(HAVE_MPI)
118
119 if (cs_glob_n_ranks > 1) {
120 MPI_Allreduce(MPI_IN_PLACE, cpt, n, CS_MPI_LNUM, MPI_MAX,
121 cs_glob_mpi_comm);
122 }
123
124#else
125
126 CS_UNUSED(cpt);
127 CS_UNUSED(n);
128
129#endif
130}
131
132/*----------------------------------------------------------------------------*/
140/*----------------------------------------------------------------------------*/
141
142inline static void
144 cs_datatype_t datatype,
145 void *val)
146{
147#if defined(HAVE_MPI)
148
149 if (cs_glob_n_ranks > 1) {
150 MPI_Allreduce(MPI_IN_PLACE, val, n, cs_datatype_to_mpi[datatype], MPI_SUM,
151 cs_glob_mpi_comm);
152 }
153
154#else
155
156 CS_UNUSED(n);
157 CS_UNUSED(datatype);
158 CS_UNUSED(val);
159
160#endif
161}
162
163/*----------------------------------------------------------------------------*/
172/*----------------------------------------------------------------------------*/
173
174inline static void
176 cs_datatype_t datatype,
177 void *val)
178{
179#if defined(HAVE_MPI)
180
181 if (cs_glob_n_ranks > 1) {
182 MPI_Allreduce(MPI_IN_PLACE, val, n, cs_datatype_to_mpi[datatype], MPI_MAX,
183 cs_glob_mpi_comm);
184 }
185
186#else
187
188 CS_UNUSED(n);
189 CS_UNUSED(datatype);
190 CS_UNUSED(val);
191
192#endif
193}
194
195/*----------------------------------------------------------------------------*/
204/*----------------------------------------------------------------------------*/
205
206inline static void
208 cs_datatype_t datatype,
209 void *val)
210{
211#if defined(HAVE_MPI)
212
213 if (cs_glob_n_ranks > 1) {
214 MPI_Allreduce(MPI_IN_PLACE, val, n, cs_datatype_to_mpi[datatype], MPI_MIN,
215 cs_glob_mpi_comm);
216 }
217
218#else
219
220 CS_UNUSED(n);
221 CS_UNUSED(datatype);
222 CS_UNUSED(val);
223
224#endif
225}
226
227/*----------------------------------------------------------------------------*/
238/*----------------------------------------------------------------------------*/
239
240inline static void
241cs_parall_bcast(int root_rank,
242 int n,
243 cs_datatype_t datatype,
244 void *val)
245{
246#if defined(HAVE_MPI)
247
248 if (cs_glob_n_ranks > 1)
249 MPI_Bcast(val, n, cs_datatype_to_mpi[datatype], root_rank,
250 cs_glob_mpi_comm);
251
252#else
253
254 CS_UNUSED(root_rank);
255 CS_UNUSED(n);
256 CS_UNUSED(datatype);
257 CS_UNUSED(val);
258
259#endif
260}
261
262/*----------------------------------------------------------------------------*/
278/*----------------------------------------------------------------------------*/
279
280void
281cs_parall_allgather_r(int n_elts,
282 int n_g_elts,
283 cs_real_t array[],
284 cs_real_t g_array[]);
285
286/*----------------------------------------------------------------------------*/
306/*----------------------------------------------------------------------------*/
307
308void
310 int n_g_elts,
311 int stride,
312 cs_real_t o_key[],
313 cs_real_t array[],
314 cs_real_t g_array[]);
315
316/*----------------------------------------------------------------------------*/
334/*----------------------------------------------------------------------------*/
335
336void
337cs_parall_gather_r(int root_rank,
338 int n_elts,
339 int n_g_elts,
340 const cs_real_t array[],
341 cs_real_t g_array[]);
342
343/*----------------------------------------------------------------------------*/
365/*----------------------------------------------------------------------------*/
366
367void
368cs_parall_gather_ordered_r(int root_rank,
369 int n_elts,
370 int n_g_elts,
371 int stride,
372 cs_real_t o_key[],
373 cs_real_t array[],
374 cs_real_t g_array[]);
375
376/*----------------------------------------------------------------------------*/
394/*----------------------------------------------------------------------------*/
395
396void
397cs_parall_scatter_r(int root_rank,
398 int n_elts,
399 int n_g_elts,
400 const cs_real_t g_array[],
401 cs_real_t array[]);
402
403/*----------------------------------------------------------------------------*/
422/*----------------------------------------------------------------------------*/
423
424void
425cs_parall_gather_f(int root_rank,
426 int n_elts,
427 int n_g_elts,
428 const float array[],
429 float g_array[]);
430
431/*----------------------------------------------------------------------------*/
450/*----------------------------------------------------------------------------*/
451
452void
453cs_parall_scatter_f(int root_rank,
454 int n_elts,
455 int n_g_elts,
456 const float g_array[],
457 float array[]);
458
459/*----------------------------------------------------------------------------*/
469/*----------------------------------------------------------------------------*/
470
471void
473 cs_real_t *max,
474 cs_real_t max_loc_vals[]);
475
476/*----------------------------------------------------------------------------*/
486/*----------------------------------------------------------------------------*/
487
488void
490 cs_real_t *min,
491 cs_real_t min_loc_vals[]);
492
493/*----------------------------------------------------------------------------*/
504/*----------------------------------------------------------------------------*/
505
506void
508 int *rank_id,
509 cs_real_t val);
510
511/*----------------------------------------------------------------------------*/
520/*----------------------------------------------------------------------------*/
521
522size_t
524
525/*----------------------------------------------------------------------------*/
535/*----------------------------------------------------------------------------*/
536
537void
538cs_parall_set_min_coll_buf_size(size_t buffer_size);
539
540/*----------------------------------------------------------------------------*/
550/*----------------------------------------------------------------------------*/
551
552inline static int
554 cs_lnum_t min_thread_elements)
555{
556#if defined(HAVE_OPENMP)
557 int n_t = omp_get_max_threads();
558 int n_t_l = n_elements / min_thread_elements;
559 if (n_t_l < n_t)
560 n_t = n_t_l;
561 if (n_t < 1)
562 n_t = 1;
563 return n_t;
564#else
565 CS_UNUSED(n_elements); /* avoid compiler warning */
566 CS_UNUSED(min_thread_elements);
567 return 1;
568#endif
569}
570
571/*----------------------------------------------------------------------------*/
584/*----------------------------------------------------------------------------*/
585
586inline static void
588 size_t type_size,
589 cs_lnum_t *s_id,
590 cs_lnum_t *e_id)
591{
592#if defined(HAVE_OPENMP)
593 const int t_id = omp_get_thread_num();
594 const int n_t = omp_get_num_threads();
595 const cs_lnum_t t_n = (n + n_t - 1) / n_t;
596 const cs_lnum_t cl_m = CS_CL_SIZE / type_size; /* Cache line multiple */
597
598 *s_id = t_id * t_n;
599 *e_id = (t_id+1) * t_n;
600 *s_id = cs_align(*s_id, cl_m);
601 *e_id = cs_align(*e_id, cl_m);
602 if (*e_id > n) *e_id = n;
603#else
604 CS_UNUSED(type_size); /* avoid compiler warning */
605 *s_id = 0;
606 *e_id = n;
607#endif
608}
609
610/*----------------------------------------------------------------------------*/
629/*----------------------------------------------------------------------------*/
630
631inline static void
633 size_t type_size,
634 cs_lnum_t *s_id,
635 cs_lnum_t *e_id)
636{
637#if defined(HAVE_OPENMP)
638 const int t_id = omp_get_thread_num();
639 const double n_t = omp_get_num_threads();
640 const cs_lnum_t cl_m = CS_CL_SIZE / type_size; /* Cache line multiple */
641
642 double r0 = (double)t_id / (double)n_t;
643 double r1 = (double)(t_id+1) / (double)n_t;
644
645 r0 = r0*r0;
646 r1 = r1*r1;
647
648 const cs_lnum_t t_0 = r0*n;
649 const cs_lnum_t t_1 = r1*n;
650
651 *s_id = t_0 * n;
652 *e_id = t_1 * n;
653 *s_id = cs_align(*s_id, cl_m);
654 *e_id = cs_align(*e_id, cl_m);
655 if (*e_id > n) *e_id = n;
656#else
657 CS_UNUSED(type_size); /* avoid compiler warning */
658 *s_id = 0;
659 *e_id = n;
660#endif
661}
662
663/*----------------------------------------------------------------------------*/
672/*----------------------------------------------------------------------------*/
673
674static inline size_t
676 size_t block_size)
677{
678 return (n % block_size) ? n/block_size + 1 : n/block_size;
679}
680
681/*----------------------------------------------------------------------------*/
682
684
685#if defined(__cplusplus)
686
687/*=============================================================================
688 * Public C++ functions
689 *============================================================================*/
690
691/*----------------------------------------------------------------------------*/
699/*----------------------------------------------------------------------------*/
700
701inline static void
702cs_parall_counter([[maybe_unused]] const cs_execution_context *ec,
703 [[maybe_unused]] cs_gnum_t cpt[],
704 [[maybe_unused]] const int n)
705{
706#if defined(HAVE_MPI)
707
708 if (ec->use_mpi()) {
709 MPI_Allreduce(MPI_IN_PLACE, cpt, n, CS_MPI_GNUM, MPI_SUM,
710 ec->comm());
711 }
712
713#endif
714}
715
716/*----------------------------------------------------------------------------*/
724/*----------------------------------------------------------------------------*/
725
726inline static void
728 [[maybe_unused]] cs_lnum_t cpt[],
729 [[maybe_unused]] const int n)
730{
731#if defined(HAVE_MPI)
732
733 if (ec->use_mpi()) {
734 MPI_Allreduce(MPI_IN_PLACE, cpt, n, CS_MPI_LNUM, MPI_MAX,
735 ec->comm());
736 }
737
738#endif
739}
740
741/*----------------------------------------------------------------------------*/
750/*----------------------------------------------------------------------------*/
751
752inline static void
753cs_parall_sum([[maybe_unused]] const cs_execution_context *ec,
754 [[maybe_unused]] int n,
755 [[maybe_unused]] cs_datatype_t datatype,
756 [[maybe_unused]] void *val)
757{
758#if defined(HAVE_MPI)
759
760 if (ec->use_mpi()) {
761 MPI_Allreduce(MPI_IN_PLACE, val, n,
762 cs_datatype_to_mpi[datatype], MPI_SUM,
763 ec->comm());
764 }
765
766#endif
767}
768
769/*----------------------------------------------------------------------------*/
779/*----------------------------------------------------------------------------*/
780
781inline static void
782cs_parall_max([[maybe_unused]] const cs_execution_context *ec,
783 [[maybe_unused]] int n,
784 [[maybe_unused]] cs_datatype_t datatype,
785 [[maybe_unused]] void *val)
786{
787#if defined(HAVE_MPI)
788
789 if (ec->use_mpi()) {
790 MPI_Allreduce(MPI_IN_PLACE, val, n,
791 cs_datatype_to_mpi[datatype], MPI_MAX,
792 ec->comm());
793 }
794
795#endif
796}
797
798/*----------------------------------------------------------------------------*/
808/*----------------------------------------------------------------------------*/
809
810inline static void
811cs_parall_min([[maybe_unused]] const cs_execution_context *ec,
812 [[maybe_unused]] int n,
813 [[maybe_unused]] cs_datatype_t datatype,
814 [[maybe_unused]] void *val)
815{
816#if defined(HAVE_MPI)
817
818 if (ec->use_mpi()) {
819 MPI_Allreduce(MPI_IN_PLACE, val, n,
820 cs_datatype_to_mpi[datatype], MPI_MIN,
821 ec->comm());
822 }
823
824#endif
825}
826
827/*----------------------------------------------------------------------------*/
842/*----------------------------------------------------------------------------*/
843
844inline static void
846 size_t type_size,
847 int t_id,
848 int n_t,
849 cs_lnum_t *s_id,
850 cs_lnum_t *e_id)
851{
852 const cs_lnum_t t_n = (n + n_t - 1) / n_t;
853 const cs_lnum_t cl_m = CS_CL_SIZE / type_size; /* Cache line multiple */
854
855 *s_id = t_id * t_n;
856 *e_id = (t_id+1) * t_n;
857 *s_id = cs_align(*s_id, cl_m);
858 *e_id = cs_align(*e_id, cl_m);
859 if (*s_id > n) *s_id = n;
860 if (*e_id > n) *e_id = n;
861}
862
863/*=============================================================================
864 * Public C++ templates
865 *============================================================================*/
866
867namespace cs {
868namespace parall {
869
870/*----------------------------------------------------------------------------*/
876/*----------------------------------------------------------------------------*/
877
878template <typename T, typename... Vals>
879static void
881(
882 const cs_execution_context *ec,
883 T& first,
884 Vals&... values
885)
886{
887#if defined(HAVE_MPI)
888 /* If no parallel ranks exit the function */
889 if (!ec->use_mpi())
890 return;
891
892 /* Count number of values */
893 constexpr int n_vals = sizeof...(Vals);
894
895 /* Set datatype for global communication */
896 cs_datatype_t datatype = cs_datatype_from_type<T>();
897
898 /* Temporary work array and parallel sum */
899 if (n_vals == 0)
900 cs_parall_sum(ec, 1, datatype, &first);
901 else {
902 /* Unpack values */
903 T *_values[] = {&values ...};
904
905 T w[n_vals + 1];
906 w[0] = first;
907 for (int i = 0; i < n_vals; i++)
908 w[i+1] = *(_values[i]);
909
910 cs_parall_sum(ec, n_vals + 1, datatype, w);
911
912 first = w[0];
913 for (int i = 0; i < n_vals; i++)
914 *(_values[i]) = w[i+1];
915 }
916#else
917 return;
918#endif
919}
920
921/*----------------------------------------------------------------------------*/
927/*----------------------------------------------------------------------------*/
928
929template <typename T, typename... Vals>
930static void
932(
933 T& first,
934 Vals&... values
935)
936{
937#if defined(HAVE_MPI)
939 sum(ec, first, values...);
940#else
941 return;
942#endif
943}
944
945/*----------------------------------------------------------------------------*/
952/*----------------------------------------------------------------------------*/
953
954template <int Stride, typename T, typename... Vals>
955static void
957(
958 const cs_execution_context *ec,
959 T first[],
960 Vals&&... values
961)
962{
963#if defined(HAVE_MPI)
964
965 if (!ec->use_mpi())
966 return;
967
968 /* Count number of values */
969 constexpr int n_vals = sizeof...(Vals);
970
971 /* Set datatype for global communication */
972 cs_datatype_t datatype = cs_datatype_from_type<T>();
973
974 /* Temporary work array and parallel sum */
975 if (n_vals == 0)
976 cs_parall_sum(ec, Stride, datatype, first);
977 else {
978 /* Unpack values */
979 T *_values[] = {values ...};
980
981 constexpr int work_size = (n_vals + 1) * Stride;
982
983 T w[work_size];
984 for (int i = 0; i < Stride; i++)
985 w[i] = first[i];
986
987 for (int i = 0; i < n_vals; i++) {
988 for (int j = 0; j < Stride; j++)
989 w[(i+1)*Stride + j] = _values[i][j];
990 }
991
992 cs_parall_sum(ec, work_size, datatype, w);
993
994 for (int i = 0; i < Stride; i++)
995 first[i] = w[i];
996
997 for (int i = 0; i < n_vals; i++) {
998 for (int j = 0; j < Stride; j++)
999 _values[i][j] = w[(i+1)*Stride + j];
1000 }
1001 }
1002
1003#endif
1004}
1005
1006/*----------------------------------------------------------------------------*/
1014/*----------------------------------------------------------------------------*/
1015
1016template <int Stride, typename T, typename... Vals>
1017static void
1019(
1020 T first[],
1021 Vals&&... values
1022)
1023{
1024#if defined(HAVE_MPI)
1025 auto *ec = cs_execution_context_glob_get();
1026 sum<Stride>(ec, first, values...);
1027#else
1028 return;
1029#endif
1030}
1031
1032/*----------------------------------------------------------------------------*/
1039/*----------------------------------------------------------------------------*/
1040
1041template <typename T, typename... Vals>
1042static void
1044(
1045 const cs_execution_context *ec,
1046 T& first,
1047 Vals&... values
1048)
1049{
1050#if defined(HAVE_MPI)
1051
1052 /* Count number of values */
1053 constexpr int n_vals = sizeof...(Vals);
1054
1055 /* Set datatype for global communication */
1056 cs_datatype_t datatype = cs_datatype_from_type<T>();
1057
1058 /* Temporary work array and parallel sum */
1059 if (n_vals == 0)
1060 cs_parall_max(ec, 1, datatype, &first);
1061 else {
1062
1063 /* Unpack values */
1064 T *_values[] = {&values ...};
1065
1066 T w[n_vals + 1];
1067 w[0] = first;
1068 for (int i = 0; i < n_vals; i++)
1069 w[i+1] = *(_values[i]);
1070
1071 cs_parall_max(ec, n_vals + 1, datatype, w);
1072
1073 first = w[0];
1074 for (int i = 0; i < n_vals; i++)
1075 *(_values[i]) = w[i+1];
1076 }
1077
1078#endif // defined(HAVE_MPI)
1079}
1080
1081/*----------------------------------------------------------------------------*/
1088/*----------------------------------------------------------------------------*/
1089
1090template <typename T, typename... Vals>
1091static void
1093(
1094 T& first,
1095 Vals&... values
1096)
1097{
1098#if defined(HAVE_MPI)
1099 auto *ec = cs_execution_context_glob_get();
1100 max(ec, first, values...);
1101#else
1102 return;
1103#endif
1104}
1105
1106/*----------------------------------------------------------------------------*/
1114/*----------------------------------------------------------------------------*/
1115
1116template <int Stride, typename T, typename... Vals>
1117static void
1119(
1120 const cs_execution_context *ec,
1121 T first[],
1122 Vals&&... values
1123)
1124{
1125#if defined(HAVE_MPI)
1126
1127 if (!ec->use_mpi())
1128 return;
1129
1130 /* Count number of values */
1131 constexpr int n_vals = sizeof...(Vals);
1132
1133 /* Set datatype for global communication */
1134 cs_datatype_t datatype = cs_datatype_from_type<T>();
1135
1136 /* Temporary work array and parallel sum */
1137 if (n_vals == 0)
1138 cs_parall_max(ec, Stride, datatype, first);
1139 else {
1140 /* Unpack values */
1141 T *_values[] = {values ...};
1142
1143 constexpr int work_size = (n_vals + 1) * Stride;
1144
1145 T w[work_size];
1146 for (int i = 0; i < Stride; i++)
1147 w[i] = first[i];
1148
1149 for (int i = 0; i < n_vals; i++)
1150 for (int j = 0; j < Stride; j++)
1151 w[(i+1)*Stride + j] = _values[i][j];
1152
1153 cs_parall_max(ec, work_size, datatype, w);
1154
1155 for (int i = 0; i < Stride; i++)
1156 first[i] = w[i];
1157
1158 for (int i = 0; i < n_vals; i++)
1159 for (int j = 0; j < Stride; j++)
1160 _values[i][j] = w[(i+1)*Stride + j];
1161 }
1162
1163#endif // defined(HAVE_MPI)
1164}
1165
1166/*----------------------------------------------------------------------------*/
1173/*----------------------------------------------------------------------------*/
1174
1175template <int Stride, typename T, typename... Vals>
1176static void
1178(
1179 T first[],
1180 Vals&&... values
1181)
1182{
1183#if defined(HAVE_MPI)
1184 auto *ec = cs_execution_context_glob_get();
1185 max<Stride>(ec, first, values...);
1186#else
1187 return;
1188#endif
1189}
1190
1191/*----------------------------------------------------------------------------*/
1198/*----------------------------------------------------------------------------*/
1199
1200template <typename T, typename... Vals>
1201static void
1203(
1204 const cs_execution_context *ec,
1205 T& first,
1206 Vals&... values
1207)
1208{
1209#if defined(HAVE_MPI)
1210
1211 if (!ec->use_mpi())
1212 return;
1213
1214 /* Count number of values */
1215 constexpr int n_vals = sizeof...(Vals);
1216
1217 /* Set datatype for global communication */
1218 cs_datatype_t datatype = cs_datatype_from_type<T>();
1219
1220 if (n_vals == 0)
1221 cs_parall_min(ec, 1, datatype, &first);
1222
1223 else {
1224 /* Temporary work array and parallel sum */
1225
1226 /* Unpack values */
1227 T *_values[] = {&values ...};
1228
1229 T w[n_vals + 1];
1230 w[0] = first;
1231 for (int i = 0; i < n_vals; i++)
1232 w[i + 1] = *(_values[i]);
1233
1234 cs_parall_min(ec, n_vals + 1, datatype, w);
1235
1236 first = w[0];
1237 for (int i = 0; i < n_vals; i++)
1238 *(_values[i]) = w[i + 1];
1239 }
1240
1241#endif // defined(HAVE_MPI)
1242}
1243
1244/*----------------------------------------------------------------------------*/
1251/*----------------------------------------------------------------------------*/
1252
1253template <typename T, typename... Vals>
1254static void
1256(
1257 T& first,
1258 Vals&... values
1259)
1260{
1261#if defined(HAVE_MPI)
1262 auto *ec = cs_execution_context_glob_get();
1263 min(ec, first, values...);
1264#else
1265 return;
1266#endif
1267}
1268
1269/*----------------------------------------------------------------------------*/
1277/*----------------------------------------------------------------------------*/
1278
1279template <int Stride, typename T, typename... Vals>
1280static void
1282(
1283 const cs_execution_context *ec,
1284 T first[],
1285 Vals&&... values
1286)
1287{
1288#if defined(HAVE_MPI)
1289
1290 if (!ec->use_mpi())
1291 return;
1292
1293 /* Count number of values */
1294 constexpr int n_vals = sizeof...(Vals);
1295
1296 /* Set datatype for global communication */
1297 cs_datatype_t datatype = cs_datatype_from_type<T>();
1298
1299 /* Temporary work array and parallel sum */
1300 if (n_vals == 0)
1301 cs_parall_min(ec, Stride, datatype, first);
1302 else {
1303 /* Unpack values */
1304 T *_values[] = {values ...};
1305
1306 constexpr int work_size = (n_vals + 1) * Stride;
1307
1308 T w[work_size];
1309 for (int i = 0; i < Stride; i++)
1310 w[i] = first[i];
1311
1312 for (int i = 0; i < n_vals; i++)
1313 for (int j = 0; j < Stride; j++)
1314 w[(i+1)*Stride + j] = _values[i][j];
1315
1316 cs_parall_min(ec, work_size, datatype, w);
1317
1318 for (int i = 0; i < Stride; i++)
1319 first[i] = w[i];
1320
1321 for (int i = 0; i < n_vals; i++)
1322 for (int j = 0; j < Stride; j++)
1323 _values[i][j] = w[(i+1)*Stride + j];
1324 }
1325
1326#endif
1327}
1328
1329/*----------------------------------------------------------------------------*/
1337/*----------------------------------------------------------------------------*/
1338
1339template <int Stride, typename T, typename... Vals>
1340static void
1342(
1343 T first[],
1344 Vals&&... values
1345)
1346{
1347#if defined(HAVE_MPI)
1348 auto *ec = cs_execution_context_glob_get();
1349 min<Stride>(ec, first, values...);
1350#else
1351 return;
1352#endif
1353}
1354
1355} /* namespace parall */
1356} /* namespace cs */
1357
1358/*----------------------------------------------------------------------------*/
1364/*----------------------------------------------------------------------------*/
1365
1366template <typename T, typename... Vals>
1367static void
1369(
1370 T& first,
1371 Vals&... values
1372)
1373{
1374#if defined(HAVE_MPI)
1375
1376 if (cs_glob_n_ranks == 1)
1377 return;
1378
1379 /* Count number of values */
1380 constexpr size_t n_vals = sizeof...(Vals);
1381
1382 /* Set datatype for global communication */
1383 cs_datatype_t datatype = cs_datatype_from_type<T>();
1384
1385 /* Temporary work array and parallel sum */
1386 if (n_vals == 0)
1387 cs_parall_sum(1, datatype, &first);
1388 else {
1389 /* Unpack values */
1390 T *_values[] = {&values ...};
1391
1392 T w[n_vals + 1];
1393 w[0] = first;
1394 for (size_t i = 0; i < n_vals; i++)
1395 w[i+1] = *(_values[i]);
1396
1397 cs_parall_sum(n_vals + 1, datatype, w);
1398
1399 first = w[0];
1400 for (size_t i = 0; i < n_vals; i++)
1401 *(_values[i]) = w[i+1];
1402 }
1403
1404#endif // defined(HAVE_MPI)
1405}
1406
1407/*----------------------------------------------------------------------------*/
1413/*----------------------------------------------------------------------------*/
1414
1415template <typename T, typename... Vals>
1416static void
1418(
1419 const cs_execution_context *ec,
1420 T& first,
1421 Vals&... values
1422)
1423{
1424#if defined(HAVE_MPI)
1425
1426 /* Count number of values */
1427 constexpr size_t n_vals = sizeof...(Vals);
1428
1429 /* Set datatype for global communication */
1430 cs_datatype_t datatype = cs_datatype_from_type<T>();
1431
1432 /* Temporary work array and parallel sum */
1433 if (n_vals == 0)
1434 cs_parall_sum(ec, 1, datatype, &first);
1435 else {
1436 /* Unpack values */
1437 T *_values[] = {&values ...};
1438
1439 T w[n_vals + 1];
1440 w[0] = first;
1441 for (size_t i = 0; i < n_vals; i++)
1442 w[i+1] = *(_values[i]);
1443
1444 cs_parall_sum(ec, n_vals + 1, datatype, w);
1445
1446 first = w[0];
1447 for (size_t i = 0; i < n_vals; i++)
1448 *(_values[i]) = w[i+1];
1449 }
1450
1451#endif // defined(HAVE_MPI)
1452}
1453
1454/*----------------------------------------------------------------------------*/
1461/*----------------------------------------------------------------------------*/
1462
1463template <int Stride, typename T, typename... Vals>
1464static void
1466(
1467 const cs_execution_context *ec,
1468 T first[],
1469 Vals&&... values
1470)
1471{
1472#if defined(HAVE_MPI)
1473 if (!ec->use_mpi())
1474 return;
1475
1476 /* Count number of values */
1477 constexpr size_t n_vals = sizeof...(Vals);
1478
1479 /* Set datatype for global communication */
1480 cs_datatype_t datatype = cs_datatype_from_type<T>();
1481
1482 /* Temporary work array and parallel sum */
1483 if (n_vals == 0)
1484 cs_parall_sum(ec, Stride, datatype, first);
1485 else {
1486 /* Unpack values */
1487 T *_values[] = {values ...};
1488
1489 constexpr size_t work_size = (n_vals + 1) * Stride;
1490
1491 T w[work_size];
1492 for (int i = 0; i < Stride; i++)
1493 w[i] = first[i];
1494
1495 for (size_t i = 0; i < n_vals; i++)
1496 for (int j = 0; j < Stride; j++)
1497 w[(i+1)*Stride + j] = _values[i][j];
1498
1499 cs_parall_sum(ec, work_size, datatype, w);
1500
1501 for (int i = 0; i < Stride; i++)
1502 first[i] = w[i];
1503
1504 for (size_t i = 0; i < n_vals; i++) {
1505 for (int j = 0; j < Stride; j++)
1506 _values[i][j] = w[(i+1)*Stride + j];
1507 }
1508 }
1509
1510#endif // defined(HAVE_MPI)
1511}
1512
1513/*----------------------------------------------------------------------------*/
1521/*----------------------------------------------------------------------------*/
1522
1523template <int Stride, typename T, typename... Vals>
1524static void
1526(
1527 T first[],
1528 Vals&&... values
1529)
1530{
1531#if defined(HAVE_MPI)
1532
1533 if (cs_glob_n_ranks == 1)
1534 return;
1535
1536 /* Count number of values */
1537 constexpr size_t n_vals = sizeof...(Vals);
1538
1539 /* Set datatype for global communication */
1540 cs_datatype_t datatype = cs_datatype_from_type<T>();
1541
1542 /* Temporary work array and parallel sum */
1543 if (n_vals == 0)
1544 cs_parall_sum(Stride, datatype, first);
1545 else {
1546 /* Unpack values */
1547 T *_values[] = {values ...};
1548
1549 constexpr size_t work_size = (n_vals + 1) * Stride;
1550
1551 T w[work_size];
1552 for (int i = 0; i < Stride; i++)
1553 w[i] = first[i];
1554
1555 for (int i = 0; i < n_vals; i++) {
1556 for (int j = 0; j < Stride; j++)
1557 w[(i+1)*Stride + j] = _values[i][j];
1558 }
1559
1560 cs_parall_sum(work_size, datatype, w);
1561
1562 for (int i = 0; i < Stride; i++)
1563 first[i] = w[i];
1564
1565 for (size_t i = 0; i < n_vals; i++) {
1566 for (int j = 0; j < Stride; j++)
1567 _values[i][j] = w[(i+1)*Stride + j];
1568 }
1569 }
1570
1571#endif // defined(HAVE_MPI)
1572}
1573
1574/*----------------------------------------------------------------------------*/
1581/*----------------------------------------------------------------------------*/
1582
1583template <typename T, typename... Vals>
1584static void
1586(
1587 T& first,
1588 Vals&... values
1589)
1590{
1591#if defined(HAVE_MPI)
1592
1593 if (cs_glob_n_ranks == 1)
1594 return;
1595
1596 /* Count number of values */
1597 constexpr size_t n_vals = sizeof...(Vals);
1598
1599 /* Set datatype for global communication */
1600 cs_datatype_t datatype = cs_datatype_from_type<T>();
1601
1602 /* Temporary work array and parallel sum */
1603 if (n_vals == 0)
1604 cs_parall_max(1, datatype, &first);
1605 else {
1606
1607 /* Unpack values */
1608 T *_values[] = {&values ...};
1609
1610 T w[n_vals + 1];
1611 w[0] = first;
1612 for (size_t i = 0; i < n_vals; i++)
1613 w[i+1] = *(_values[i]);
1614
1615 cs_parall_max(n_vals + 1, datatype, w);
1616
1617 first = w[0];
1618 for (size_t i = 0; i < n_vals; i++)
1619 *(_values[i]) = w[i+1];
1620 }
1621
1622#endif // defined(HAVE_MPI)
1623}
1624
1625/*----------------------------------------------------------------------------*/
1631/*----------------------------------------------------------------------------*/
1632
1633template <typename T, typename... Vals>
1634static void
1636(
1637 const cs_execution_context *ec,
1638 T& first,
1639 Vals&... values
1640)
1641{
1642#if defined(HAVE_MPI)
1643
1644 /* Count number of values */
1645 constexpr size_t n_vals = sizeof...(Vals);
1646
1647 /* Set datatype for global communication */
1648 cs_datatype_t datatype = cs_datatype_from_type<T>();
1649
1650 /* Temporary work array and parallel sum */
1651 if (n_vals == 0)
1652 cs_parall_max(ec, 1, datatype, &first);
1653 else {
1654
1655 /* Unpack values */
1656 T *_values[] = {&values ...};
1657
1658 T w[n_vals + 1];
1659 w[0] = first;
1660 for (size_t i = 0; i < n_vals; i++)
1661 w[i+1] = *(_values[i]);
1662
1663 cs_parall_max(ec, n_vals + 1, datatype, w);
1664
1665 first = w[0];
1666 for (size_t i = 0; i < n_vals; i++)
1667 *(_values[i]) = w[i+1];
1668 }
1669
1670#endif // defined(HAVE_MPI)
1671}
1672
1673/*----------------------------------------------------------------------------*/
1681/*----------------------------------------------------------------------------*/
1682
1683template <int Stride, typename T, typename... Vals>
1684static void
1686(
1687 T first[],
1688 Vals&&... values
1689)
1690{
1691#if defined(HAVE_MPI)
1692
1693 if (cs_glob_n_ranks == 1)
1694 return;
1695
1696 /* Count number of values */
1697 constexpr size_t n_vals = sizeof...(Vals);
1698
1699 /* Set datatype for global communication */
1700 cs_datatype_t datatype = cs_datatype_from_type<T>();
1701
1702 /* Temporary work array and parallel sum */
1703 if (n_vals == 0)
1704 cs_parall_max(Stride, datatype, first);
1705 else {
1706 /* Unpack values */
1707 T *_values[] = {values ...};
1708
1709 constexpr size_t work_size = (n_vals + 1) * Stride;
1710
1711 T w[work_size];
1712 for (int i = 0; i < Stride; i++)
1713 w[i] = first[i];
1714
1715 for (size_t i = 0; i < n_vals; i++)
1716 for (int j = 0; j < Stride; j++)
1717 w[(i+1)*Stride + j] = _values[i][j];
1718
1719 cs_parall_max(work_size, datatype, w);
1720
1721 for (int i = 0; i < Stride; i++)
1722 first[i] = w[i];
1723
1724 for (size_t i = 0; i < n_vals; i++)
1725 for (int j = 0; j < Stride; j++)
1726 _values[i][j] = w[(i+1)*Stride + j];
1727 }
1728
1729#endif // defined(HAVE_MPI)
1730}
1731
1732/*----------------------------------------------------------------------------*/
1739/*----------------------------------------------------------------------------*/
1740
1741template <int Stride, typename T, typename... Vals>
1742static void
1744(
1745 const cs_execution_context *ec,
1746 T first[],
1747 Vals&&... values
1748)
1749{
1750#if defined(HAVE_MPI)
1751
1752 /* Count number of values */
1753 constexpr size_t n_vals = sizeof...(Vals);
1754
1755 /* Set datatype for global communication */
1756 cs_datatype_t datatype = cs_datatype_from_type<T>();
1757
1758 /* Temporary work array and parallel sum */
1759 if (n_vals == 0)
1760 cs_parall_max(ec, Stride, datatype, first);
1761 else {
1762 /* Unpack values */
1763 T *_values[] = {values ...};
1764
1765 constexpr size_t work_size = (n_vals + 1) * Stride;
1766
1767 T w[work_size];
1768 for (int i = 0; i < Stride; i++)
1769 w[i] = first[i];
1770
1771 for (size_t i = 0; i < n_vals; i++)
1772 for (int j = 0; j < Stride; j++)
1773 w[(i+1)*Stride + j] = _values[i][j];
1774
1775 cs_parall_max(ec, work_size, datatype, w);
1776
1777 for (int i = 0; i < Stride; i++)
1778 first[i] = w[i];
1779
1780 for (size_t i = 0; i < n_vals; i++)
1781 for (int j = 0; j < Stride; j++)
1782 _values[i][j] = w[(i+1)*Stride + j];
1783 }
1784
1785#endif // defined(HAVE_MPI)
1786}
1787
1788/*----------------------------------------------------------------------------*/
1795/*----------------------------------------------------------------------------*/
1796
1797template <typename T, typename... Vals>
1798static void
1800(
1801 T& first,
1802 Vals&... values
1803)
1804{
1805#if defined(HAVE_MPI)
1806
1807 if (cs_glob_n_ranks == 1)
1808 return;
1809
1810 /* Count number of values */
1811 constexpr size_t n_vals = sizeof...(Vals);
1812
1813 /* Set datatype for global communication */
1814 cs_datatype_t datatype = cs_datatype_from_type<T>();
1815
1816 if (n_vals == 0)
1817 cs_parall_min(1, datatype, &first);
1818
1819 else {
1820 /* Temporary work array and parallel sum */
1821
1822 /* Unpack values */
1823 T *_values[] = {&values ...};
1824
1825 T w[n_vals + 1];
1826 w[0] = first;
1827 for (size_t i = 0; i < n_vals; i++)
1828 w[i + 1] = *(_values[i]);
1829
1830 cs_parall_min(n_vals + 1, datatype, w);
1831
1832 first = w[0];
1833 for (size_t i = 0; i < n_vals; i++)
1834 *(_values[i]) = w[i + 1];
1835 }
1836
1837#endif // defined(HAVE_MPI)
1838}
1839
1840/*----------------------------------------------------------------------------*/
1846/*----------------------------------------------------------------------------*/
1847
1848template <typename T, typename... Vals>
1849static void
1851(
1852 const cs_execution_context *ec,
1853 T& first,
1854 Vals&... values
1855)
1856{
1857#if defined(HAVE_MPI)
1858
1859 /* Count number of values */
1860 constexpr size_t n_vals = sizeof...(Vals);
1861
1862 /* Set datatype for global communication */
1863 cs_datatype_t datatype = cs_datatype_from_type<T>();
1864
1865 if (n_vals == 0)
1866 cs_parall_min(ec, 1, datatype, &first);
1867
1868 else {
1869 /* Temporary work array and parallel sum */
1870
1871 /* Unpack values */
1872 T *_values[] = {&values ...};
1873
1874 T w[n_vals + 1];
1875 w[0] = first;
1876 for (size_t i = 0; i < n_vals; i++)
1877 w[i + 1] = *(_values[i]);
1878
1879 cs_parall_min(ec, n_vals + 1, datatype, w);
1880
1881 first = w[0];
1882 for (size_t i = 0; i < n_vals; i++)
1883 *(_values[i]) = w[i + 1];
1884 }
1885
1886#endif // defined(HAVE_MPI)
1887}
1888
1889/*----------------------------------------------------------------------------*/
1897/*----------------------------------------------------------------------------*/
1898
1899template <int Stride, typename T, typename... Vals>
1900static void
1902(
1903 T first[],
1904 Vals&&... values
1905)
1906{
1907#if defined(HAVE_MPI)
1908
1909 if (cs_glob_n_ranks == 1)
1910 return;
1911
1912 /* Count number of values */
1913 constexpr size_t n_vals = sizeof...(Vals);
1914
1915 /* Set datatype for global communication */
1916 cs_datatype_t datatype = cs_datatype_from_type<T>();
1917
1918 /* Temporary work array and parallel sum */
1919 if (n_vals == 0)
1920 cs_parall_min(Stride, datatype, first);
1921 else {
1922 /* Unpack values */
1923 T *_values[] = {values ...};
1924
1925 constexpr size_t work_size = (n_vals + 1) * Stride;
1926
1927 T w[work_size];
1928 for (int i = 0; i < Stride; i++)
1929 w[i] = first[i];
1930
1931 for (size_t i = 0; i < n_vals; i++)
1932 for (int j = 0; j < Stride; j++)
1933 w[(i+1)*Stride + j] = _values[i][j];
1934
1935 cs_parall_min(work_size, datatype, w);
1936
1937 for (int i = 0; i < Stride; i++)
1938 first[i] = w[i];
1939
1940 for (size_t i = 0; i < n_vals; i++)
1941 for (int j = 0; j < Stride; j++)
1942 _values[i][j] = w[(i+1)*Stride + j];
1943 }
1944
1945#endif
1946}
1947
1948/*----------------------------------------------------------------------------*/
1955/*----------------------------------------------------------------------------*/
1956
1957template <int Stride, typename T, typename... Vals>
1958static void
1960(
1961 const cs_execution_context *ec,
1962 T first[],
1963 Vals&&... values
1964)
1965{
1966#if defined(HAVE_MPI)
1967
1968 /* Count number of values */
1969 constexpr size_t n_vals = sizeof...(Vals);
1970
1971 /* Set datatype for global communication */
1972 cs_datatype_t datatype = cs_datatype_from_type<T>();
1973
1974 /* Temporary work array and parallel sum */
1975 if (n_vals == 0)
1976 cs_parall_min(ec, Stride, datatype, first);
1977
1978 else {
1979 /* Unpack values */
1980 T *_values[] = {values ...};
1981
1982 constexpr size_t work_size = (n_vals + 1) * Stride;
1983
1984 T w[work_size];
1985 for (int i = 0; i < Stride; i++)
1986 w[i] = first[i];
1987
1988 for (size_t i = 0; i < n_vals; i++)
1989 for (int j = 0; j < Stride; j++)
1990 w[(i+1)*Stride + j] = _values[i][j];
1991
1992 cs_parall_min(ec, work_size, datatype, w);
1993
1994 for (int i = 0; i < Stride; i++)
1995 first[i] = w[i];
1996
1997 for (size_t i = 0; i < n_vals; i++)
1998 for (int j = 0; j < Stride; j++)
1999 _values[i][j] = w[(i+1)*Stride + j];
2000 }
2001#endif // defined(HAVE_MPI)
2002}
2003
2004#endif //__cplusplus
2005
2006/*----------------------------------------------------------------------------*/
2007
2008#endif /* __CS_PARALL_H__ */
Definition: cs_execution_context.h:61
bool use_mpi() const
Does the execution context uses MPI parallelism ?
Definition: cs_execution_context.h:128
int cs_glob_n_ranks
Definition: cs_defs.cpp:175
cs_datatype_t
Definition: cs_defs.h:315
#define BEGIN_C_DECLS
Definition: cs_defs.h:554
double cs_real_t
Floating-point value.
Definition: cs_defs.h:357
unsigned cs_gnum_t
global mesh entity number
Definition: cs_defs.h:342
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:669
#define CS_UNUSED(x)
Definition: cs_defs.h:543
#define END_C_DECLS
Definition: cs_defs.h:555
int cs_lnum_t
local mesh entity id
Definition: cs_defs.h:350
#define CS_CL_SIZE
Definition: cs_defs.h:513
const cs_execution_context * cs_execution_context_glob_get(void)
Get the global execution context.
Definition: cs_execution_context.cpp:61
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.cpp:531
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:241
static void cs_parall_min_strided(T first[], Vals &&... values)
Minimum values of a given datatype on all default communicator processes.
Definition: cs_parall.h:1902
void cs_parall_set_min_coll_buf_size(size_t buffer_size)
Define minimum recommended scatter or gather buffer size.
Definition: cs_parall.cpp:854
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.cpp:598
static void cs_parall_sum_strided(const cs_execution_context *ec, T first[], Vals &&... values)
Sum strided-values of a given datatype over a communicator.
Definition: cs_parall.h:1466
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.cpp:356
static void cs_parall_thread_range_upper(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 for upper triangular matrix elements.
Definition: cs_parall.h:632
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:175
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:114
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.cpp:410
static int cs_parall_n_threads(cs_lnum_t n_elements, cs_lnum_t min_thread_elements)
Compute recommended number of threads for a section.
Definition: cs_parall.h:553
static void cs_parall_max_scalars(T &first, Vals &... values)
Maximum values of a given datatype on all default communicator processes.
Definition: cs_parall.h:1586
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:86
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.cpp:647
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:143
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.cpp:487
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.cpp:777
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:587
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.cpp:317
static void cs_parall_sum_scalars(T &first, Vals &... values)
Sum values of a given datatype on all default communicator processes.
Definition: cs_parall.h:1369
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:207
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.cpp:712
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.cpp:279
size_t cs_parall_get_min_coll_buf_size(void)
Return minimum recommended scatter or gather buffer size.
Definition: cs_parall.cpp:832
static size_t cs_parall_block_count(size_t n, size_t block_size)
Compute number of blocks needed for a given array and block sizes.
Definition: cs_parall.h:675
static void cs_parall_min_scalars(T &first, Vals &... values)
Minimum values of a given datatype on all default communicator processes.
Definition: cs_parall.h:1800
static void cs_parall_max_strided(T first[], Vals &&... values)
Maximum values of a given datatype on all default communicator processes.
Definition: cs_parall.h:1686
cs_e2n_sum_t
Definition: cs_parall.h:52
@ CS_E2N_SUM_SCATTER_ATOMIC
Definition: cs_parall.h:57
@ CS_E2N_SUM_SCATTER
Definition: cs_parall.h:54
@ CS_E2N_SUM_GATHER
Definition: cs_parall.h:59
cs_e2n_sum_t cs_glob_e2n_sum_type
static void sum(const cs_execution_context *ec, T &first, Vals &... values)
Sum values of a given datatype over a given communicator.
Definition: cs_parall.h:881
static void max(const cs_execution_context *ec, T &first, Vals &... values)
Maximum values of a given datatype on a given communicator processes.
Definition: cs_parall.h:1044
static void max(T first[], Vals &&... values)
Maximum values of a given datatype on all default communicator processes.
Definition: cs_parall.h:1178
static void min(const cs_execution_context *ec, T &first, Vals &... values)
Minimum values of a given datatype on a given communicator processes.
Definition: cs_parall.h:1203
static void min(T first[], Vals &&... values)
Minimum values of a given datatype on all default communicator processes.
Definition: cs_parall.h:1342
Definition: cs_array.h:1098
Definition: parall.f90:26