UPstreamReduce.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "Pstream.H"
29#include "PstreamReduceOps.H"
30#include "UPstreamWrapping.H"
31
32#include <mpi.h>
33#include <cinttypes>
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37// Specialisations for bool
38
39void Foam::reduce
40(
41 bool& value,
42 const andOp<bool>&,
43 const int tag, /* (unused) */
44 const label comm
45)
46{
47 // This can also work:
48 // PstreamDetail::allReduce(&value, 1, MPI_BYTE, MPI_BAND, comm);
49 PstreamDetail::allReduce(&value, 1, MPI_C_BOOL, MPI_LAND, comm);
50}
51
52
53void Foam::reduce
54(
55 bool& value,
56 const orOp<bool>&,
57 const int tag, /* (unused) */
58 const label comm
59)
60{
61 // This can also work:
62 // PstreamDetail::allReduce(&value, 1, MPI_BYTE, MPI_BOR, comm);
63 PstreamDetail::allReduce(&value, 1, MPI_C_BOOL, MPI_LOR, comm);
64}
65
66
67// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68
69// Specialisations for common reduction types
70
71#undef Pstream_CommonReductions
72#define Pstream_CommonReductions(Native, TaggedType) \
73 \
74void Foam::reduce \
75( \
76 Native& value, \
77 const minOp<Native>&, \
78 const int tag, /* (unused) */ \
79 const label comm \
80) \
81{ \
82 PstreamDetail::allReduce<Native> \
83 ( \
84 &value, 1, TaggedType, MPI_MIN, comm \
85 ); \
86} \
87 \
88void Foam::reduce \
89( \
90 Native& value, \
91 const maxOp<Native>&, \
92 const int tag, /* (unused) */ \
93 const label comm \
94) \
95{ \
96 PstreamDetail::allReduce<Native> \
97 ( \
98 &value, 1, TaggedType, MPI_MAX, comm \
99 ); \
100} \
101 \
102void Foam::reduce \
103( \
104 Native& value, \
105 const sumOp<Native>&, \
106 const int tag, /* (unused) */ \
107 const label comm \
108) \
109{ \
110 PstreamDetail::allReduce<Native> \
111 ( \
112 &value, 1, TaggedType, MPI_SUM, comm \
113 ); \
114} \
115 \
116void Foam::reduce \
117( \
118 Native values[], \
119 const int size, \
120 const sumOp<Native>&, \
121 const int tag, /* (unused) */ \
122 const label comm \
123) \
124{ \
125 PstreamDetail::allReduce<Native> \
126 ( \
127 values, size, TaggedType, MPI_SUM, comm \
128 ); \
129} \
130
131
132Pstream_CommonReductions(int32_t, MPI_INT32_T);
133Pstream_CommonReductions(int64_t, MPI_INT64_T);
134Pstream_CommonReductions(uint32_t, MPI_UINT32_T);
135Pstream_CommonReductions(uint64_t, MPI_UINT64_T);
136Pstream_CommonReductions(float, MPI_FLOAT);
137Pstream_CommonReductions(double, MPI_DOUBLE);
138
139#undef Pstream_CommonReductions
140
141
142// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
143
144// Specialisations for floating-point types
145
146#undef Pstream_FloatReductions
147#define Pstream_FloatReductions(Native, TaggedType) \
148 \
149void Foam::sumReduce \
150( \
151 Native& value, \
152 label& count, \
153 const int tag, /* (unused) */ \
154 const label comm \
155) \
156{ \
157 if (UPstream::parRun()) \
158 { \
159 Native values[2]; \
160 values[0] = value; \
161 values[1] = static_cast<Native>(count); \
162 \
163 PstreamDetail::allReduce<Native> \
164 ( \
165 values, 2, TaggedType, MPI_SUM, comm \
166 ); \
167 \
168 value = values[0]; \
169 count = static_cast<label>(values[1]); \
170 } \
171} \
172 \
173void Foam::reduce \
174( \
175 Native& value, \
176 const sumOp<Native>&, \
177 const int tag, /* (unused) */ \
178 const label comm, \
179 label& requestID \
180) \
181{ \
182 PstreamDetail::allReduce<Native> \
183 ( \
184 &value, 1, TaggedType, MPI_SUM, comm, &requestID \
185 ); \
186} \
187 \
188void Foam::reduce \
189( \
190 Native values[], \
191 const int size, \
192 const sumOp<Native>&, \
193 const int tag, /* (unused) */ \
194 const label comm, \
195 label& requestID \
196) \
197{ \
198 PstreamDetail::allReduce<Native> \
199 ( \
200 values, size, TaggedType, MPI_SUM, comm, &requestID \
201 ); \
202}
203
204
205Pstream_FloatReductions(float, MPI_FLOAT);
206Pstream_FloatReductions(double, MPI_DOUBLE);
207
208#undef Pstream_FloatReductions
209
210
211// ************************************************************************* //
Inter-processor communication reduction functions.
#define Pstream_FloatReductions(Native, TaggedType)
#define Pstream_CommonReductions(Native, TaggedType)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)