LvArray
math.hpp
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2021, Lawrence Livermore National Security, LLC and LvArray contributors.
3  * All rights reserved.
4  * See the LICENSE file for details.
5  * SPDX-License-Identifier: (BSD-3-Clause)
6  */
7 
13 #pragma once
14 
15 // Source includes
16 #include "LvArrayConfig.hpp"
17 #include "Macros.hpp"
18 
19 // System includes
20 #include <cmath>
21 #include <type_traits>
22 
23 #if defined( LVARRAY_USE_CUDA )
24  #include <cuda_fp16.h>
25 #endif
26 
27 namespace LvArray
28 {
29 
33 namespace math
34 {
35 
36 namespace internal
37 {
38 
47 template< typename T, typename U >
48 LVARRAY_HOST_DEVICE inline constexpr
49 T convert( T const, U const u )
50 { return u; }
51 
57 template< typename T >
58 LVARRAY_HOST_DEVICE inline constexpr
59 int numValues( T const )
60 { return 1; }
61 
66 template< typename T >
67 struct SingleType
68 {
70  using type = T;
71 };
72 
78 template< typename T >
79 LVARRAY_HOST_DEVICE inline constexpr
81 { return x; }
82 
88 template< typename T >
89 LVARRAY_HOST_DEVICE inline constexpr
91 { return x; }
92 
99 template< typename T >
100 LVARRAY_HOST_DEVICE inline constexpr
101 T lessThan( T const x, T const y )
102 { return __hlt( x, y ); }
103 
104 #if defined( LVARRAY_USE_CUDA )
105 
112 template< typename U >
113 LVARRAY_HOST_DEVICE inline constexpr
114 __half convert( __half const, U const u )
115 { return __float2half_rn( u ); }
116 
124 template< typename U >
125 LVARRAY_HOST_DEVICE inline constexpr
126 __half2 convert( __half2 const, U const u )
127 { return __float2half2_rn( u ); }
128 
134 LVARRAY_HOST_DEVICE inline
135 __half2 convert( __half2 const, __half const u )
136 {
137 #if defined( __CUDA_ARCH__ )
138  return __half2half2( u );
139 #else
140  return __float2half2_rn( u );
141 #endif
142 }
143 
153 template< typename U, typename V >
154 LVARRAY_HOST_DEVICE inline constexpr
155 __half2 convert( __half2 const, U const u, V const v )
156 { return __floats2half2_rn( u, v ); }
157 
164 LVARRAY_HOST_DEVICE inline
165 __half2 convert( __half2 const, __half const u, __half const v )
166 {
167 #if defined( __CUDA_ARCH__ )
168  return __halves2half2( u, v );
169 #else
170  return __floats2half2_rn( u, v );
171 #endif
172 }
173 
178 LVARRAY_HOST_DEVICE inline constexpr
179 int numValues( __half2 const & )
180 { return 2; }
181 
185 template<>
186 struct SingleType< __half2 >
187 {
189  using type = __half;
190 };
191 
196 LVARRAY_DEVICE inline
197 __half getFirst( __half2 const x )
198 { return __low2half( x ); }
199 
204 LVARRAY_DEVICE inline
205 __half getSecond( __half2 const x )
206 { return __high2half( x ); }
207 
213 LVARRAY_DEVICE inline
214 __half lessThan( __half const x, __half const y )
215 { return __hlt( x, y ); }
216 
222 LVARRAY_DEVICE inline
223 __half2 lessThan( __half2 const x, __half2 const y )
224 { return __hlt2( x, y ); }
225 
226 #endif
227 
228 } // namespace internal
229 
233 
240 template< typename T >
241 LVARRAY_HOST_DEVICE inline constexpr
243 { return internal::numValues( T() ); }
244 
249 template< typename T >
251 
260 template< typename T, typename U >
261 LVARRAY_HOST_DEVICE inline constexpr
262 T convert( U const u )
263 { return internal::convert( T(), u ); }
264 
275 template< typename T, typename U, typename V >
276 LVARRAY_HOST_DEVICE inline constexpr
277 T convert( U const u, V const v )
278 { return internal::convert( T(), u, v ); }
279 
286 template< typename T >
287 LVARRAY_DEVICE inline
289 { return internal::getFirst( x ); }
290 
297 template< typename T >
298 LVARRAY_DEVICE inline
300 { return internal::getSecond( x ); }
301 
308 template< typename T >
309 LVARRAY_HOST_DEVICE inline constexpr
310 std::enable_if_t< std::is_arithmetic< T >::value, T >
311 max( T const a, T const b )
312 {
313 #if defined(__CUDA_ARCH__)
314  return ::max( a, b );
315 #else
316  return std::max( a, b );
317 #endif
318 }
319 
320 #if defined( LVARRAY_USE_CUDA )
321 
323 LVARRAY_DEVICE inline
324 __half max( __half const a, __half const b )
325 {
326 #if CUDART_VERSION > 11000 && (__CUDA_ARCH__ >= 800 || !defined(__CUDA_ARCH__))
327  return __hmax( a, b );
328 #else
329  return a > b ? a : b;
330 #endif
331 }
332 
334 LVARRAY_DEVICE inline
335 __half2 max( __half2 const a, __half2 const b )
336 {
337 #if CUDART_VERSION > 11000 && (__CUDA_ARCH__ >= 800 || !defined(__CUDA_ARCH__))
338  return __hmax2( a, b );
339 #else
340  __half2 const aFactor = __hge2( a, b );
341  __half2 const bFactor = convert< __half2 >( 1 ) - aFactor;
342  return a * aFactor + bFactor * b;
343 #endif
344 }
345 
346 #endif
347 
348 
355 template< typename T >
356 LVARRAY_HOST_DEVICE inline constexpr
357 std::enable_if_t< std::is_arithmetic< T >::value, T >
358 min( T const a, T const b )
359 {
360 #if defined(__CUDA_ARCH__)
361  return ::min( a, b );
362 #else
363  return std::min( a, b );
364 #endif
365 }
366 
367 #if defined( LVARRAY_USE_CUDA )
368 
370 LVARRAY_DEVICE inline
371 __half min( __half const a, __half const b )
372 {
373 #if CUDART_VERSION > 11000 && (__CUDA_ARCH__ >= 800 || !defined(__CUDA_ARCH__))
374  return __hmin( a, b );
375 #else
376  return a < b ? a : b;
377 #endif
378 }
379 
381 LVARRAY_DEVICE inline
382 __half2 min( __half2 const a, __half2 const b )
383 {
384 #if CUDART_VERSION > 11000 && (__CUDA_ARCH__ >= 800 || !defined(__CUDA_ARCH__))
385  return __hmin2( a, b );
386 #else
387  __half2 const aFactor = __hle2( a, b );
388  __half2 const bFactor = convert< __half2 >( 1 ) - aFactor;
389  return a * aFactor + bFactor * b;
390 #endif
391 }
392 
393 #endif
394 
400 template< typename T >
401 LVARRAY_HOST_DEVICE inline constexpr
402 T abs( T const x )
403 {
404 #if defined(__CUDA_ARCH__)
405  return ::abs( x );
406 #else
407  return std::abs( x );
408 #endif
409 }
410 
411 #if defined( LVARRAY_USE_CUDA )
412 
414 LVARRAY_DEVICE inline
415 __half abs( __half const x )
416 {
417 #if CUDART_VERSION > 11000
418  return __habs( x );
419 #else
420  return x > __half( 0 ) ? x : -x;
421 #endif
422 }
423 
425 LVARRAY_DEVICE inline
426 __half2 abs( __half2 const x )
427 {
428 #if CUDART_VERSION > 11000
429  return __habs2( x );
430 #else
431  return x - __hle2( x, convert< __half2 >( 0 ) ) * ( x + x );
432 #endif
433 }
434 
435 #endif
436 
442 template< typename T >
443 LVARRAY_HOST_DEVICE inline constexpr
444 T square( T const x )
445 { return x * x; }
446 
448 
452 
460 LVARRAY_HOST_DEVICE inline
461 float sqrt( float const x )
462 {
463 #if defined(__CUDA_ARCH__)
464  return ::sqrtf( x );
465 #else
466  return std::sqrt( x );
467 #endif
468 }
469 
471 template< typename T >
472 LVARRAY_HOST_DEVICE inline
473 double sqrt( T const x )
474 {
475 #if defined(__CUDA_ARCH__)
476  return ::sqrt( double( x ) );
477 #else
478  return std::sqrt( x );
479 #endif
480 }
481 
482 #if defined( LVARRAY_USE_CUDA )
483 
485 LVARRAY_DEVICE inline
486 __half sqrt( __half const x )
487 { return ::hsqrt( x ); }
488 
490 LVARRAY_DEVICE inline
491 __half2 sqrt( __half2 const x )
492 { return ::h2sqrt( x ); }
493 
494 #endif
495 
502 LVARRAY_HOST_DEVICE inline
503 float invSqrt( float const x )
504 {
505 #if defined(__CUDA_ARCH__)
506  return ::rsqrtf( x );
507 #else
508  return 1 / std::sqrt( x );
509 #endif
510 }
511 
513 template< typename T >
514 LVARRAY_HOST_DEVICE inline
515 double invSqrt( T const x )
516 {
517 #if defined( __CUDA_ARCH__ )
518  return ::rsqrt( double( x ) );
519 #else
520  return 1 / std::sqrt( x );
521 #endif
522 }
523 
524 #if defined( LVARRAY_USE_CUDA )
525 
527 LVARRAY_DEVICE inline
528 __half invSqrt( __half const x )
529 { return ::hrsqrt( x ); }
530 
532 LVARRAY_DEVICE inline
533 __half2 invSqrt( __half2 const x )
534 { return ::h2rsqrt( x ); }
535 
536 #endif
537 
539 
543 
551 LVARRAY_HOST_DEVICE inline
552 float sin( float const theta )
553 {
554 #if defined(__CUDA_ARCH__)
555  return ::sinf( theta );
556 #else
557  return std::sin( theta );
558 #endif
559 }
560 
562 template< typename T >
563 LVARRAY_HOST_DEVICE inline
564 double sin( T const theta )
565 {
566 #if defined(__CUDA_ARCH__)
567  return ::sin( double( theta ) );
568 #else
569  return std::sin( theta );
570 #endif
571 }
572 
573 #if defined( LVARRAY_USE_CUDA )
574 
576 LVARRAY_DEVICE inline
577 __half sin( __half const theta )
578 { return ::hsin( theta ); }
579 
581 LVARRAY_DEVICE inline
582 __half2 sin( __half2 const theta )
583 { return ::h2sin( theta ); }
584 
585 #endif
586 
593 LVARRAY_HOST_DEVICE inline
594 float cos( float const theta )
595 {
596 #if defined(__CUDA_ARCH__)
597  return ::cosf( theta );
598 #else
599  return std::cos( theta );
600 #endif
601 }
602 
604 template< typename T >
605 LVARRAY_HOST_DEVICE inline
606 double cos( T const theta )
607 {
608 #if defined(__CUDA_ARCH__)
609  return ::cos( double( theta ) );
610 #else
611  return std::cos( theta );
612 #endif
613 }
614 
615 #if defined( LVARRAY_USE_CUDA )
616 
618 LVARRAY_DEVICE inline
619 __half cos( __half const theta )
620 { return ::hcos( theta ); }
621 
623 LVARRAY_DEVICE inline
624 __half2 cos( __half2 const theta )
625 { return ::h2cos( theta ); }
626 
627 #endif
628 
635 LVARRAY_HOST_DEVICE inline
636 void sincos( float const theta, float & sinTheta, float & cosTheta )
637 {
638 #if defined(__CUDA_ARCH__)
639  ::sincos( theta, &sinTheta, &cosTheta );
640 #else
641  sinTheta = std::sin( theta );
642  cosTheta = std::cos( theta );
643 #endif
644 }
645 
647 template< typename T >
648 LVARRAY_HOST_DEVICE inline
649 void sincos( double const theta, double & sinTheta, double & cosTheta )
650 {
651 #if defined(__CUDA_ARCH__)
652  ::sincos( theta, &sinTheta, &cosTheta );
653 #else
654  sinTheta = std::sin( theta );
655  cosTheta = std::cos( theta );
656 #endif
657 }
658 
660 template< typename T >
661 LVARRAY_HOST_DEVICE inline
662 void sincos( T const theta, double & sinTheta, double & cosTheta )
663 {
664 #if defined(__CUDA_ARCH__)
665  double s, c;
666  ::sincos( theta, &s, &c );
667  sinTheta = s;
668  cosTheta = c;
669 #else
670  sinTheta = std::sin( theta );
671  cosTheta = std::cos( theta );
672 #endif
673 }
674 
675 #if defined( LVARRAY_USE_CUDA )
676 
678 LVARRAY_DEVICE inline
679 void sincos( __half const theta, __half & sinTheta, __half & cosTheta )
680 {
681  sinTheta = ::hsin( theta );
682  cosTheta = ::hcos( theta );
683 }
684 
686 LVARRAY_DEVICE inline
687 void sincos( __half2 const theta, __half2 & sinTheta, __half2 & cosTheta )
688 {
689  sinTheta = ::h2sin( theta );
690  cosTheta = ::h2cos( theta );
691 }
692 
693 #endif
694 
701 LVARRAY_HOST_DEVICE inline
702 float tan( float const theta )
703 {
704 #if defined(__CUDA_ARCH__)
705  return ::tanf( theta );
706 #else
707  return std::tan( theta );
708 #endif
709 }
710 
712 template< typename T >
713 LVARRAY_HOST_DEVICE inline
714 double tan( T const theta )
715 {
716 #if defined(__CUDA_ARCH__)
717  return ::tan( double( theta ) );
718 #else
719  return std::tan( theta );
720 #endif
721 }
722 
723 #if defined( LVARRAY_USE_CUDA )
724 
726 LVARRAY_DEVICE inline
727 __half tan( __half const theta )
728 {
729  __half s, c;
730  sincos( theta, s, c );
731  return s / c;
732 }
733 
735 LVARRAY_DEVICE inline
736 __half2 tan( __half2 const theta )
737 {
738  __half2 s, c;
739  sincos( theta, s, c );
740  return s / c;
741 }
742 
743 #endif
744 
746 
750 
752 namespace internal
753 {
754 
766 template< typename T >
767 LVARRAY_DEVICE inline
768 T asinImpl( T const x )
769 {
770  T const negate = lessThan( x, math::convert< T >( 0 ) );
771  T const absX = abs( x );
772 
773  T ret = math::convert< T >( -0.0187293 ) * absX + math::convert< T >( 0.0742610 );
774  ret = ret * absX - math::convert< T >( 0.2121144 );
775  ret = ret * absX + math::convert< T >( 1.5707288 );
776  ret = math::convert< T >( 3.14159265358979 * 0.5 ) - ret * sqrt( math::convert< T >( 1 ) - absX );
777  ret = ret - negate * ( ret + ret );
778  T const smallAngle = lessThan( absX, math::convert< T >( 1.7e-1 ) );
779  return smallAngle * x + ( math::convert< T >( 1 ) - smallAngle ) * ret;
780 }
781 
788 template< typename T >
789 LVARRAY_DEVICE inline
790 T acosImpl( T const x )
791 {
792  T const negate = lessThan( x, math::convert< T >( 0 ) );
793  T const absX = abs( x );
794 
795  T ret = math::convert< T >( -0.0187293 ) * absX + math::convert< T >( 0.0742610 );
796  ret = ret * absX - math::convert< T >( 0.2121144 );
797  ret = ret * absX + math::convert< T >( 1.5707288 );
798  ret = ret * sqrt( math::convert< T >( 1 ) - absX );
799  ret = ret - negate * ( ret + ret );
800  return negate * math::convert< T >( 3.14159265358979 ) + ret;
801 }
802 
810 template< typename T >
811 LVARRAY_DEVICE inline
812 T atan2Impl( T const y, T const x )
813 {
814  T const absX = abs( x );
815  T const absY = abs( y );
816  T const ratio = min( absX, absY ) / max( absX, absY );
817  T const ratio2 = ratio * ratio;
818 
819  T ret = math::convert< T >( -0.013480470 ) * ratio2 + math::convert< T >( 0.057477314 );
820  ret = ret * ratio2 - math::convert< T >( 0.121239071 );
821  ret = ret * ratio2 + math::convert< T >( 0.195635925 );
822  ret = ret * ratio2 - math::convert< T >( 0.332994597 );
823  ret = ret * ratio2 + math::convert< T >( 0.999995630 );
824  ret = ret * ratio;
825 
826  // For single values the following works out to:
827  // ret = absX < absY ? 1.570796327 - ret : ret;
828  // ret = x < 0 ? 3.141592654 - ret : ret;
829  // ret = y < 0 ? -ret : ret;
830  ret = internal::lessThan( absX, absY ) * ( math::convert< T >( 1.570796327 ) - ret - ret ) + ret;
831  ret = internal::lessThan( x, math::convert< T >( 0 ) ) * ( math::convert< T >( 3.141592654 ) - ret - ret ) + ret;
832  ret = ret - internal::lessThan( y, math::convert< T >( 0 ) ) * ( ret + ret );
833 
834  return ret;
835 }
836 
837 } // namespace internal
838 
845 LVARRAY_HOST_DEVICE inline
846 float asin( float const x )
847 {
848 #if defined(__CUDA_ARCH__)
849  return ::asinf( x );
850 #else
851  return std::asin( x );
852 #endif
853 }
854 
856 template< typename T >
857 LVARRAY_HOST_DEVICE inline
858 double asin( T const x )
859 {
860 #if defined(__CUDA_ARCH__)
861  return ::asin( double( x ) );
862 #else
863  return std::asin( x );
864 #endif
865 }
866 
867 #if defined( LVARRAY_USE_CUDA )
868 
870 LVARRAY_DEVICE inline
871 __half asin( __half const x )
872 { return internal::asinImpl( x ); }
873 
875 LVARRAY_DEVICE inline
876 __half2 asin( __half2 const x )
877 { return internal::asinImpl( x ); }
878 
879 #endif
880 
887 LVARRAY_HOST_DEVICE inline
888 float acos( float const x )
889 {
890 #if defined(__CUDA_ARCH__)
891  return ::acosf( x );
892 #else
893  return std::acos( x );
894 #endif
895 }
896 
898 template< typename T >
899 LVARRAY_HOST_DEVICE inline
900 double acos( T const x )
901 {
902 #if defined(__CUDA_ARCH__)
903  return ::acos( double( x ) );
904 #else
905  return std::acos( x );
906 #endif
907 }
908 
909 #if defined( LVARRAY_USE_CUDA )
910 
912 LVARRAY_DEVICE inline
913 __half acos( __half const x )
914 { return internal::acosImpl( x ); }
915 
917 LVARRAY_DEVICE inline
918 __half2 acos( __half2 const x )
919 { return internal::acosImpl( x ); }
920 
921 #endif
922 
930 LVARRAY_HOST_DEVICE inline
931 float atan2( float const y, float const x )
932 {
933 #if defined(__CUDA_ARCH__)
934  return ::atan2f( y, x );
935 #else
936  return std::atan2( y, x );
937 #endif
938 }
939 
941 template< typename T >
942 LVARRAY_HOST_DEVICE inline
943 double atan2( T const y, T const x )
944 {
945 #if defined(__CUDA_ARCH__)
946  return ::atan2( double( y ), double( x ) );
947 #else
948  return std::atan2( y, x );
949 #endif
950 }
951 
952 #if defined( LVARRAY_USE_CUDA )
953 
955 LVARRAY_DEVICE inline
956 __half atan2( __half const y, __half const x )
957 { return internal::atan2Impl( y, x ); }
958 
960 LVARRAY_DEVICE inline
961 __half2 atan2( __half2 const y, __half2 const x )
962 { return internal::atan2Impl( y, x ); }
963 
964 #endif
965 
967 
971 
979 LVARRAY_HOST_DEVICE inline
980 float exp( float const x )
981 {
982 #if defined(__CUDA_ARCH__)
983  return ::expf( x );
984 #else
985  return std::exp( x );
986 #endif
987 }
988 
990 template< typename T >
991 LVARRAY_HOST_DEVICE inline
992 double exp( T const x )
993 {
994 #if defined(__CUDA_ARCH__)
995  return ::exp( double( x ) );
996 #else
997  return std::exp( x );
998 #endif
999 }
1000 
1001 #if defined( LVARRAY_USE_CUDA )
1002 
1004 LVARRAY_DEVICE inline
1005 __half exp( __half const x )
1006 { return ::hexp( x ); }
1007 
1009 LVARRAY_DEVICE inline
1010 __half2 exp( __half2 const x )
1011 { return ::h2exp( x ); }
1012 
1013 #endif
1014 
1021 LVARRAY_HOST_DEVICE inline
1022 float log( float const x )
1023 {
1024 #if defined(__CUDA_ARCH__)
1025  return ::logf( x );
1026 #else
1027  return std::log( x );
1028 #endif
1029 }
1030 
1032 template< typename T >
1033 LVARRAY_HOST_DEVICE inline
1034 double log( T const x )
1035 {
1036 #if defined(__CUDA_ARCH__)
1037  return ::log( double( x ) );
1038 #else
1039  return std::log( x );
1040 #endif
1041 }
1042 
1043 #if defined( LVARRAY_USE_CUDA )
1044 
1046 LVARRAY_DEVICE inline
1047 __half log( __half const x )
1048 { return ::hlog( x ); }
1049 
1051 LVARRAY_DEVICE inline
1052 __half2 log( __half2 const x )
1053 { return ::h2log( x ); }
1054 
1055 #endif
1056 
1058 
1059 } // namespace math
1060 } // namespace LvArray
LVARRAY_HOST_DEVICE double exp(T const x)
Definition: math.hpp:992
LVARRAY_HOST_DEVICE constexpr T convert(U const u)
Convert u to type.
Definition: math.hpp:262
LVARRAY_HOST_DEVICE double tan(T const theta)
Definition: math.hpp:714
LVARRAY_HOST_DEVICE float asin(float const x)
Definition: math.hpp:846
LVARRAY_HOST_DEVICE double sqrt(T const x)
Definition: math.hpp:473
LVARRAY_HOST_DEVICE float acos(float const x)
Definition: math.hpp:888
LVARRAY_DEVICE T acosImpl(T const x)
Definition: math.hpp:790
LVARRAY_HOST_DEVICE double acos(T const x)
Definition: math.hpp:900
LVARRAY_HOST_DEVICE float cos(float const theta)
Definition: math.hpp:594
LVARRAY_HOST_DEVICE constexpr std::enable_if_t< std::is_arithmetic< T >::value, T > min(T const a, T const b)
Definition: math.hpp:358
LVARRAY_DEVICE T asinImpl(T const x)
Definition: math.hpp:768
LVARRAY_HOST_DEVICE constexpr T convert(U const u, V const v)
Convert u and v to a dual type.
Definition: math.hpp:277
The type of a single value of type T.
Definition: math.hpp:67
LVARRAY_HOST_DEVICE float log(float const x)
Definition: math.hpp:1022
LVARRAY_HOST_DEVICE float sin(float const theta)
Definition: math.hpp:552
LVARRAY_HOST_DEVICE void sincos(float const theta, float &sinTheta, float &cosTheta)
Compute the sine and cosine of theta.
Definition: math.hpp:636
LVARRAY_HOST_DEVICE double asin(T const x)
Definition: math.hpp:858
T type
An alias for T.
Definition: math.hpp:70
LVARRAY_HOST_DEVICE double cos(T const theta)
Definition: math.hpp:606
#define LVARRAY_DEVICE
Mark a function for only device usage.
Definition: Macros.hpp:552
LVARRAY_HOST_DEVICE constexpr T lessThan(T const x, T const y)
Definition: math.hpp:101
LVARRAY_DEVICE SingleType< T > getFirst(T const x)
Definition: math.hpp:288
LVARRAY_HOST_DEVICE float sqrt(float const x)
Definition: math.hpp:461
The top level namespace.
Definition: Array.hpp:24
LVARRAY_DEVICE SingleType< T > getSecond(T const x)
Definition: math.hpp:299
Contains a bunch of macro definitions.
LVARRAY_HOST_DEVICE constexpr T abs(T const x)
Definition: math.hpp:402
LVARRAY_HOST_DEVICE float atan2(float const y, float const x)
Definition: math.hpp:931
LVARRAY_DEVICE T atan2Impl(T const y, T const x)
Definition: math.hpp:812
LVARRAY_HOST_DEVICE float exp(float const x)
Definition: math.hpp:980
LVARRAY_HOST_DEVICE float tan(float const theta)
Definition: math.hpp:702
LVARRAY_HOST_DEVICE constexpr std::enable_if_t< std::is_arithmetic< T >::value, T > max(T const a, T const b)
Definition: math.hpp:311
LVARRAY_HOST_DEVICE double atan2(T const y, T const x)
Definition: math.hpp:943
typename internal::SingleType< T >::type SingleType
The type of a single value of type T.
Definition: math.hpp:250
LVARRAY_HOST_DEVICE constexpr int numValues()
Return the number of values stored in type.
Definition: math.hpp:242
LVARRAY_HOST_DEVICE double log(T const x)
Definition: math.hpp:1034
LVARRAY_HOST_DEVICE float invSqrt(float const x)
Definition: math.hpp:503
LVARRAY_HOST_DEVICE double sin(T const theta)
Definition: math.hpp:564
#define LVARRAY_HOST_DEVICE
Mark a function for both host and device usage.
Definition: Macros.hpp:549
LVARRAY_HOST_DEVICE constexpr T square(T const x)
Definition: math.hpp:444