UPstreamWrapping.H
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) 2012-2016 OpenFOAM Foundation
9  Copyright (C) 2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Namespace
28  Foam::PstreamDetail
29 
30 Description
31  Some implementation details for Pstream and/or MPI.
32 
33 InNamespace
34  Foam::PstreamDetail
35 
36 Description
37  Functions to wrap MPI_Bcast, MPI_Allreduce, MPI_Iallreduce etc.
38 
39 SourceFiles
40  UPstreamWrappingTemplates.C
41 
42 \*---------------------------------------------------------------------------*/
43 
44 #ifndef Foam_UPstreamWrapping_H
45 #define Foam_UPstreamWrapping_H
46 
47 #include "UPstream.H"
48 #include <mpi.h>
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 namespace Foam
53 {
54 namespace PstreamDetail
55 {
56 
57 // MPI_Bcast, using root=0
58 template<class Type>
59 void broadcast0
60 (
61  Type* values,
62  int count,
63  MPI_Datatype datatype,
64  const label comm
65 );
66 
67 // MPI_Reduce, using root=0
68 template<class Type>
69 void reduce0
70 (
71  Type* values,
72  int count,
73  MPI_Datatype datatype,
74  MPI_Op optype,
75  const label comm
76 );
77 
78 // MPI_Allreduce or MPI_Iallreduce
79 template<class Type>
80 void allReduce
81 (
82  Type* values,
83  int count,
84  MPI_Datatype datatype,
85  MPI_Op optype,
86  const label comm, // Communicator
87  label* requestID = nullptr // Non-null for MPI_Iallreduce
88 );
89 
90 
91 // MPI_Alltoall or MPI_Ialltoall with one element per rank
92 template<class Type>
93 void allToAll
94 (
95  const UList<Type>& sendData,
96  UList<Type>& recvData,
97  MPI_Datatype datatype,
98  const label comm, // Communicator
99  label* requestID = nullptr // Non-null for MPI_Ialltoall
100 );
101 
102 
103 // MPI_Alltoallv or MPI_Ialltoallv
104 template<class Type>
105 void allToAllv
106 (
107  const Type* sendData,
108  const UList<int>& sendCounts,
109  const UList<int>& sendOffsets,
110 
111  Type* recvData,
112  const UList<int>& recvCounts,
113  const UList<int>& recvOffsets,
114 
115  MPI_Datatype datatype,
116  const label comm, // Communicator
117  label* requestID = nullptr // Non-null for MPI_Ialltoallv
118 );
119 
120 
121 // MPI_Gather or MPI_Igather
122 template<class Type>
123 void gather
124 (
125  const Type* sendData,
126  int sendCount,
127 
128  Type* recvData, // Ignored on non-root rank
129  int recvCount, // Ignored on non-root rank
130 
131  MPI_Datatype datatype, // The send/recv data type
132  const label comm, // Communicator
133  label* requestID = nullptr // Non-null for MPI_Igather
134 );
135 
136 
137 // MPI_Scatter or MPI_Iscatter
138 template<class Type>
139 void scatter
140 (
141  const Type* sendData, // Ignored on non-root rank
142  int sendCount, // Ignored on non-root rank
143 
144  Type* recvData,
145  int recvCount,
146 
147  MPI_Datatype datatype, // The send/recv data type
148  const label comm, // Communicator
149  label* requestID = nullptr // Non-null for MPI_Iscatter
150 );
151 
152 
153 // MPI_Gatherv or MPI_Igatherv
154 template<class Type>
155 void gatherv
156 (
157  const Type* sendData,
158  int sendCount, // Ignored on master if recvCounts[0] == 0
159 
160  Type* recvData, // Ignored on non-root rank
161  const UList<int>& recvCounts, // Ignored on non-root rank
162  const UList<int>& recvOffsets, // Ignored on non-root rank
163 
164  MPI_Datatype datatype, // The send/recv data type
165  const label comm, // Communicator
166  label* requestID = nullptr // Non-null for MPI_Igatherv
167 );
168 
169 
170 // MPI_Scatterv or MPI_Iscatterv
171 template<class Type>
172 void scatterv
173 (
174  const Type* sendData, // Ignored on non-root rank
175  const UList<int>& sendCounts, // Ignored on non-root rank
176  const UList<int>& sendOffsets, // Ignored on non-root rank
177 
178  Type* recvData,
179  int recvCount,
180 
181  MPI_Datatype datatype, // The send/recv data type
182  const label comm, // Communicator
183  label* requestID = nullptr // Non-null for MPI_Igatherv
184 );
185 
186 
187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 
189 } // End namespace PstreamDetail
190 } // End namespace Foam
191 
192 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
193 
194 #ifdef NoRepository
195  #include "UPstreamWrappingTemplates.C"
196 #endif
197 
198 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
199 
200 #endif
201 
202 // ************************************************************************* //
void scatterv(const Type *sendData, const UList< int > &sendCounts, const UList< int > &sendOffsets, Type *recvData, int recvCount, MPI_Datatype datatype, const label comm, label *requestID=nullptr)
void gatherv(const Type *sendData, int sendCount, Type *recvData, const UList< int > &recvCounts, const UList< int > &recvOffsets, MPI_Datatype datatype, const label comm, label *requestID=nullptr)
void allReduce(Type *values, int count, MPI_Datatype datatype, MPI_Op optype, const label comm, label *requestID=nullptr)
void broadcast0(Type *values, int count, MPI_Datatype datatype, const label comm)
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:164
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of &#39;true&#39; entries.
Definition: BitOps.H:73
void gather(const Type *sendData, int sendCount, Type *recvData, int recvCount, MPI_Datatype datatype, const label comm, label *requestID=nullptr)
void allToAll(const UList< Type > &sendData, UList< Type > &recvData, MPI_Datatype datatype, const label comm, label *requestID=nullptr)
void allToAllv(const Type *sendData, const UList< int > &sendCounts, const UList< int > &sendOffsets, Type *recvData, const UList< int > &recvCounts, const UList< int > &recvOffsets, MPI_Datatype datatype, const label comm, label *requestID=nullptr)
void reduce0(Type *values, int count, MPI_Datatype datatype, MPI_Op optype, const label comm)
void scatter(const Type *sendData, int sendCount, Type *recvData, int recvCount, MPI_Datatype datatype, const label comm, label *requestID=nullptr)
Namespace for OpenFOAM.