dune-common 2.9.0
loop.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightInfo: Copyright (C) DUNE Project contributors, see file LICENSE.md in module root
2// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3#ifndef DUNE_COMMON_SIMD_LOOP_HH
4#define DUNE_COMMON_SIMD_LOOP_HH
5
6#include <array>
7#include <cmath>
8#include <cstddef>
9#include <cstdlib>
10#include <ostream>
11
12#include <dune/common/math.hh>
15
16namespace Dune {
17
18
19/*
20 * silence warnings from GCC about using integer operands on a bool
21 * (when instantiated for T=bool)
22 */
23#if __GNUC__ >= 7
24# pragma GCC diagnostic push
25# pragma GCC diagnostic ignored "-Wbool-operation"
26# pragma GCC diagnostic ignored "-Wint-in-bool-context"
27# define GCC_WARNING_DISABLED
28#endif
29
30/*
31 * silence warnings from Clang about using bitwise operands on
32 * a bool (when instantiated for T=bool)
33 */
34#ifdef __clang__
35#if __has_warning("-Wbitwise-instead-of-logical")
36# pragma clang diagnostic push
37# pragma clang diagnostic ignored "-Wbitwise-instead-of-logical"
38# define CLANG_WARNING_DISABLED
39#endif
40#endif
41
42/*
43 * Introduce a simd pragma if OpenMP is available in standard version >= 4
44 */
45#if _OPENMP >= 201307
46 #define DUNE_PRAGMA_OMP_SIMD _Pragma("omp simd")
47#else
48 #define DUNE_PRAGMA_OMP_SIMD
49#endif
50
51
63 template<class T, std::size_t S, std::size_t A = 0>
64 class alignas(A==0?alignof(T):A) LoopSIMD : public std::array<T,S> {
65
66 public:
67
68 //default constructor
70 assert(reinterpret_cast<uintptr_t>(this) % std::min(alignof(LoopSIMD<T,S,A>),alignof(std::max_align_t)) == 0);
71 }
72
73 // broadcast constructor initializing the content with a given value
75 this->fill(i);
76 }
77
78 template<std::size_t OA>
79 explicit LoopSIMD(const LoopSIMD<T,S,OA>& other)
80 : std::array<T,S>(other)
81 {
82 assert(reinterpret_cast<uintptr_t>(this) % std::min(alignof(LoopSIMD<T,S,A>),alignof(std::max_align_t)) == 0);
83 }
84
85 /*
86 * Definition of basic operators
87 */
88
89 //Prefix operators
90#define DUNE_SIMD_LOOP_PREFIX_OP(SYMBOL) \
91 auto operator SYMBOL() { \
92 DUNE_PRAGMA_OMP_SIMD \
93 for(std::size_t i=0; i<S; i++){ \
94 SYMBOL(*this)[i]; \
95 } \
96 return *this; \
97 } \
98 static_assert(true, "expecting ;")
99
102#undef DUNE_SIMD_LOOP_PREFIX_OP
103
104 //Unary operators
105#define DUNE_SIMD_LOOP_UNARY_OP(SYMBOL) \
106 auto operator SYMBOL() const { \
107 LoopSIMD<T,S,A> out; \
108 DUNE_PRAGMA_OMP_SIMD \
109 for(std::size_t i=0; i<S; i++){ \
110 out[i] = SYMBOL((*this)[i]); \
111 } \
112 return out; \
113 } \
114 static_assert(true, "expecting ;")
115
119
120 auto operator!() const {
123 for(std::size_t i=0; i<S; i++){
124 out[i] = !((*this)[i]);
125 }
126 return out;
127 }
128#undef DUNE_SIMD_LOOP_UNARY_OP
129
130 //Postfix operators
131#define DUNE_SIMD_LOOP_POSTFIX_OP(SYMBOL) \
132 auto operator SYMBOL(int){ \
133 LoopSIMD<T,S,A> out = *this; \
134 SYMBOL(*this); \
135 return out; \
136 } \
137 static_assert(true, "expecting ;")
138
141#undef DUNE_SIMD_LOOP_POSTFIX_OP
142
143 //Assignment operators
144#define DUNE_SIMD_LOOP_ASSIGNMENT_OP(SYMBOL) \
145 auto operator SYMBOL(const Simd::Scalar<T> s) { \
146 DUNE_PRAGMA_OMP_SIMD \
147 for(std::size_t i=0; i<S; i++){ \
148 (*this)[i] SYMBOL s; \
149 } \
150 return *this; \
151 } \
152 \
153 auto operator SYMBOL(const LoopSIMD<T,S,A> &v) { \
154 DUNE_PRAGMA_OMP_SIMD \
155 for(std::size_t i=0; i<S; i++){ \
156 (*this)[i] SYMBOL v[i]; \
157 } \
158 return *this; \
159 } \
160 static_assert(true, "expecting ;")
161
172#undef DUNE_SIMD_LOOP_ASSIGNMENT_OP
173 };
174
175 //Arithmetic operators
176#define DUNE_SIMD_LOOP_BINARY_OP(SYMBOL) \
177 template<class T, std::size_t S, std::size_t A> \
178 auto operator SYMBOL(const LoopSIMD<T,S,A> &v, const Simd::Scalar<T> s) { \
179 LoopSIMD<T,S,A> out; \
180 DUNE_PRAGMA_OMP_SIMD \
181 for(std::size_t i=0; i<S; i++){ \
182 out[i] = v[i] SYMBOL s; \
183 } \
184 return out; \
185 } \
186 template<class T, std::size_t S, std::size_t A> \
187 auto operator SYMBOL(const Simd::Scalar<T> s, const LoopSIMD<T,S,A> &v) { \
188 LoopSIMD<T,S,A> out; \
189 DUNE_PRAGMA_OMP_SIMD \
190 for(std::size_t i=0; i<S; i++){ \
191 out[i] = s SYMBOL v[i]; \
192 } \
193 return out; \
194 } \
195 template<class T, std::size_t S, std::size_t A> \
196 auto operator SYMBOL(const LoopSIMD<T,S,A> &v, \
197 const LoopSIMD<T,S,A> &w) { \
198 LoopSIMD<T,S,A> out; \
199 DUNE_PRAGMA_OMP_SIMD \
200 for(std::size_t i=0; i<S; i++){ \
201 out[i] = v[i] SYMBOL w[i]; \
202 } \
203 return out; \
204 } \
205 static_assert(true, "expecting ;")
206
212
216
217#undef DUNE_SIMD_LOOP_BINARY_OP
218
219 //Bitshift operators
220#define DUNE_SIMD_LOOP_BITSHIFT_OP(SYMBOL) \
221 template<class T, std::size_t S, std::size_t A, class U> \
222 auto operator SYMBOL(const LoopSIMD<T,S,A> &v, const U s) { \
223 LoopSIMD<T,S,A> out; \
224 DUNE_PRAGMA_OMP_SIMD \
225 for(std::size_t i=0; i<S; i++){ \
226 out[i] = v[i] SYMBOL s; \
227 } \
228 return out; \
229 } \
230 template<class T, std::size_t S, std::size_t A, class U, std::size_t AU> \
231 auto operator SYMBOL(const LoopSIMD<T,S,A> &v, \
232 const LoopSIMD<U,S,AU> &w) { \
233 LoopSIMD<T,S,A> out; \
234 DUNE_PRAGMA_OMP_SIMD \
235 for(std::size_t i=0; i<S; i++){ \
236 out[i] = v[i] SYMBOL w[i]; \
237 } \
238 return out; \
239 } \
240 static_assert(true, "expecting ;")
241
244
245#undef DUNE_SIMD_LOOP_BITSHIFT_OP
246
247 //Comparison operators
248#define DUNE_SIMD_LOOP_COMPARISON_OP(SYMBOL) \
249 template<class T, std::size_t S, std::size_t A, class U> \
250 auto operator SYMBOL(const LoopSIMD<T,S,A> &v, const U s) { \
251 Simd::Mask<LoopSIMD<T,S,A>> out; \
252 DUNE_PRAGMA_OMP_SIMD \
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::Scalar<T> s, const LoopSIMD<T,S,A> &v) { \
260 Simd::Mask<LoopSIMD<T,S,A>> out; \
261 DUNE_PRAGMA_OMP_SIMD \
262 for(std::size_t i=0; i<S; i++){ \
263 out[i] = s SYMBOL v[i]; \
264 } \
265 return out; \
266 } \
267 template<class T, std::size_t S, std::size_t A> \
268 auto operator SYMBOL(const LoopSIMD<T,S,A> &v, \
269 const LoopSIMD<T,S,A> &w) { \
270 Simd::Mask<LoopSIMD<T,S,A>> out; \
271 DUNE_PRAGMA_OMP_SIMD \
272 for(std::size_t i=0; i<S; i++){ \
273 out[i] = v[i] SYMBOL w[i]; \
274 } \
275 return out; \
276 } \
277 static_assert(true, "expecting ;")
278
285#undef DUNE_SIMD_LOOP_COMPARISON_OP
286
287 //Boolean operators
288#define DUNE_SIMD_LOOP_BOOLEAN_OP(SYMBOL) \
289 template<class T, std::size_t S, std::size_t A> \
290 auto operator SYMBOL(const LoopSIMD<T,S,A> &v, const Simd::Scalar<T> s) { \
291 Simd::Mask<LoopSIMD<T,S,A>> out; \
292 DUNE_PRAGMA_OMP_SIMD \
293 for(std::size_t i=0; i<S; i++){ \
294 out[i] = v[i] SYMBOL s; \
295 } \
296 return out; \
297 } \
298 template<class T, std::size_t S, std::size_t A> \
299 auto operator SYMBOL(const Simd::Mask<T> s, const LoopSIMD<T,S,A> &v) { \
300 Simd::Mask<LoopSIMD<T,S,A>> out; \
301 DUNE_PRAGMA_OMP_SIMD \
302 for(std::size_t i=0; i<S; i++){ \
303 out[i] = s SYMBOL v[i]; \
304 } \
305 return out; \
306 } \
307 template<class T, std::size_t S, std::size_t A> \
308 auto operator SYMBOL(const LoopSIMD<T,S,A> &v, \
309 const LoopSIMD<T,S,A> &w) { \
310 Simd::Mask<LoopSIMD<T,S,A>> out; \
311 DUNE_PRAGMA_OMP_SIMD \
312 for(std::size_t i=0; i<S; i++){ \
313 out[i] = v[i] SYMBOL w[i]; \
314 } \
315 return out; \
316 } \
317 static_assert(true, "expecting ;")
318
321#undef DUNE_SIMD_LOOP_BOOLEAN_OP
322
323 //prints a given LoopSIMD
324 template<class T, std::size_t S, std::size_t A>
325 std::ostream& operator<< (std::ostream &os, const LoopSIMD<T,S,A> &v) {
326 os << "[";
327 for(std::size_t i=0; i<S-1; i++) {
328 os << v[i] << ", ";
329 }
330 os << v[S-1] << "]";
331 return os;
332 }
333
334 namespace Simd {
335 namespace Overloads {
336 /*
337 * Implementation/Overloads of the functions needed for
338 * SIMD-interface-compatibility
339 */
340
341 //Implementation of SIMD-interface-types
342 template<class T, std::size_t S, std::size_t A>
343 struct ScalarType<LoopSIMD<T,S,A>> {
345 };
346
347 template<class U, class T, std::size_t S, std::size_t A>
348 struct RebindType<U, LoopSIMD<T,S,A>> {
350 };
351
352 //Implementation of SIMD-interface-functionality
353 template<class T, std::size_t S, std::size_t A>
354 struct LaneCount<LoopSIMD<T,S,A>> : index_constant<S*lanes<T>()> {};
355
356 template<class T, std::size_t S, std::size_t A>
357 auto lane(ADLTag<5>, std::size_t l, LoopSIMD<T,S,A> &&v)
358 -> decltype(std::move(Simd::lane(l%lanes<T>(), v[l/lanes<T>()])))
359 {
360 return std::move(Simd::lane(l%lanes<T>(), v[l/lanes<T>()]));
361 }
362
363 template<class T, std::size_t S, std::size_t A>
364 auto lane(ADLTag<5>, std::size_t l, const LoopSIMD<T,S,A> &v)
365 -> decltype(Simd::lane(l%lanes<T>(), v[l/lanes<T>()]))
366 {
367 return Simd::lane(l%lanes<T>(), v[l/lanes<T>()]);
368 }
369
370 template<class T, std::size_t S, std::size_t A>
371 auto lane(ADLTag<5>, std::size_t l, LoopSIMD<T,S,A> &v)
372 -> decltype(Simd::lane(l%lanes<T>(), v[l/lanes<T>()]))
373 {
374 return Simd::lane(l%lanes<T>(), v[l/lanes<T>()]);
375 }
376
377 template<class T, std::size_t S, std::size_t AM, std::size_t AD>
379 LoopSIMD<T,S,AD> ifTrue, LoopSIMD<T,S,AD> ifFalse) {
381 for(std::size_t i=0; i<S; i++) {
382 out[i] = Simd::cond(mask[i], ifTrue[i], ifFalse[i]);
383 }
384 return out;
385 }
386
387 template<class M, class T, std::size_t S, std::size_t A>
388 auto cond(ADLTag<5, std::is_same<bool, Simd::Scalar<M> >::value
389 && Simd::lanes<M>() == Simd::lanes<LoopSIMD<T,S,A> >()>,
390 M mask, LoopSIMD<T,S,A> ifTrue, LoopSIMD<T,S,A> ifFalse)
391 {
392 LoopSIMD<T,S,A> out;
393 for(auto l : range(Simd::lanes(mask)))
394 Simd::lane(l, out) = Simd::lane(l, mask) ? Simd::lane(l, ifTrue) : Simd::lane(l, ifFalse);
395 return out;
396 }
397
398 template<class M, std::size_t S, std::size_t A>
400 bool out = false;
401 for(std::size_t i=0; i<S; i++) {
402 out |= Simd::anyTrue(mask[i]);
403 }
404 return out;
405 }
406
407 template<class M, std::size_t S, std::size_t A>
409 bool out = true;
410 for(std::size_t i=0; i<S; i++) {
411 out &= Simd::allTrue(mask[i]);
412 }
413 return out;
414 }
415
416 template<class M, std::size_t S, std::size_t A>
418 bool out = false;
419 for(std::size_t i=0; i<S; i++) {
420 out |= Simd::anyFalse(mask[i]);
421 }
422 return out;
423 }
424
425 template<class M, std::size_t S, std::size_t A>
427 bool out = true;
428 for(std::size_t i=0; i<S; i++) {
429 out &= Simd::allFalse(mask[i]);
430 }
431 return out;
432 }
433 } //namespace Overloads
434
435 } //namespace Simd
436
437
438 /*
439 * Overloads the unary cmath-operations. Operations requiring
440 * or returning more than one argument are not supported.
441 * Due to inconsistency with the return values, cmath-operations
442 * on integral types are also not supported-
443 */
444
445#define DUNE_SIMD_LOOP_CMATH_UNARY_OP(expr) \
446 template<class T, std::size_t S, std::size_t A, typename Sfinae = \
447 typename std::enable_if_t<!std::is_integral<Simd::Scalar<T>>::value> > \
448 auto expr(const LoopSIMD<T,S,A> &v) { \
449 using std::expr; \
450 LoopSIMD<T,S,A> out; \
451 for(std::size_t i=0; i<S; i++) { \
452 out[i] = expr(v[i]); \
453 } \
454 return out; \
455 } \
456 static_assert(true, "expecting ;")
457
458#define DUNE_SIMD_LOOP_CMATH_UNARY_OP_WITH_RETURN(expr, returnType) \
459 template<class T, std::size_t S, std::size_t A, typename Sfinae = \
460 typename std::enable_if_t<!std::is_integral<Simd::Scalar<T>>::value> > \
461 auto expr(const LoopSIMD<T,S,A> &v) { \
462 using std::expr; \
463 LoopSIMD<returnType,S> out; \
464 for(std::size_t i=0; i<S; i++) { \
465 out[i] = expr(v[i]); \
466 } \
467 return out; \
468 } \
469 static_assert(true, "expecting ;")
470
483
493
496
501
512
515
516#undef DUNE_SIMD_LOOP_CMATH_UNARY_OP
517#undef DUNE_SIMD_LOOP_CMATH_UNARY_OP_WITH_RETURN
518
519
520 /* not implemented cmath-functions:
521 * atan2
522 * frexp, idexp
523 * modf
524 * scalbn, scalbln
525 * pow
526 * hypot
527 * remainder, remquo
528 * copysign
529 * nan
530 * nextafter, nexttoward
531 * fdim, fmax, fmin
532 */
533
534 /*
535 * Overloads specific functions usually provided by the std library
536 * More overloads will be provided should the need arise.
537 */
538
539#define DUNE_SIMD_LOOP_STD_UNARY_OP(expr) \
540 template<class T, std::size_t S, std::size_t A> \
541 auto expr(const LoopSIMD<T,S,A> &v) { \
542 using std::expr; \
543 LoopSIMD<T,S,A> out; \
544 for(std::size_t i=0; i<S; i++) { \
545 out[i] = expr(v[i]); \
546 } \
547 return out; \
548 } \
549 \
550 template<class T, std::size_t S, std::size_t A> \
551 auto expr(const LoopSIMD<std::complex<T>,S,A> &v) { \
552 using std::expr; \
553 LoopSIMD<T,S,A> out; \
554 for(std::size_t i=0; i<S; i++) { \
555 out[i] = expr(v[i]); \
556 } \
557 return out; \
558 } \
559 static_assert(true, "expecting ;")
560
563
564#undef DUNE_SIMD_LOOP_STD_UNARY_OP
565
566#define DUNE_SIMD_LOOP_STD_BINARY_OP(expr) \
567 template<class T, std::size_t S, std::size_t A> \
568 auto expr(const LoopSIMD<T,S,A> &v, const LoopSIMD<T,S,A> &w) { \
569 using std::expr; \
570 LoopSIMD<T,S,A> out; \
571 for(std::size_t i=0; i<S; i++) { \
572 out[i] = expr(v[i],w[i]); \
573 } \
574 return out; \
575 } \
576 static_assert(true, "expecting ;")
577
580
581#undef DUNE_SIMD_LOOP_STD_BINARY_OP
582
583 namespace MathOverloads {
584 template<class T, std::size_t S, std::size_t A>
587 for(auto l : range(S))
588 out[l] = Dune::isNaN(v[l]);
589 return out;
590 }
591
592 template<class T, std::size_t S, std::size_t A>
595 for(auto l : range(S))
596 out[l] = Dune::isInf(v[l]);
597 return out;
598 }
599
600 template<class T, std::size_t S, std::size_t A>
603 for(auto l : range(S))
604 out[l] = Dune::isFinite(v[l]);
605 return out;
606 }
607 } //namespace MathOverloads
608
609 template<class T, std::size_t S, std::size_t A>
610 struct IsNumber<LoopSIMD<T,S,A>> :
611 public std::integral_constant<bool, IsNumber<T>::value>{
612 };
613
614#ifdef CLANG_WARNING_DISABLED
615# pragma clang diagnostic pop
616# undef CLANG_WARNING_DISABLED
617#endif
618
619#ifdef GCC_WARNING_DISABLED
620# pragma GCC diagnostic pop
621# undef GCC_WARNING_DISABLED
622#endif
623
624} //namespace Dune
625
626#endif
Traits for type conversions and type information.
#define DUNE_PRAGMA_OMP_SIMD
Definition: loop.hh:48
Some useful basic math stuff.
std::integral_constant< std::size_t, i > index_constant
An index constant with value i.
Definition: indices.hh:30
static StaticIntegralRange< T, to, from > range(std::integral_constant< T, from >, std::integral_constant< T, to >) noexcept
Definition: rangeutilities.hh:300
std::ostream & operator<<(std::ostream &s, const bigunsignedint< k > &x)
Definition: bigunsignedint.hh:278
I round(const T &val, typename EpsilonType< T >::Type epsilon)
round using epsilon
Definition: float_cmp.cc:311
I trunc(const T &val, typename EpsilonType< T >::Type epsilon)
truncate using epsilon
Definition: float_cmp.cc:407
bool anyTrue(const Mask &mask)
Whether any entry is true
Definition: simd/interface.hh:429
V cond(M &&mask, const V &ifTrue, const V &ifFalse)
Like the ?: operator.
Definition: simd/interface.hh:386
bool allTrue(const Mask &mask)
Whether all entries are true
Definition: simd/interface.hh:439
bool anyFalse(const Mask &mask)
Whether any entry is false
Definition: simd/interface.hh:449
constexpr std::size_t lanes()
Number of lanes in a SIMD type.
Definition: simd/interface.hh:305
decltype(auto) lane(std::size_t l, V &&v)
Extract an element of a SIMD type.
Definition: simd/interface.hh:324
Rebind< bool, V > Mask
Mask type type of some SIMD type.
Definition: simd/interface.hh:289
bool allFalse(const Mask &mask)
Whether all entries are false
Definition: simd/interface.hh:459
typename Overloads::ScalarType< std::decay_t< V > >::type Scalar
Element type of some SIMD type.
Definition: simd/interface.hh:235
Mask< V > mask(ADLTag< 0, std::is_same< V, Mask< V > >::value >, const V &v)
implements Simd::mask()
Definition: defaults.hh:153
bool allFalse(ADLTag< 0 >, const Mask &mask)
implements Simd::allFalse()
Definition: defaults.hh:124
bool allTrue(ADLTag< 0 >, const Mask &mask)
implements Simd::allTrue()
Definition: defaults.hh:104
bool anyFalse(ADLTag< 0 >, const Mask &mask)
implements Simd::anyFalse()
Definition: defaults.hh:114
STL namespace.
Dune namespace.
Definition: alignedallocator.hh:13
DUNE_SIMD_LOOP_CMATH_UNARY_OP_WITH_RETURN(ilogb, int)
DUNE_SIMD_LOOP_STD_BINARY_OP(max)
DUNE_SIMD_LOOP_BOOLEAN_OP && DUNE_SIMD_LOOP_BOOLEAN_OP(||);template< class T, std::size_t S, std::size_t A > std::ostream &operator<<(std::ostream &os, const LoopSIMD< T, S, A > &v
Definition: loop.hh:320
DUNE_SIMD_LOOP_BINARY_OP(+)
DUNE_SIMD_LOOP_CMATH_UNARY_OP(cos)
DUNE_SIMD_LOOP_COMPARISON_OP(<)
DUNE_SIMD_LOOP_BITSHIFT_OP(<<)
DUNE_SIMD_LOOP_STD_UNARY_OP(real)
auto min(const AlignedNumber< T, align > &a, const AlignedNumber< T, align > &b)
Definition: debugalign.hh:447
auto max(const AlignedNumber< T, align > &a, const AlignedNumber< T, align > &b)
Definition: debugalign.hh:425
bool anyTrue(ADLTag< 5 >, const AlignedNumber< bool, align > &mask)
Definition: debugalign.hh:543
const AlignedNumber< T, align > & cond(ADLTag< 5 >, AlignedNumber< bool, align > mask, const AlignedNumber< T, align > &ifTrue, const AlignedNumber< T, align > &ifFalse)
Definition: debugalign.hh:535
T & lane(ADLTag< 5 >, std::size_t l, AlignedNumber< T, align > &v)
Definition: debugalign.hh:520
bool isNaN(const FieldVector< K, SIZE > &b, PriorityTag< 2 >, ADLTag)
Definition: fvector.hh:604
bool isInf(const FieldVector< K, SIZE > &b, PriorityTag< 2 >, ADLTag)
Definition: fvector.hh:594
auto isFinite(const LoopSIMD< T, S, A > &v, PriorityTag< 3 >, ADLTag)
Definition: loop.hh:601
auto isNaN(const LoopSIMD< T, S, A > &v, PriorityTag< 3 >, ADLTag)
Definition: loop.hh:585
auto isFinite(const FieldVector< K, SIZE > &b, PriorityTag< 2 >, ADLTag)
Definition: fvector.hh:584
auto isInf(const LoopSIMD< T, S, A > &v, PriorityTag< 3 >, ADLTag)
Definition: loop.hh:593
Whether this type acts as a scalar in the context of (hierarchically blocked) containers.
Definition: typetraits.hh:194
Tag to make sure the functions in this namespace can be found by ADL.
Definition: math.hh:230
Tag used to force late-binding lookup in Dune::Simd::Overloads.
Definition: base.hh:182
should have a member type type
Definition: standard.hh:60
should have a member type type
Definition: standard.hh:67
should be derived from a Dune::index_constant
Definition: standard.hh:74
Definition: loop.hh:64
LoopSIMD(Simd::Scalar< T > i)
Definition: loop.hh:74
DUNE_SIMD_LOOP_PREFIX_OP(++)
auto operator!() const
Definition: loop.hh:120
DUNE_SIMD_LOOP_POSTFIX_OP(--)
DUNE_SIMD_LOOP_ASSIGNMENT_OP * DUNE_SIMD_LOOP_ASSIGNMENT_OP(/=);DUNE_SIMD_LOOP_ASSIGNMENT_OP(%=
DUNE_SIMD_LOOP_ASSIGNMENT_OP & DUNE_SIMD_LOOP_ASSIGNMENT_OP(|=);DUNE_SIMD_LOOP_ASSIGNMENT_OP(^=
DUNE_SIMD_LOOP_UNARY_OP(-)
DUNE_SIMD_LOOP_PREFIX_OP(--)
DUNE_SIMD_LOOP_POSTFIX_OP(++)
DUNE_SIMD_LOOP_ASSIGNMENT_OP(-=)
LoopSIMD(const LoopSIMD< T, S, OA > &other)
Definition: loop.hh:79
LoopSIMD()
Definition: loop.hh:69
DUNE_SIMD_LOOP_UNARY_OP(~)
DUNE_SIMD_LOOP_ASSIGNMENT_OP(+=)
DUNE_SIMD_LOOP_UNARY_OP(+)
DUNE_SIMD_LOOP_ASSIGNMENT_OP(<<=)
DUNE_SIMD_LOOP_ASSIGNMENT_OP(> >=)
Simd::Scalar< T > type
Definition: loop.hh:344
Helper class for tagging priorities.
Definition: typeutilities.hh:73
Include file for users of the SIMD abstraction layer.