DUNE PDELab (2.8)

loop.hh
1#ifndef DUNE_COMMON_SIMD_LOOP_HH
2#define DUNE_COMMON_SIMD_LOOP_HH
3
4#include <array>
5#include <cmath>
6#include <cstddef>
7#include <cstdlib>
8#include <ostream>
9
10#include <dune/common/math.hh>
13
14namespace Dune {
15
16/*
17 * silence warnings from GCC about using integer operands on a bool
18 * (when instantiated for T=bool)
19 */
20#if __GNUC__ >= 7
21# pragma GCC diagnostic push
22# pragma GCC diagnostic ignored "-Wbool-operation"
23# pragma GCC diagnostic ignored "-Wint-in-bool-context"
24#endif
25
37 template<class T, std::size_t S, std::size_t A = 0>
38 class alignas(A==0?alignof(T):A) LoopSIMD : public std::array<T,S> {
39
40 public:
41
42 //default constructor
43 LoopSIMD() {
44 assert(reinterpret_cast<uintptr_t>(this) % std::min(alignof(LoopSIMD<T,S,A>),alignof(std::max_align_t)) == 0);
45 }
46
47 // broadcast constructor initializing the content with a given value
49 this->fill(i);
50 }
51
52 template<std::size_t OA>
53 explicit LoopSIMD(const LoopSIMD<T,S,OA>& other)
54 : std::array<T,S>(other)
55 {
56 assert(reinterpret_cast<uintptr_t>(this) % std::min(alignof(LoopSIMD<T,S,A>),alignof(std::max_align_t)) == 0);
57 }
58
59 /*
60 * Definition of basic operators
61 */
62
63 //Prefix operators
64#define DUNE_SIMD_LOOP_PREFIX_OP(SYMBOL) \
65 auto operator SYMBOL() { \
66 for(std::size_t i=0; i<S; i++){ \
67 SYMBOL(*this)[i]; \
68 } \
69 return *this; \
70 } \
71 static_assert(true, "expecting ;")
72
73 DUNE_SIMD_LOOP_PREFIX_OP(++);
74 DUNE_SIMD_LOOP_PREFIX_OP(--);
75#undef DUNE_SIMD_LOOP_PREFIX_OP
76
77 //Unary operators
78#define DUNE_SIMD_LOOP_UNARY_OP(SYMBOL) \
79 auto operator SYMBOL() const { \
80 LoopSIMD<T,S,A> out; \
81 for(std::size_t i=0; i<S; i++){ \
82 out[i] = SYMBOL((*this)[i]); \
83 } \
84 return out; \
85 } \
86 static_assert(true, "expecting ;")
87
88 DUNE_SIMD_LOOP_UNARY_OP(+);
89 DUNE_SIMD_LOOP_UNARY_OP(-);
90 DUNE_SIMD_LOOP_UNARY_OP(~);
91
92 auto operator!() const {
94 for(std::size_t i=0; i<S; i++){
95 out[i] = !((*this)[i]);
96 }
97 return out;
98 }
99#undef DUNE_SIMD_LOOP_UNARY_OP
100
101 //Postfix operators
102#define DUNE_SIMD_LOOP_POSTFIX_OP(SYMBOL) \
103 auto operator SYMBOL(int){ \
104 LoopSIMD<T,S,A> out = *this; \
105 SYMBOL(*this); \
106 return out; \
107 } \
108 static_assert(true, "expecting ;")
109
110 DUNE_SIMD_LOOP_POSTFIX_OP(++);
111 DUNE_SIMD_LOOP_POSTFIX_OP(--);
112#undef DUNE_SIMD_LOOP_POSTFIX_OP
113
114 //Assignment operators
115#define DUNE_SIMD_LOOP_ASSIGNMENT_OP(SYMBOL) \
116 auto operator SYMBOL(const Simd::Scalar<T> s) { \
117 for(std::size_t i=0; i<S; i++){ \
118 (*this)[i] SYMBOL s; \
119 } \
120 return *this; \
121 } \
122 \
123 auto operator SYMBOL(const LoopSIMD<T,S,A> &v) { \
124 for(std::size_t i=0; i<S; i++){ \
125 (*this)[i] SYMBOL v[i]; \
126 } \
127 return *this; \
128 } \
129 static_assert(true, "expecting ;")
130
131 DUNE_SIMD_LOOP_ASSIGNMENT_OP(+=);
132 DUNE_SIMD_LOOP_ASSIGNMENT_OP(-=);
133 DUNE_SIMD_LOOP_ASSIGNMENT_OP(*=);
134 DUNE_SIMD_LOOP_ASSIGNMENT_OP(/=);
135 DUNE_SIMD_LOOP_ASSIGNMENT_OP(%=);
136 DUNE_SIMD_LOOP_ASSIGNMENT_OP(<<=);
137 DUNE_SIMD_LOOP_ASSIGNMENT_OP(>>=);
138 DUNE_SIMD_LOOP_ASSIGNMENT_OP(&=);
139 DUNE_SIMD_LOOP_ASSIGNMENT_OP(|=);
140 DUNE_SIMD_LOOP_ASSIGNMENT_OP(^=);
141#undef DUNE_SIMD_LOOP_ASSIGNMENT_OP
142 };
143
144 //Arithmetic operators
145#define DUNE_SIMD_LOOP_BINARY_OP(SYMBOL) \
146 template<class T, std::size_t S, std::size_t A> \
147 auto operator SYMBOL(const LoopSIMD<T,S,A> &v, const Simd::Scalar<T> s) { \
148 LoopSIMD<T,S,A> out; \
149 for(std::size_t i=0; i<S; i++){ \
150 out[i] = v[i] SYMBOL s; \
151 } \
152 return out; \
153 } \
154 template<class T, std::size_t S, std::size_t A> \
155 auto operator SYMBOL(const Simd::Scalar<T> s, const LoopSIMD<T,S,A> &v) { \
156 LoopSIMD<T,S,A> out; \
157 for(std::size_t i=0; i<S; i++){ \
158 out[i] = s SYMBOL v[i]; \
159 } \
160 return out; \
161 } \
162 template<class T, std::size_t S, std::size_t A> \
163 auto operator SYMBOL(const LoopSIMD<T,S,A> &v, \
164 const LoopSIMD<T,S,A> &w) { \
165 LoopSIMD<T,S,A> out; \
166 for(std::size_t i=0; i<S; i++){ \
167 out[i] = v[i] SYMBOL w[i]; \
168 } \
169 return out; \
170 } \
171 static_assert(true, "expecting ;")
172
173 DUNE_SIMD_LOOP_BINARY_OP(+);
174 DUNE_SIMD_LOOP_BINARY_OP(-);
175 DUNE_SIMD_LOOP_BINARY_OP(*);
176 DUNE_SIMD_LOOP_BINARY_OP(/);
177 DUNE_SIMD_LOOP_BINARY_OP(%);
178
179 DUNE_SIMD_LOOP_BINARY_OP(&);
180 DUNE_SIMD_LOOP_BINARY_OP(|);
181 DUNE_SIMD_LOOP_BINARY_OP(^);
182
183#undef DUNE_SIMD_LOOP_BINARY_OP
184
185 //Bitshift operators
186#define DUNE_SIMD_LOOP_BITSHIFT_OP(SYMBOL) \
187 template<class T, std::size_t S, std::size_t A, class U> \
188 auto operator SYMBOL(const LoopSIMD<T,S,A> &v, const U s) { \
189 LoopSIMD<T,S,A> out; \
190 for(std::size_t i=0; i<S; i++){ \
191 out[i] = v[i] SYMBOL s; \
192 } \
193 return out; \
194 } \
195 template<class T, std::size_t S, std::size_t A, class U, std::size_t AU> \
196 auto operator SYMBOL(const LoopSIMD<T,S,A> &v, \
197 const LoopSIMD<U,S,AU> &w) { \
198 LoopSIMD<T,S,A> out; \
199 for(std::size_t i=0; i<S; i++){ \
200 out[i] = v[i] SYMBOL w[i]; \
201 } \
202 return out; \
203 } \
204 static_assert(true, "expecting ;")
205
206 DUNE_SIMD_LOOP_BITSHIFT_OP(<<);
207 DUNE_SIMD_LOOP_BITSHIFT_OP(>>);
208
209#undef DUNE_SIMD_LOOP_BITSHIFT_OP
210
211 //Comparison operators
212#define DUNE_SIMD_LOOP_COMPARISON_OP(SYMBOL) \
213 template<class T, std::size_t S, std::size_t A, class U> \
214 auto operator SYMBOL(const LoopSIMD<T,S,A> &v, const U s) { \
215 Simd::Mask<LoopSIMD<T,S,A>> out; \
216 for(std::size_t i=0; i<S; i++){ \
217 out[i] = v[i] SYMBOL s; \
218 } \
219 return out; \
220 } \
221 template<class T, std::size_t S, std::size_t A> \
222 auto operator SYMBOL(const Simd::Scalar<T> s, const LoopSIMD<T,S,A> &v) { \
223 Simd::Mask<LoopSIMD<T,S,A>> out; \
224 for(std::size_t i=0; i<S; i++){ \
225 out[i] = s SYMBOL v[i]; \
226 } \
227 return out; \
228 } \
229 template<class T, std::size_t S, std::size_t A> \
230 auto operator SYMBOL(const LoopSIMD<T,S,A> &v, \
231 const LoopSIMD<T,S,A> &w) { \
232 Simd::Mask<LoopSIMD<T,S,A>> out; \
233 for(std::size_t i=0; i<S; i++){ \
234 out[i] = v[i] SYMBOL w[i]; \
235 } \
236 return out; \
237 } \
238 static_assert(true, "expecting ;")
239
240 DUNE_SIMD_LOOP_COMPARISON_OP(<);
241 DUNE_SIMD_LOOP_COMPARISON_OP(>);
242 DUNE_SIMD_LOOP_COMPARISON_OP(<=);
243 DUNE_SIMD_LOOP_COMPARISON_OP(>=);
244 DUNE_SIMD_LOOP_COMPARISON_OP(==);
245 DUNE_SIMD_LOOP_COMPARISON_OP(!=);
246#undef DUNE_SIMD_LOOP_COMPARISON_OP
247
248 //Boolean operators
249#define DUNE_SIMD_LOOP_BOOLEAN_OP(SYMBOL) \
250 template<class T, std::size_t S, std::size_t A> \
251 auto operator SYMBOL(const LoopSIMD<T,S,A> &v, const Simd::Scalar<T> s) { \
252 Simd::Mask<LoopSIMD<T,S,A>> out; \
253 for(std::size_t i=0; i<S; i++){ \
254 out[i] = v[i] SYMBOL s; \
255 } \
256 return out; \
257 } \
258 template<class T, std::size_t S, std::size_t A> \
259 auto operator SYMBOL(const Simd::Mask<T> s, const LoopSIMD<T,S,A> &v) { \
260 Simd::Mask<LoopSIMD<T,S,A>> out; \
261 for(std::size_t i=0; i<S; i++){ \
262 out[i] = s SYMBOL v[i]; \
263 } \
264 return out; \
265 } \
266 template<class T, std::size_t S, std::size_t A> \
267 auto operator SYMBOL(const LoopSIMD<T,S,A> &v, \
268 const LoopSIMD<T,S,A> &w) { \
269 Simd::Mask<LoopSIMD<T,S,A>> out; \
270 for(std::size_t i=0; i<S; i++){ \
271 out[i] = v[i] SYMBOL w[i]; \
272 } \
273 return out; \
274 } \
275 static_assert(true, "expecting ;")
276
277 DUNE_SIMD_LOOP_BOOLEAN_OP(&&);
278 DUNE_SIMD_LOOP_BOOLEAN_OP(||);
279#undef DUNE_SIMD_LOOP_BOOLEAN_OP
280
281 //prints a given LoopSIMD
282 template<class T, std::size_t S, std::size_t A>
283 std::ostream& operator<< (std::ostream &os, const LoopSIMD<T,S,A> &v) {
284 os << "[";
285 for(std::size_t i=0; i<S-1; i++) {
286 os << v[i] << ", ";
287 }
288 os << v[S-1] << "]";
289 return os;
290 }
291
292 namespace Simd {
293 namespace Overloads {
294 /*
295 * Implementation/Overloads of the functions needed for
296 * SIMD-interface-compatibility
297 */
298
299 //Implementation of SIMD-interface-types
300 template<class T, std::size_t S, std::size_t A>
301 struct ScalarType<LoopSIMD<T,S,A>> {
302 using type = Simd::Scalar<T>;
303 };
304
305 template<class U, class T, std::size_t S, std::size_t A>
306 struct RebindType<U, LoopSIMD<T,S,A>> {
307 using type = LoopSIMD<Simd::Rebind<U, T>,S,A>;
308 };
309
310 //Implementation of SIMD-interface-functionality
311 template<class T, std::size_t S, std::size_t A>
312 struct LaneCount<LoopSIMD<T,S,A>> : index_constant<S*lanes<T>()> {};
313
314 template<class T, std::size_t S, std::size_t A>
315 auto lane(ADLTag<5>, std::size_t l, LoopSIMD<T,S,A> &&v)
316 -> decltype(std::move(Simd::lane(l%lanes<T>(), v[l/lanes<T>()])))
317 {
318 return std::move(Simd::lane(l%lanes<T>(), v[l/lanes<T>()]));
319 }
320
321 template<class T, std::size_t S, std::size_t A>
322 auto lane(ADLTag<5>, std::size_t l, const LoopSIMD<T,S,A> &v)
323 -> decltype(Simd::lane(l%lanes<T>(), v[l/lanes<T>()]))
324 {
325 return Simd::lane(l%lanes<T>(), v[l/lanes<T>()]);
326 }
327
328 template<class T, std::size_t S, std::size_t A>
329 auto lane(ADLTag<5>, std::size_t l, LoopSIMD<T,S,A> &v)
330 -> decltype(Simd::lane(l%lanes<T>(), v[l/lanes<T>()]))
331 {
332 return Simd::lane(l%lanes<T>(), v[l/lanes<T>()]);
333 }
334
335 template<class T, std::size_t S, std::size_t AM, std::size_t AD>
336 auto cond(ADLTag<5>, Simd::Mask<LoopSIMD<T,S,AM>> mask,
337 LoopSIMD<T,S,AD> ifTrue, LoopSIMD<T,S,AD> ifFalse) {
338 LoopSIMD<T,S,AD> out;
339 for(std::size_t i=0; i<S; i++) {
340 out[i] = Simd::cond(mask[i], ifTrue[i], ifFalse[i]);
341 }
342 return out;
343 }
344
345 template<class M, class T, std::size_t S, std::size_t A>
346 auto cond(ADLTag<5, std::is_same<bool, Simd::Scalar<M> >::value
347 && Simd::lanes<M>() == Simd::lanes<LoopSIMD<T,S,A> >()>,
348 M mask, LoopSIMD<T,S,A> ifTrue, LoopSIMD<T,S,A> ifFalse)
349 {
350 LoopSIMD<T,S,A> out;
351 for(auto l : range(Simd::lanes(mask)))
352 Simd::lane(l, out) = Simd::lane(l, mask) ? Simd::lane(l, ifTrue) : Simd::lane(l, ifFalse);
353 return out;
354 }
355
356 template<class M, std::size_t S, std::size_t A>
357 bool anyTrue(ADLTag<5>, LoopSIMD<M,S,A> mask) {
358 bool out = false;
359 for(std::size_t i=0; i<S; i++) {
360 out |= Simd::anyTrue(mask[i]);
361 }
362 return out;
363 }
364
365 template<class M, std::size_t S, std::size_t A>
366 bool allTrue(ADLTag<5>, LoopSIMD<M,S,A> mask) {
367 bool out = true;
368 for(std::size_t i=0; i<S; i++) {
369 out &= Simd::allTrue(mask[i]);
370 }
371 return out;
372 }
373
374 template<class M, std::size_t S, std::size_t A>
375 bool anyFalse(ADLTag<5>, LoopSIMD<M,S,A> mask) {
376 bool out = false;
377 for(std::size_t i=0; i<S; i++) {
378 out |= Simd::anyFalse(mask[i]);
379 }
380 return out;
381 }
382
383 template<class M, std::size_t S, std::size_t A>
384 bool allFalse(ADLTag<5>, LoopSIMD<M,S,A> mask) {
385 bool out = true;
386 for(std::size_t i=0; i<S; i++) {
387 out &= Simd::allFalse(mask[i]);
388 }
389 return out;
390 }
391 } //namespace Overloads
392
393 } //namespace Simd
394
395
396 /*
397 * Overloads the unary cmath-operations. Operations requiring
398 * or returning more than one argument are not supported.
399 * Due to inconsistency with the return values, cmath-operations
400 * on integral types are also not supported-
401 */
402
403#define DUNE_SIMD_LOOP_CMATH_UNARY_OP(expr) \
404 template<class T, std::size_t S, std::size_t A, typename Sfinae = \
405 typename std::enable_if_t<!std::is_integral<Simd::Scalar<T>>::value> > \
406 auto expr(const LoopSIMD<T,S,A> &v) { \
407 using std::expr; \
408 LoopSIMD<T,S,A> out; \
409 for(std::size_t i=0; i<S; i++) { \
410 out[i] = expr(v[i]); \
411 } \
412 return out; \
413 } \
414 static_assert(true, "expecting ;")
415
416#define DUNE_SIMD_LOOP_CMATH_UNARY_OP_WITH_RETURN(expr, returnType) \
417 template<class T, std::size_t S, std::size_t A, typename Sfinae = \
418 typename std::enable_if_t<!std::is_integral<Simd::Scalar<T>>::value> > \
419 auto expr(const LoopSIMD<T,S,A> &v) { \
420 using std::expr; \
421 LoopSIMD<returnType,S> out; \
422 for(std::size_t i=0; i<S; i++) { \
423 out[i] = expr(v[i]); \
424 } \
425 return out; \
426 } \
427 static_assert(true, "expecting ;")
428
429 DUNE_SIMD_LOOP_CMATH_UNARY_OP(cos);
430 DUNE_SIMD_LOOP_CMATH_UNARY_OP(sin);
431 DUNE_SIMD_LOOP_CMATH_UNARY_OP(tan);
432 DUNE_SIMD_LOOP_CMATH_UNARY_OP(acos);
433 DUNE_SIMD_LOOP_CMATH_UNARY_OP(asin);
434 DUNE_SIMD_LOOP_CMATH_UNARY_OP(atan);
435 DUNE_SIMD_LOOP_CMATH_UNARY_OP(cosh);
436 DUNE_SIMD_LOOP_CMATH_UNARY_OP(sinh);
437 DUNE_SIMD_LOOP_CMATH_UNARY_OP(tanh);
438 DUNE_SIMD_LOOP_CMATH_UNARY_OP(acosh);
439 DUNE_SIMD_LOOP_CMATH_UNARY_OP(asinh);
440 DUNE_SIMD_LOOP_CMATH_UNARY_OP(atanh);
441
442 DUNE_SIMD_LOOP_CMATH_UNARY_OP(exp);
443 DUNE_SIMD_LOOP_CMATH_UNARY_OP(log);
444 DUNE_SIMD_LOOP_CMATH_UNARY_OP(log10);
445 DUNE_SIMD_LOOP_CMATH_UNARY_OP(exp2);
446 DUNE_SIMD_LOOP_CMATH_UNARY_OP(expm1);
447 DUNE_SIMD_LOOP_CMATH_UNARY_OP_WITH_RETURN(ilogb, int);
448 DUNE_SIMD_LOOP_CMATH_UNARY_OP(log1p);
449 DUNE_SIMD_LOOP_CMATH_UNARY_OP(log2);
450 DUNE_SIMD_LOOP_CMATH_UNARY_OP(logb);
451
452 DUNE_SIMD_LOOP_CMATH_UNARY_OP(sqrt);
453 DUNE_SIMD_LOOP_CMATH_UNARY_OP(cbrt);
454
455 DUNE_SIMD_LOOP_CMATH_UNARY_OP(erf);
456 DUNE_SIMD_LOOP_CMATH_UNARY_OP(erfc);
457 DUNE_SIMD_LOOP_CMATH_UNARY_OP(tgamma);
458 DUNE_SIMD_LOOP_CMATH_UNARY_OP(lgamma);
459
460 DUNE_SIMD_LOOP_CMATH_UNARY_OP(ceil);
461 DUNE_SIMD_LOOP_CMATH_UNARY_OP(floor);
462 DUNE_SIMD_LOOP_CMATH_UNARY_OP(trunc);
463 DUNE_SIMD_LOOP_CMATH_UNARY_OP(round);
464 DUNE_SIMD_LOOP_CMATH_UNARY_OP_WITH_RETURN(lround, long);
465 DUNE_SIMD_LOOP_CMATH_UNARY_OP_WITH_RETURN(llround, long long);
466 DUNE_SIMD_LOOP_CMATH_UNARY_OP(rint);
467 DUNE_SIMD_LOOP_CMATH_UNARY_OP_WITH_RETURN(lrint, long);
468 DUNE_SIMD_LOOP_CMATH_UNARY_OP_WITH_RETURN(llrint, long long);
469 DUNE_SIMD_LOOP_CMATH_UNARY_OP(nearbyint);
470
471 DUNE_SIMD_LOOP_CMATH_UNARY_OP(fabs);
472 DUNE_SIMD_LOOP_CMATH_UNARY_OP(abs);
473
474#undef DUNE_SIMD_LOOP_CMATH_UNARY_OP
475#undef DUNE_SIMD_LOOP_CMATH_UNARY_OP_WITH_RETURN
476
477
478 /* not implemented cmath-functions:
479 * atan2
480 * frexp, idexp
481 * modf
482 * scalbn, scalbln
483 * pow
484 * hypot
485 * remainder, remquo
486 * copysign
487 * nan
488 * nextafter, nexttoward
489 * fdim, fmax, fmin
490 */
491
492 /*
493 * Overloads specific functions usually provided by the std library
494 * More overloads will be provided should the need arise.
495 */
496
497#define DUNE_SIMD_LOOP_STD_UNARY_OP(expr) \
498 template<class T, std::size_t S, std::size_t A> \
499 auto expr(const LoopSIMD<T,S,A> &v) { \
500 using std::expr; \
501 LoopSIMD<T,S,A> out; \
502 for(std::size_t i=0; i<S; i++) { \
503 out[i] = expr(v[i]); \
504 } \
505 return out; \
506 } \
507 \
508 template<class T, std::size_t S, std::size_t A> \
509 auto expr(const LoopSIMD<std::complex<T>,S,A> &v) { \
510 using std::expr; \
511 LoopSIMD<T,S,A> out; \
512 for(std::size_t i=0; i<S; i++) { \
513 out[i] = expr(v[i]); \
514 } \
515 return out; \
516 } \
517 static_assert(true, "expecting ;")
518
519 DUNE_SIMD_LOOP_STD_UNARY_OP(real);
520 DUNE_SIMD_LOOP_STD_UNARY_OP(imag);
521
522#undef DUNE_SIMD_LOOP_STD_UNARY_OP
523
524#define DUNE_SIMD_LOOP_STD_BINARY_OP(expr) \
525 template<class T, std::size_t S, std::size_t A> \
526 auto expr(const LoopSIMD<T,S,A> &v, const LoopSIMD<T,S,A> &w) { \
527 using std::expr; \
528 LoopSIMD<T,S,A> out; \
529 for(std::size_t i=0; i<S; i++) { \
530 out[i] = expr(v[i],w[i]); \
531 } \
532 return out; \
533 } \
534 static_assert(true, "expecting ;")
535
536 DUNE_SIMD_LOOP_STD_BINARY_OP(max);
537 DUNE_SIMD_LOOP_STD_BINARY_OP(min);
538
539#undef DUNE_SIMD_LOOP_STD_BINARY_OP
540
541 namespace MathOverloads {
542 template<class T, std::size_t S, std::size_t A>
543 auto isNaN(const LoopSIMD<T,S,A> &v, PriorityTag<3>, ADLTag) {
544 Simd::Mask<LoopSIMD<T,S,A>> out;
545 for(auto l : range(S))
546 out[l] = Dune::isNaN(v[l]);
547 return out;
548 }
549
550 template<class T, std::size_t S, std::size_t A>
551 auto isInf(const LoopSIMD<T,S,A> &v, PriorityTag<3>, ADLTag) {
552 Simd::Mask<LoopSIMD<T,S,A>> out;
553 for(auto l : range(S))
554 out[l] = Dune::isInf(v[l]);
555 return out;
556 }
557
558 template<class T, std::size_t S, std::size_t A>
559 auto isFinite(const LoopSIMD<T,S,A> &v, PriorityTag<3>, ADLTag) {
560 Simd::Mask<LoopSIMD<T,S,A>> out;
561 for(auto l : range(S))
562 out[l] = Dune::isFinite(v[l]);
563 return out;
564 }
565 } //namepace MathOverloads
566
567 template<class T, std::size_t S, std::size_t A>
568 struct IsNumber<LoopSIMD<T,S,A>> :
569 public std::integral_constant<bool, IsNumber<T>::value>{
570 };
571
572#if __GNUC__ >= 7
573# pragma GCC diagnostic pop
574#endif
575
576} //namespace Dune
577
578#endif
Definition: loop.hh:38
Traits for type conversions and type information.
std::integral_constant< std::size_t, i > index_constant
An index constant with value i.
Definition: indices.hh:28
I round(const T &val, typename EpsilonType< T >::Type epsilon)
round using epsilon
Definition: float_cmp.cc:309
I trunc(const T &val, typename EpsilonType< T >::Type epsilon)
truncate using epsilon
Definition: float_cmp.cc:405
Mask< V > mask(ADLTag< 0, std::is_same< V, Mask< V > >::value >, const V &v)
implements Simd::mask()
Definition: defaults.hh:151
bool allFalse(ADLTag< 0 >, const Mask &mask)
implements Simd::allFalse()
Definition: defaults.hh:122
bool allTrue(ADLTag< 0 >, const Mask &mask)
implements Simd::allTrue()
Definition: defaults.hh:102
bool anyFalse(ADLTag< 0 >, const Mask &mask)
implements Simd::anyFalse()
Definition: defaults.hh:112
auto min(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::min()
Definition: defaults.hh:87
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:79
bool anyTrue(const Mask &mask)
Whether any entry is true
Definition: interface.hh:427
V cond(M &&mask, const V &ifTrue, const V &ifFalse)
Like the ?: operator.
Definition: interface.hh:384
bool allTrue(const Mask &mask)
Whether all entries are true
Definition: interface.hh:437
bool anyFalse(const Mask &mask)
Whether any entry is false
Definition: interface.hh:447
constexpr std::size_t lanes()
Number of lanes in a SIMD type.
Definition: interface.hh:303
decltype(auto) lane(std::size_t l, V &&v)
Extract an element of a SIMD type.
Definition: interface.hh:322
Rebind< bool, V > Mask
Mask type type of some SIMD type.
Definition: interface.hh:287
bool allFalse(const Mask &mask)
Whether all entries are false
Definition: interface.hh:457
typename Overloads::ScalarType< std::decay_t< V > >::type Scalar
Element type of some SIMD type.
Definition: interface.hh:233
Some useful basic math stuff.
Dune namespace.
Definition: alignedallocator.hh:11
const T1 cond(bool b, const T1 &v1, const T2 &v2)
conditional evaluate
Definition: conditional.hh:26
Include file for users of the SIMD abstraction layer.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 24, 22:29, 2024)