DUNE PDELab (2.8)

mpitraits.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_MPITRAITS_HH
4#define DUNE_MPITRAITS_HH
5
16#if HAVE_MPI
17
18#include <cstddef>
19#include <cstdint>
20#include <type_traits>
21#include <utility>
22#include <complex>
23
24#include <mpi.h>
25
26namespace Dune
27{
37 template<typename T>
38 struct MPITraits
39 {
40 private:
41 MPITraits(){}
42 MPITraits(const MPITraits&){}
43 static MPI_Datatype datatype;
44 static MPI_Datatype vectortype;
45 public:
46 static inline MPI_Datatype getType()
47 {
48 if(datatype==MPI_DATATYPE_NULL) {
49 MPI_Type_contiguous(sizeof(T),MPI_BYTE,&datatype);
50 MPI_Type_commit(&datatype);
51 }
52 return datatype;
53 }
54 static constexpr bool is_intrinsic = false;
55 };
56 template<class T>
57 MPI_Datatype MPITraits<T>::datatype = MPI_DATATYPE_NULL;
58
59#ifndef DOXYGEN
60
61 // A Macro for defining traits for the primitive data types
62#define ComposeMPITraits(p,m) \
63 template<> \
64 struct MPITraits<p>{ \
65 static inline MPI_Datatype getType(){ \
66 return m; \
67 } \
68 static constexpr bool is_intrinsic = true; \
69 }
70
71 ComposeMPITraits(char, MPI_CHAR);
72 ComposeMPITraits(unsigned char,MPI_UNSIGNED_CHAR);
73 ComposeMPITraits(short,MPI_SHORT);
74 ComposeMPITraits(unsigned short,MPI_UNSIGNED_SHORT);
75 ComposeMPITraits(int,MPI_INT);
76 ComposeMPITraits(unsigned int,MPI_UNSIGNED);
77 ComposeMPITraits(long,MPI_LONG);
78 ComposeMPITraits(unsigned long,MPI_UNSIGNED_LONG);
79 ComposeMPITraits(float,MPI_FLOAT);
80 ComposeMPITraits(double,MPI_DOUBLE);
81 ComposeMPITraits(long double,MPI_LONG_DOUBLE);
82 ComposeMPITraits(std::complex<double>, MPI_CXX_DOUBLE_COMPLEX);
83 ComposeMPITraits(std::complex<long double>, MPI_CXX_LONG_DOUBLE_COMPLEX);
84 ComposeMPITraits(std::complex<float>, MPI_CXX_FLOAT_COMPLEX);
85
86
87#undef ComposeMPITraits
88
89 template<class K, int n> class FieldVector;
90
91 template<class K, int n>
92 struct MPITraits<FieldVector<K,n> >
93 {
94 static MPI_Datatype datatype;
95 static MPI_Datatype vectortype;
96
97 static inline MPI_Datatype getType()
98 {
99 if(datatype==MPI_DATATYPE_NULL) {
100 MPI_Type_contiguous(n, MPITraits<K>::getType(), &vectortype);
101 MPI_Type_commit(&vectortype);
102 FieldVector<K,n> fvector;
103 MPI_Aint base;
104 MPI_Aint displ;
105 MPI_Get_address(&fvector, &base);
106 MPI_Get_address(&(fvector[0]), &displ);
107 displ -= base;
108 int length[1]={1};
109
110 MPI_Type_create_struct(1, length, &displ, &vectortype, &datatype);
111 MPI_Type_commit(&datatype);
112 }
113 return datatype;
114 }
115
116 };
117
118 template<class K, int n>
119 MPI_Datatype MPITraits<FieldVector<K,n> >::datatype = MPI_DATATYPE_NULL;
120 template<class K, int n>
121 MPI_Datatype MPITraits<FieldVector<K,n> >::vectortype = {MPI_DATATYPE_NULL};
122
123
124 template<int k>
125 class bigunsignedint;
126
127 template<int k>
128 struct MPITraits<bigunsignedint<k> >
129 {
130 static MPI_Datatype datatype;
131 static MPI_Datatype vectortype;
132
133 static inline MPI_Datatype getType()
134 {
135 if(datatype==MPI_DATATYPE_NULL) {
136 MPI_Type_contiguous(bigunsignedint<k>::n, MPITraits<std::uint16_t>::getType(),
137 &vectortype);
138 //MPI_Type_commit(&vectortype);
139 bigunsignedint<k> data;
140 MPI_Aint base;
141 MPI_Aint displ;
142 MPI_Get_address(&data, &base);
143 MPI_Get_address(&(data.digit), &displ);
144 displ -= base;
145 int length[1]={1};
146 MPI_Type_create_struct(1, length, &displ, &vectortype, &datatype);
147 MPI_Type_commit(&datatype);
148 }
149 return datatype;
150 }
151 };
152}
153
154namespace Dune
155{
156 template<int k>
157 MPI_Datatype MPITraits<bigunsignedint<k> >::datatype = MPI_DATATYPE_NULL;
158 template<int k>
159 MPI_Datatype MPITraits<bigunsignedint<k> >::vectortype = MPI_DATATYPE_NULL;
160
161 template<typename T1, typename T2>
162 struct MPITraits<std::pair<T1,T2 > >
163 {
164 public:
165 inline static MPI_Datatype getType();
166 private:
167 static MPI_Datatype type;
168 };
169 template<typename T1, typename T2>
170 MPI_Datatype MPITraits<std::pair<T1,T2> >::getType()
171 {
172 if(type==MPI_DATATYPE_NULL) {
173 int length[2] = {1, 1};
174 MPI_Aint disp[2];
175 MPI_Datatype types[2] = {MPITraits<T1>::getType(),
176 MPITraits<T2>::getType()};
177
178 using Pair = std::pair<T1, T2>;
179 static_assert(std::is_standard_layout<Pair>::value, "offsetof() is only defined for standard layout types");
180 disp[0] = offsetof(Pair, first);
181 disp[1] = offsetof(Pair, second);
182
183 MPI_Datatype tmp;
184 MPI_Type_create_struct(2, length, disp, types, &tmp);
185
186 MPI_Type_create_resized(tmp, 0, sizeof(Pair), &type);
187 MPI_Type_commit(&type);
188
189 MPI_Type_free(&tmp);
190 }
191 return type;
192 }
193
194 template<typename T1, typename T2>
195 MPI_Datatype MPITraits<std::pair<T1,T2> >::type=MPI_DATATYPE_NULL;
196
197#endif // !DOXYGEN
198
199} // namespace Dune
200
201#endif // HAVE_MPI
202
205#endif
vector space out of a tensor product of fields.
Definition: fvector.hh:95
Dune namespace.
Definition: alignedallocator.hh:11
STL namespace.
A traits class describing the mapping of types onto MPI_Datatypes.
Definition: mpitraits.hh:39
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 24, 22:29, 2024)