blob: ce66c2c7868b3fe6e27135bea79b083380f97bfc [file] [log] [blame]
Austin Schuh189376f2018-12-20 22:11:15 +11001// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2015-2016 Gael Guennebaud <gael.guennebaud@inria.fr>
5//
6// This Source Code Form is subject to the terms of the Mozilla
7// Public License v. 2.0. If a copy of the MPL was not distributed
8// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10// workaround issue between gcc >= 4.7 and cuda 5.5
11#if (defined __GNUC__) && (__GNUC__>4 || __GNUC_MINOR__>=7)
12 #undef _GLIBCXX_ATOMIC_BUILTINS
13 #undef _GLIBCXX_USE_INT128
14#endif
15
16#define EIGEN_TEST_NO_LONGDOUBLE
17#define EIGEN_TEST_NO_COMPLEX
18#define EIGEN_TEST_FUNC cuda_basic
19#define EIGEN_DEFAULT_DENSE_INDEX_TYPE int
20
21#include <math_constants.h>
22#include <cuda.h>
23#include "main.h"
24#include "cuda_common.h"
25
26// Check that dense modules can be properly parsed by nvcc
27#include <Eigen/Dense>
28
29// struct Foo{
30// EIGEN_DEVICE_FUNC
31// void operator()(int i, const float* mats, float* vecs) const {
32// using namespace Eigen;
33// // Matrix3f M(data);
34// // Vector3f x(data+9);
35// // Map<Vector3f>(data+9) = M.inverse() * x;
36// Matrix3f M(mats+i/16);
37// Vector3f x(vecs+i*3);
38// // using std::min;
39// // using std::sqrt;
40// Map<Vector3f>(vecs+i*3) << x.minCoeff(), 1, 2;// / x.dot(x);//(M.inverse() * x) / x.x();
41// //x = x*2 + x.y() * x + x * x.maxCoeff() - x / x.sum();
42// }
43// };
44
45template<typename T>
46struct coeff_wise {
47 EIGEN_DEVICE_FUNC
48 void operator()(int i, const typename T::Scalar* in, typename T::Scalar* out) const
49 {
50 using namespace Eigen;
51 T x1(in+i);
52 T x2(in+i+1);
53 T x3(in+i+2);
54 Map<T> res(out+i*T::MaxSizeAtCompileTime);
55
56 res.array() += (in[0] * x1 + x2).array() * x3.array();
57 }
58};
59
60template<typename T>
61struct replicate {
62 EIGEN_DEVICE_FUNC
63 void operator()(int i, const typename T::Scalar* in, typename T::Scalar* out) const
64 {
65 using namespace Eigen;
66 T x1(in+i);
67 int step = x1.size() * 4;
68 int stride = 3 * step;
69
70 typedef Map<Array<typename T::Scalar,Dynamic,Dynamic> > MapType;
71 MapType(out+i*stride+0*step, x1.rows()*2, x1.cols()*2) = x1.replicate(2,2);
72 MapType(out+i*stride+1*step, x1.rows()*3, x1.cols()) = in[i] * x1.colwise().replicate(3);
73 MapType(out+i*stride+2*step, x1.rows(), x1.cols()*3) = in[i] * x1.rowwise().replicate(3);
74 }
75};
76
77template<typename T>
78struct redux {
79 EIGEN_DEVICE_FUNC
80 void operator()(int i, const typename T::Scalar* in, typename T::Scalar* out) const
81 {
82 using namespace Eigen;
83 int N = 10;
84 T x1(in+i);
85 out[i*N+0] = x1.minCoeff();
86 out[i*N+1] = x1.maxCoeff();
87 out[i*N+2] = x1.sum();
88 out[i*N+3] = x1.prod();
89 out[i*N+4] = x1.matrix().squaredNorm();
90 out[i*N+5] = x1.matrix().norm();
91 out[i*N+6] = x1.colwise().sum().maxCoeff();
92 out[i*N+7] = x1.rowwise().maxCoeff().sum();
93 out[i*N+8] = x1.matrix().colwise().squaredNorm().sum();
94 }
95};
96
97template<typename T1, typename T2>
98struct prod_test {
99 EIGEN_DEVICE_FUNC
100 void operator()(int i, const typename T1::Scalar* in, typename T1::Scalar* out) const
101 {
102 using namespace Eigen;
103 typedef Matrix<typename T1::Scalar, T1::RowsAtCompileTime, T2::ColsAtCompileTime> T3;
104 T1 x1(in+i);
105 T2 x2(in+i+1);
106 Map<T3> res(out+i*T3::MaxSizeAtCompileTime);
107 res += in[i] * x1 * x2;
108 }
109};
110
111template<typename T1, typename T2>
112struct diagonal {
113 EIGEN_DEVICE_FUNC
114 void operator()(int i, const typename T1::Scalar* in, typename T1::Scalar* out) const
115 {
116 using namespace Eigen;
117 T1 x1(in+i);
118 Map<T2> res(out+i*T2::MaxSizeAtCompileTime);
119 res += x1.diagonal();
120 }
121};
122
123template<typename T>
124struct eigenvalues {
125 EIGEN_DEVICE_FUNC
126 void operator()(int i, const typename T::Scalar* in, typename T::Scalar* out) const
127 {
128 using namespace Eigen;
129 typedef Matrix<typename T::Scalar, T::RowsAtCompileTime, 1> Vec;
130 T M(in+i);
131 Map<Vec> res(out+i*Vec::MaxSizeAtCompileTime);
132 T A = M*M.adjoint();
133 SelfAdjointEigenSolver<T> eig;
134 eig.computeDirect(M);
135 res = eig.eigenvalues();
136 }
137};
138
139void test_cuda_basic()
140{
141 ei_test_init_cuda();
142
143 int nthreads = 100;
144 Eigen::VectorXf in, out;
145
146 #ifndef __CUDA_ARCH__
147 int data_size = nthreads * 512;
148 in.setRandom(data_size);
149 out.setRandom(data_size);
150 #endif
151
152 CALL_SUBTEST( run_and_compare_to_cuda(coeff_wise<Vector3f>(), nthreads, in, out) );
153 CALL_SUBTEST( run_and_compare_to_cuda(coeff_wise<Array44f>(), nthreads, in, out) );
154
155 CALL_SUBTEST( run_and_compare_to_cuda(replicate<Array4f>(), nthreads, in, out) );
156 CALL_SUBTEST( run_and_compare_to_cuda(replicate<Array33f>(), nthreads, in, out) );
157
158 CALL_SUBTEST( run_and_compare_to_cuda(redux<Array4f>(), nthreads, in, out) );
159 CALL_SUBTEST( run_and_compare_to_cuda(redux<Matrix3f>(), nthreads, in, out) );
160
161 CALL_SUBTEST( run_and_compare_to_cuda(prod_test<Matrix3f,Matrix3f>(), nthreads, in, out) );
162 CALL_SUBTEST( run_and_compare_to_cuda(prod_test<Matrix4f,Vector4f>(), nthreads, in, out) );
163
164 CALL_SUBTEST( run_and_compare_to_cuda(diagonal<Matrix3f,Vector3f>(), nthreads, in, out) );
165 CALL_SUBTEST( run_and_compare_to_cuda(diagonal<Matrix4f,Vector4f>(), nthreads, in, out) );
166
167 CALL_SUBTEST( run_and_compare_to_cuda(eigenvalues<Matrix3f>(), nthreads, in, out) );
168 CALL_SUBTEST( run_and_compare_to_cuda(eigenvalues<Matrix2f>(), nthreads, in, out) );
169
170}