blob: 7e41e9e2077cb3020d8b4ace5e7b1a535898e73a [file] [log] [blame]
Gael Guennebaud8cef5412008-06-21 17:28:07 +00001// This file is part of Eigen, a lightweight C++ template library
Benoit Jacob6347b1d2009-05-22 20:25:33 +02002// for linear algebra.
Gael Guennebaud8cef5412008-06-21 17:28:07 +00003//
Gael Guennebaud28e64b02010-06-24 23:21:58 +02004// Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud@inria.fr>
Gael Guennebaud8cef5412008-06-21 17:28:07 +00005//
Benoit Jacob69124cf2012-07-13 14:42:47 -04006// 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/.
Gael Guennebaud19bc4182013-02-25 01:17:44 +01009
Gael Guennebaud8cef5412008-06-21 17:28:07 +000010#include "main.h"
Gael Guennebaud8cef5412008-06-21 17:28:07 +000011
Benoit Jacobe2775862010-04-28 18:51:38 -040012template<typename ArrayType> void array(const ArrayType& m)
Gael Guennebaud8cef5412008-06-21 17:28:07 +000013{
Hauke Heibelf1679c72010-06-20 17:37:56 +020014 typedef typename ArrayType::Index Index;
Benoit Jacobe2775862010-04-28 18:51:38 -040015 typedef typename ArrayType::Scalar Scalar;
Gael Guennebaudd476cad2016-06-25 10:12:06 +020016 typedef typename ArrayType::RealScalar RealScalar;
Benoit Jacobe2775862010-04-28 18:51:38 -040017 typedef Array<Scalar, ArrayType::RowsAtCompileTime, 1> ColVectorType;
18 typedef Array<Scalar, 1, ArrayType::ColsAtCompileTime> RowVectorType;
Gael Guennebaud8cef5412008-06-21 17:28:07 +000019
Hauke Heibelf1679c72010-06-20 17:37:56 +020020 Index rows = m.rows();
Hauke Heibelf578dc72010-12-16 17:34:13 +010021 Index cols = m.cols();
Gael Guennebaud8cef5412008-06-21 17:28:07 +000022
Benoit Jacobe2775862010-04-28 18:51:38 -040023 ArrayType m1 = ArrayType::Random(rows, cols),
24 m2 = ArrayType::Random(rows, cols),
Gael Guennebaud8cef5412008-06-21 17:28:07 +000025 m3(rows, cols);
Christoph Hertzberg3be9f5c2015-04-16 13:25:20 +020026 ArrayType m4 = m1; // copy constructor
27 VERIFY_IS_APPROX(m1, m4);
Gael Guennebaud8cef5412008-06-21 17:28:07 +000028
Gael Guennebaud627595a2009-06-10 11:20:30 +020029 ColVectorType cv1 = ColVectorType::Random(rows);
30 RowVectorType rv1 = RowVectorType::Random(cols);
31
Benoit Jacob47160402010-10-25 10:15:22 -040032 Scalar s1 = internal::random<Scalar>(),
Hauke Heibel8cb3e362012-03-02 16:27:27 +010033 s2 = internal::random<Scalar>();
34
Gael Guennebaud59dc1da2008-09-03 17:16:28 +000035 // scalar addition
Gael Guennebaudc70d5422010-01-18 22:54:20 +010036 VERIFY_IS_APPROX(m1 + s1, s1 + m1);
Benoit Jacobe2775862010-04-28 18:51:38 -040037 VERIFY_IS_APPROX(m1 + s1, ArrayType::Constant(rows,cols,s1) + m1);
Gael Guennebaudc70d5422010-01-18 22:54:20 +010038 VERIFY_IS_APPROX(s1 - m1, (-m1)+s1 );
Benoit Jacobe2775862010-04-28 18:51:38 -040039 VERIFY_IS_APPROX(m1 - s1, m1 - ArrayType::Constant(rows,cols,s1));
40 VERIFY_IS_APPROX(s1 - m1, ArrayType::Constant(rows,cols,s1) - m1);
41 VERIFY_IS_APPROX((m1*Scalar(2)) - s2, (m1+m1) - ArrayType::Constant(rows,cols,s2) );
Gael Guennebaud8cef5412008-06-21 17:28:07 +000042 m3 = m1;
Gael Guennebaudc70d5422010-01-18 22:54:20 +010043 m3 += s2;
44 VERIFY_IS_APPROX(m3, m1 + s2);
Gael Guennebaud8cef5412008-06-21 17:28:07 +000045 m3 = m1;
Gael Guennebaudc70d5422010-01-18 22:54:20 +010046 m3 -= s1;
Hauke Heibelf578dc72010-12-16 17:34:13 +010047 VERIFY_IS_APPROX(m3, m1 - s1);
48
49 // scalar operators via Maps
50 m3 = m1;
51 ArrayType::Map(m1.data(), m1.rows(), m1.cols()) -= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
52 VERIFY_IS_APPROX(m1, m3 - m2);
53
54 m3 = m1;
55 ArrayType::Map(m1.data(), m1.rows(), m1.cols()) += ArrayType::Map(m2.data(), m2.rows(), m2.cols());
56 VERIFY_IS_APPROX(m1, m3 + m2);
57
58 m3 = m1;
59 ArrayType::Map(m1.data(), m1.rows(), m1.cols()) *= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
60 VERIFY_IS_APPROX(m1, m3 * m2);
61
62 m3 = m1;
63 m2 = ArrayType::Random(rows,cols);
64 m2 = (m2==0).select(1,m2);
65 ArrayType::Map(m1.data(), m1.rows(), m1.cols()) /= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
66 VERIFY_IS_APPROX(m1, m3 / m2);
Gael Guennebaud05ad0832008-07-19 13:03:23 +000067
Gael Guennebaud59dc1da2008-09-03 17:16:28 +000068 // reductions
Gael Guennebaud42af5872013-02-23 22:58:14 +010069 VERIFY_IS_APPROX(m1.abs().colwise().sum().sum(), m1.abs().sum());
70 VERIFY_IS_APPROX(m1.abs().rowwise().sum().sum(), m1.abs().sum());
71 using std::abs;
72 VERIFY_IS_MUCH_SMALLER_THAN(abs(m1.colwise().sum().sum() - m1.sum()), m1.abs().sum());
73 VERIFY_IS_MUCH_SMALLER_THAN(abs(m1.rowwise().sum().sum() - m1.sum()), m1.abs().sum());
Gael Guennebaud28e139a2013-02-23 23:06:45 +010074 if (!internal::isMuchSmallerThan(abs(m1.sum() - (m1+m2).sum()), m1.abs().sum(), test_precision<Scalar>()))
Hauke Heibel83e3c452010-12-16 18:53:02 +010075 VERIFY_IS_NOT_APPROX(((m1+m2).rowwise().sum()).sum(), m1.sum());
Gael Guennebaud66e99ab2016-06-06 15:11:41 +020076 VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(internal::scalar_sum_op<Scalar,Scalar>()));
Gael Guennebaud627595a2009-06-10 11:20:30 +020077
78 // vector-wise ops
79 m3 = m1;
80 VERIFY_IS_APPROX(m3.colwise() += cv1, m1.colwise() + cv1);
81 m3 = m1;
82 VERIFY_IS_APPROX(m3.colwise() -= cv1, m1.colwise() - cv1);
83 m3 = m1;
84 VERIFY_IS_APPROX(m3.rowwise() += rv1, m1.rowwise() + rv1);
85 m3 = m1;
86 VERIFY_IS_APPROX(m3.rowwise() -= rv1, m1.rowwise() - rv1);
Gael Guennebaud0a18eec2014-09-19 13:25:28 +020087
88 // Conversion from scalar
89 VERIFY_IS_APPROX((m3 = s1), ArrayType::Constant(rows,cols,s1));
90 VERIFY_IS_APPROX((m3 = 1), ArrayType::Constant(rows,cols,1));
91 VERIFY_IS_APPROX((m3.topLeftCorner(rows,cols) = 1), ArrayType::Constant(rows,cols,1));
92 typedef Array<Scalar,
93 ArrayType::RowsAtCompileTime==Dynamic?2:ArrayType::RowsAtCompileTime,
94 ArrayType::ColsAtCompileTime==Dynamic?2:ArrayType::ColsAtCompileTime,
95 ArrayType::Options> FixedArrayType;
96 FixedArrayType f1(s1);
97 VERIFY_IS_APPROX(f1, FixedArrayType::Constant(s1));
98 FixedArrayType f2(numext::real(s1));
99 VERIFY_IS_APPROX(f2, FixedArrayType::Constant(numext::real(s1)));
100 FixedArrayType f3((int)100*numext::real(s1));
101 VERIFY_IS_APPROX(f3, FixedArrayType::Constant((int)100*numext::real(s1)));
102 f1.setRandom();
103 FixedArrayType f4(f1.data());
104 VERIFY_IS_APPROX(f4, f1);
105
Gael Guennebaudd476cad2016-06-25 10:12:06 +0200106 // pow
107 VERIFY_IS_APPROX(m1.pow(2), m1.square());
108 VERIFY_IS_APPROX(pow(m1,2), m1.square());
109 VERIFY_IS_APPROX(m1.pow(3), m1.cube());
110 VERIFY_IS_APPROX(pow(m1,3), m1.cube());
111 VERIFY_IS_APPROX((-m1).pow(3), -m1.cube());
112 VERIFY_IS_APPROX(pow(2*m1,3), 8*m1.cube());
113 ArrayType exponents = ArrayType::Constant(rows, cols, RealScalar(2));
114 VERIFY_IS_APPROX(Eigen::pow(m1,exponents), m1.square());
115 VERIFY_IS_APPROX(m1.pow(exponents), m1.square());
116 VERIFY_IS_APPROX(Eigen::pow(2*m1,exponents), 4*m1.square());
117 VERIFY_IS_APPROX((2*m1).pow(exponents), 4*m1.square());
118 VERIFY_IS_APPROX(Eigen::pow(m1,2*exponents), m1.square().square());
119 VERIFY_IS_APPROX(m1.pow(2*exponents), m1.square().square());
120 VERIFY_IS_APPROX(pow(m1(0,0), exponents), ArrayType::Constant(rows,cols,m1(0,0)*m1(0,0)));
121
Gael Guennebaud0a18eec2014-09-19 13:25:28 +0200122 // Check possible conflicts with 1D ctor
123 typedef Array<Scalar, Dynamic, 1> OneDArrayType;
124 OneDArrayType o1(rows);
125 VERIFY(o1.size()==rows);
126 OneDArrayType o4((int)rows);
127 VERIFY(o4.size()==rows);
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000128}
129
Benoit Jacobe2775862010-04-28 18:51:38 -0400130template<typename ArrayType> void comparisons(const ArrayType& m)
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000131{
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100132 using std::abs;
Hauke Heibelf1679c72010-06-20 17:37:56 +0200133 typedef typename ArrayType::Index Index;
Benoit Jacobe2775862010-04-28 18:51:38 -0400134 typedef typename ArrayType::Scalar Scalar;
Benoit Jacob874ff5a2009-01-26 17:56:04 +0000135 typedef typename NumTraits<Scalar>::Real RealScalar;
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000136
Hauke Heibelf1679c72010-06-20 17:37:56 +0200137 Index rows = m.rows();
138 Index cols = m.cols();
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000139
Benoit Jacob47160402010-10-25 10:15:22 -0400140 Index r = internal::random<Index>(0, rows-1),
141 c = internal::random<Index>(0, cols-1);
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000142
Benoit Jacobe2775862010-04-28 18:51:38 -0400143 ArrayType m1 = ArrayType::Random(rows, cols),
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200144 m2 = ArrayType::Random(rows, cols),
145 m3(rows, cols),
146 m4 = m1;
147
148 m4 = (m4.abs()==Scalar(0)).select(1,m4);
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000149
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100150 VERIFY(((m1 + Scalar(1)) > m1).all());
151 VERIFY(((m1 - Scalar(1)) < m1).all());
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000152 if (rows*cols>1)
153 {
154 m3 = m1;
155 m3(r,c) += 1;
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100156 VERIFY(! (m1 < m3).all() );
157 VERIFY(! (m1 > m3).all() );
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000158 }
Christoph Hertzberg494fa992015-05-07 17:28:40 +0200159 VERIFY(!(m1 > m2 && m1 < m2).any());
160 VERIFY((m1 <= m2 || m1 >= m2).all());
Gael Guennebaud627595a2009-06-10 11:20:30 +0200161
Christoph Hertzberg36052c42013-10-17 15:37:29 +0200162 // comparisons array to scalar
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100163 VERIFY( (m1 != (m1(r,c)+1) ).any() );
Christoph Hertzberg36052c42013-10-17 15:37:29 +0200164 VERIFY( (m1 > (m1(r,c)-1) ).any() );
165 VERIFY( (m1 < (m1(r,c)+1) ).any() );
166 VERIFY( (m1 == m1(r,c) ).any() );
167
168 // comparisons scalar to array
169 VERIFY( ( (m1(r,c)+1) != m1).any() );
170 VERIFY( ( (m1(r,c)-1) < m1).any() );
171 VERIFY( ( (m1(r,c)+1) > m1).any() );
172 VERIFY( ( m1(r,c) == m1).any() );
Gael Guennebaud627595a2009-06-10 11:20:30 +0200173
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000174 // test Select
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100175 VERIFY_IS_APPROX( (m1<m2).select(m1,m2), m1.cwiseMin(m2) );
176 VERIFY_IS_APPROX( (m1>m2).select(m1,m2), m1.cwiseMax(m2) );
Gael Guennebaud9f795582009-12-16 19:18:40 +0100177 Scalar mid = (m1.cwiseAbs().minCoeff() + m1.cwiseAbs().maxCoeff())/Scalar(2);
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000178 for (int j=0; j<cols; ++j)
179 for (int i=0; i<rows; ++i)
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100180 m3(i,j) = abs(m1(i,j))<mid ? 0 : m1(i,j);
Benoit Jacobe2775862010-04-28 18:51:38 -0400181 VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid))
182 .select(ArrayType::Zero(rows,cols),m1), m3);
Gael Guennebaude14aa8c2008-09-03 17:56:06 +0000183 // shorter versions:
Benoit Jacobe2775862010-04-28 18:51:38 -0400184 VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid))
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000185 .select(0,m1), m3);
Benoit Jacobe2775862010-04-28 18:51:38 -0400186 VERIFY_IS_APPROX( (m1.abs()>=ArrayType::Constant(rows,cols,mid))
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000187 .select(m1,0), m3);
Gael Guennebaude14aa8c2008-09-03 17:56:06 +0000188 // even shorter version:
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100189 VERIFY_IS_APPROX( (m1.abs()<mid).select(0,m1), m3);
Gael Guennebaud627595a2009-06-10 11:20:30 +0200190
Gael Guennebaud56c7e162009-01-24 15:22:44 +0000191 // count
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100192 VERIFY(((m1.abs()+1)>RealScalar(0.1)).count() == rows*cols);
Benoit Jacobaaaade42010-05-30 16:00:58 -0400193
Gael Guennebaud35c11582011-05-31 22:17:34 +0200194 // and/or
195 VERIFY( (m1<RealScalar(0) && m1>RealScalar(0)).count() == 0);
196 VERIFY( (m1<RealScalar(0) || m1>=RealScalar(0)).count() == rows*cols);
197 RealScalar a = m1.abs().mean();
198 VERIFY( (m1<-a || m1>a).count() == (m1.abs()>a).count());
199
Benoit Jacobaaaade42010-05-30 16:00:58 -0400200 typedef Array<typename ArrayType::Index, Dynamic, 1> ArrayOfIndices;
Gael Guennebaud4ebb8042010-06-02 09:45:57 +0200201
Gael Guennebaud9f795582009-12-16 19:18:40 +0100202 // TODO allows colwise/rowwise for array
Benoit Jacobaaaade42010-05-30 16:00:58 -0400203 VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).colwise().count(), ArrayOfIndices::Constant(cols,rows).transpose());
204 VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).rowwise().count(), ArrayOfIndices::Constant(rows, cols));
Benoit Jacobe8009992008-11-03 22:47:00 +0000205}
206
Benoit Jacobe2775862010-04-28 18:51:38 -0400207template<typename ArrayType> void array_real(const ArrayType& m)
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100208{
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100209 using std::abs;
Gael Guennebaud9b9177f2013-07-05 13:35:34 +0200210 using std::sqrt;
Hauke Heibelf1679c72010-06-20 17:37:56 +0200211 typedef typename ArrayType::Index Index;
Benoit Jacobe2775862010-04-28 18:51:38 -0400212 typedef typename ArrayType::Scalar Scalar;
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100213 typedef typename NumTraits<Scalar>::Real RealScalar;
214
Hauke Heibelf1679c72010-06-20 17:37:56 +0200215 Index rows = m.rows();
216 Index cols = m.cols();
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100217
Benoit Jacobe2775862010-04-28 18:51:38 -0400218 ArrayType m1 = ArrayType::Random(rows, cols),
Gael Guennebaudc695cbf2013-06-24 13:33:44 +0200219 m2 = ArrayType::Random(rows, cols),
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200220 m3(rows, cols),
221 m4 = m1;
Benoit Steiner83149622015-12-10 13:13:45 -0800222
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200223 m4 = (m4.abs()==Scalar(0)).select(1,m4);
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100224
Gael Guennebaud9b1ad5e2012-03-09 12:08:06 +0100225 Scalar s1 = internal::random<Scalar>();
226
Deanna Hood41b717d2015-03-18 03:11:03 +1000227 // these tests are mostly to check possible compilation issues with free-functions.
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100228 VERIFY_IS_APPROX(m1.sin(), sin(m1));
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100229 VERIFY_IS_APPROX(m1.cos(), cos(m1));
Deanna Hoodf89fcef2015-03-11 13:13:30 +1000230 VERIFY_IS_APPROX(m1.tan(), tan(m1));
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100231 VERIFY_IS_APPROX(m1.asin(), asin(m1));
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100232 VERIFY_IS_APPROX(m1.acos(), acos(m1));
Jitse Niesende150b12014-06-19 15:12:33 +0100233 VERIFY_IS_APPROX(m1.atan(), atan(m1));
Deanna Hoodf89fcef2015-03-11 13:13:30 +1000234 VERIFY_IS_APPROX(m1.sinh(), sinh(m1));
235 VERIFY_IS_APPROX(m1.cosh(), cosh(m1));
236 VERIFY_IS_APPROX(m1.tanh(), tanh(m1));
Gael Guennebaud13950562016-05-20 14:58:19 +0200237#if EIGEN_HAS_C99_MATH
Eugene Brevdofa4f9332015-12-07 15:24:49 -0800238 VERIFY_IS_APPROX(m1.lgamma(), lgamma(m1));
Eugene Brevdof7362772015-12-24 21:15:38 -0800239 VERIFY_IS_APPROX(m1.digamma(), digamma(m1));
Eugene Brevdofa4f9332015-12-07 15:24:49 -0800240 VERIFY_IS_APPROX(m1.erf(), erf(m1));
241 VERIFY_IS_APPROX(m1.erfc(), erfc(m1));
Benoit Steiner83149622015-12-10 13:13:45 -0800242#endif // EIGEN_HAS_C99_MATH
Deanna Hooda5e49972015-03-11 08:56:42 +1000243 VERIFY_IS_APPROX(m1.arg(), arg(m1));
Deanna Hood31fdd672015-03-11 06:39:23 +1000244 VERIFY_IS_APPROX(m1.round(), round(m1));
245 VERIFY_IS_APPROX(m1.floor(), floor(m1));
246 VERIFY_IS_APPROX(m1.ceil(), ceil(m1));
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200247 VERIFY((m1.isNaN() == (Eigen::isnan)(m1)).all());
248 VERIFY((m1.isInf() == (Eigen::isinf)(m1)).all());
249 VERIFY((m1.isFinite() == (Eigen::isfinite)(m1)).all());
Deanna Hood41b717d2015-03-18 03:11:03 +1000250 VERIFY_IS_APPROX(m1.inverse(), inverse(m1));
251 VERIFY_IS_APPROX(m1.abs(), abs(m1));
252 VERIFY_IS_APPROX(m1.abs2(), abs2(m1));
253 VERIFY_IS_APPROX(m1.square(), square(m1));
254 VERIFY_IS_APPROX(m1.cube(), cube(m1));
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100255 VERIFY_IS_APPROX(cos(m1+RealScalar(3)*m2), cos((m1+RealScalar(3)*m2).eval()));
Gael Guennebaud3f32f5e2015-11-27 16:27:53 +0100256 VERIFY_IS_APPROX(m1.sign(), sign(m1));
Gael Guennebaud4ebb8042010-06-02 09:45:57 +0200257
Deanna Hood41b717d2015-03-18 03:11:03 +1000258
259 // avoid NaNs with abs() so verification doesn't fail
260 m3 = m1.abs();
261 VERIFY_IS_APPROX(m3.sqrt(), sqrt(abs(m1)));
Benoit Steinere25e3a02015-12-03 18:16:35 -0800262 VERIFY_IS_APPROX(m3.rsqrt(), Scalar(1)/sqrt(abs(m1)));
Deanna Hood41b717d2015-03-18 03:11:03 +1000263 VERIFY_IS_APPROX(m3.log(), log(m3));
Gael Guennebaud89099b02016-06-01 17:00:08 +0200264 VERIFY_IS_APPROX(m3.log1p(), log1p(m3));
Deanna Hood41b717d2015-03-18 03:11:03 +1000265 VERIFY_IS_APPROX(m3.log10(), log10(m3));
266
267
268 VERIFY((!(m1>m2) == (m1<=m2)).all());
269
270 VERIFY_IS_APPROX(sin(m1.asin()), m1);
271 VERIFY_IS_APPROX(cos(m1.acos()), m1);
272 VERIFY_IS_APPROX(tan(m1.atan()), m1);
273 VERIFY_IS_APPROX(sinh(m1), 0.5*(exp(m1)-exp(-m1)));
274 VERIFY_IS_APPROX(cosh(m1), 0.5*(exp(m1)+exp(-m1)));
275 VERIFY_IS_APPROX(tanh(m1), (0.5*(exp(m1)-exp(-m1)))/(0.5*(exp(m1)+exp(-m1))));
Gael Guennebaud2e0ece72015-10-06 17:22:12 +0200276 VERIFY_IS_APPROX(arg(m1), ((m1<0).template cast<Scalar>())*std::acos(-1.0));
Deanna Hood41b717d2015-03-18 03:11:03 +1000277 VERIFY((round(m1) <= ceil(m1) && round(m1) >= floor(m1)).all());
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200278 VERIFY((Eigen::isnan)((m1*0.0)/0.0).all());
279 VERIFY((Eigen::isinf)(m4/0.0).all());
280 VERIFY(((Eigen::isfinite)(m1) && (!(Eigen::isfinite)(m1*0.0/0.0)) && (!(Eigen::isfinite)(m4/0.0))).all());
Deanna Hood41b717d2015-03-18 03:11:03 +1000281 VERIFY_IS_APPROX(inverse(inverse(m1)),m1);
282 VERIFY((abs(m1) == m1 || abs(m1) == -m1).all());
283 VERIFY_IS_APPROX(m3, sqrt(abs2(m1)));
Gael Guennebaud3f32f5e2015-11-27 16:27:53 +0100284 VERIFY_IS_APPROX( m1.sign(), -(-m1).sign() );
285 VERIFY_IS_APPROX( m1*m1.sign(),m1.abs());
286 VERIFY_IS_APPROX(m1.sign() * m1.abs(), m1);
Benoit Jacob5d63d2c2010-04-28 22:42:34 -0400287
Gael Guennebaud62670c82013-06-10 23:40:56 +0200288 VERIFY_IS_APPROX(numext::abs2(numext::real(m1)) + numext::abs2(numext::imag(m1)), numext::abs2(m1));
289 VERIFY_IS_APPROX(numext::abs2(real(m1)) + numext::abs2(imag(m1)), numext::abs2(m1));
Benoit Jacob5d63d2c2010-04-28 22:42:34 -0400290 if(!NumTraits<Scalar>::IsComplex)
Gael Guennebaud62670c82013-06-10 23:40:56 +0200291 VERIFY_IS_APPROX(numext::real(m1), m1);
Gael Guennebaud4ebb8042010-06-02 09:45:57 +0200292
Jitse Niesen6fecb6f2014-02-24 14:10:17 +0000293 // shift argument of logarithm so that it is not zero
294 Scalar smallNumber = NumTraits<Scalar>::dummy_precision();
Deanna Hood41b717d2015-03-18 03:11:03 +1000295 VERIFY_IS_APPROX((m3 + smallNumber).log() , log(abs(m1) + smallNumber));
Gael Guennebaud89099b02016-06-01 17:00:08 +0200296 VERIFY_IS_APPROX((m3 + smallNumber + 1).log() , log1p(abs(m1) + smallNumber));
Gael Guennebaud4ebb8042010-06-02 09:45:57 +0200297
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100298 VERIFY_IS_APPROX(m1.exp() * m2.exp(), exp(m1+m2));
299 VERIFY_IS_APPROX(m1.exp(), exp(m1));
300 VERIFY_IS_APPROX(m1.exp() / m2.exp(),(m1-m2).exp());
Gael Guennebaud575ac542010-06-19 23:17:07 +0200301
Gael Guennebaud575ac542010-06-19 23:17:07 +0200302 VERIFY_IS_APPROX(m3.pow(RealScalar(0.5)), m3.sqrt());
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100303 VERIFY_IS_APPROX(pow(m3,RealScalar(0.5)), m3.sqrt());
Benoit Steinere25e3a02015-12-03 18:16:35 -0800304
305 VERIFY_IS_APPROX(m3.pow(RealScalar(-0.5)), m3.rsqrt());
306 VERIFY_IS_APPROX(pow(m3,RealScalar(-0.5)), m3.rsqrt());
307
Deanna Hood41b717d2015-03-18 03:11:03 +1000308 VERIFY_IS_APPROX(log10(m3), log(m3)/log(10));
Hauke Heibel81c13362012-03-07 08:58:42 +0100309
310 // scalar by array division
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100311 const RealScalar tiny = sqrt(std::numeric_limits<RealScalar>::epsilon());
Hauke Heibeldd9365e2012-03-09 14:04:13 +0100312 s1 += Scalar(tiny);
313 m1 += ArrayType::Constant(rows,cols,Scalar(tiny));
Hauke Heibel81c13362012-03-07 08:58:42 +0100314 VERIFY_IS_APPROX(s1/m1, s1 * m1.inverse());
Eugene Brevdof7362772015-12-24 21:15:38 -0800315
Gael Guennebaudccb408e2016-05-19 18:34:41 +0200316
317
Gael Guennebaud13950562016-05-20 14:58:19 +0200318#if EIGEN_HAS_C99_MATH
Eugene Brevdob1a9afe2016-03-13 15:45:34 -0700319 // check special functions (comparing against numpy implementation)
Gael Guennebaudccb408e2016-05-19 18:34:41 +0200320 if (!NumTraits<Scalar>::IsComplex)
321 {
Eugene Brevdo7ea35bf2016-03-03 19:39:41 -0800322
Eugene Brevdob1a9afe2016-03-13 15:45:34 -0700323 {
324 // Test various propreties of igamma & igammac. These are normalized
325 // gamma integrals where
326 // igammac(a, x) = Gamma(a, x) / Gamma(a)
327 // igamma(a, x) = gamma(a, x) / Gamma(a)
328 // where Gamma and gamma are considered the standard unnormalized
329 // upper and lower incomplete gamma functions, respectively.
330 ArrayType a = m1.abs() + 2;
331 ArrayType x = m2.abs() + 2;
332 ArrayType zero = ArrayType::Zero(rows, cols);
333 ArrayType one = ArrayType::Constant(rows, cols, Scalar(1.0));
334 ArrayType a_m1 = a - one;
335 ArrayType Gamma_a_x = Eigen::igammac(a, x) * a.lgamma().exp();
336 ArrayType Gamma_a_m1_x = Eigen::igammac(a_m1, x) * a_m1.lgamma().exp();
337 ArrayType gamma_a_x = Eigen::igamma(a, x) * a.lgamma().exp();
338 ArrayType gamma_a_m1_x = Eigen::igamma(a_m1, x) * a_m1.lgamma().exp();
339
340 // Gamma(a, 0) == Gamma(a)
341 VERIFY_IS_APPROX(Eigen::igammac(a, zero), one);
342
343 // Gamma(a, x) + gamma(a, x) == Gamma(a)
344 VERIFY_IS_APPROX(Gamma_a_x + gamma_a_x, a.lgamma().exp());
345
346 // Gamma(a, x) == (a - 1) * Gamma(a-1, x) + x^(a-1) * exp(-x)
347 VERIFY_IS_APPROX(Gamma_a_x, (a - 1) * Gamma_a_m1_x + x.pow(a-1) * (-x).exp());
348
349 // gamma(a, x) == (a - 1) * gamma(a-1, x) - x^(a-1) * exp(-x)
350 VERIFY_IS_APPROX(gamma_a_x, (a - 1) * gamma_a_m1_x - x.pow(a-1) * (-x).exp());
351 }
352
353 // Check exact values of igamma and igammac against a third party calculation.
Eugene Brevdo0b9e0ab2016-03-04 21:12:10 -0800354 Scalar a_s[] = {Scalar(0), Scalar(1), Scalar(1.5), Scalar(4), Scalar(0.0001), Scalar(1000.5)};
355 Scalar x_s[] = {Scalar(0), Scalar(1), Scalar(1.5), Scalar(4), Scalar(0.0001), Scalar(1000.5)};
Eugene Brevdo7ea35bf2016-03-03 19:39:41 -0800356
357 // location i*6+j corresponds to a_s[i], x_s[j].
358 Scalar nan = std::numeric_limits<Scalar>::quiet_NaN();
Eugene Brevdo0b9e0ab2016-03-04 21:12:10 -0800359 Scalar igamma_s[][6] = {{0.0, nan, nan, nan, nan, nan},
360 {0.0, 0.6321205588285578, 0.7768698398515702,
361 0.9816843611112658, 9.999500016666262e-05, 1.0},
362 {0.0, 0.4275932955291202, 0.608374823728911,
363 0.9539882943107686, 7.522076445089201e-07, 1.0},
364 {0.0, 0.01898815687615381, 0.06564245437845008,
365 0.5665298796332909, 4.166333347221828e-18, 1.0},
366 {0.0, 0.9999780593618628, 0.9999899967080838,
367 0.9999996219837988, 0.9991370418689945, 1.0},
368 {0.0, 0.0, 0.0, 0.0, 0.0, 0.5042041932513908}};
369 Scalar igammac_s[][6] = {{nan, nan, nan, nan, nan, nan},
370 {1.0, 0.36787944117144233, 0.22313016014842982,
371 0.018315638888734182, 0.9999000049998333, 0.0},
372 {1.0, 0.5724067044708798, 0.3916251762710878,
373 0.04601170568923136, 0.9999992477923555, 0.0},
374 {1.0, 0.9810118431238462, 0.9343575456215499,
375 0.4334701203667089, 1.0, 0.0},
376 {1.0, 2.1940638138146658e-05, 1.0003291916285e-05,
377 3.7801620118431334e-07, 0.0008629581310054535,
378 0.0},
379 {1.0, 1.0, 1.0, 1.0, 1.0, 0.49579580674813944}};
Eugene Brevdo7ea35bf2016-03-03 19:39:41 -0800380 for (int i = 0; i < 6; ++i) {
381 for (int j = 0; j < 6; ++j) {
Eugene Brevdo0b9e0ab2016-03-04 21:12:10 -0800382 if ((std::isnan)(igamma_s[i][j])) {
383 VERIFY((std::isnan)(numext::igamma(a_s[i], x_s[j])));
384 } else {
385 VERIFY_IS_APPROX(numext::igamma(a_s[i], x_s[j]), igamma_s[i][j]);
386 }
387
388 if ((std::isnan)(igammac_s[i][j])) {
389 VERIFY((std::isnan)(numext::igammac(a_s[i], x_s[j])));
390 } else {
391 VERIFY_IS_APPROX(numext::igammac(a_s[i], x_s[j]), igammac_s[i][j]);
392 }
Eugene Brevdo7ea35bf2016-03-03 19:39:41 -0800393 }
394 }
Eugene Brevdof7362772015-12-24 21:15:38 -0800395 }
Eugene Brevdo14897602015-12-24 21:28:18 -0800396#endif // EIGEN_HAS_C99_MATH
Eugene Brevdof7362772015-12-24 21:15:38 -0800397
Gael Guennebaudc695cbf2013-06-24 13:33:44 +0200398 // check inplace transpose
399 m3 = m1;
400 m3.transposeInPlace();
Deanna Hood41b717d2015-03-18 03:11:03 +1000401 VERIFY_IS_APPROX(m3, m1.transpose());
Gael Guennebaudc695cbf2013-06-24 13:33:44 +0200402 m3.transposeInPlace();
Deanna Hood41b717d2015-03-18 03:11:03 +1000403 VERIFY_IS_APPROX(m3, m1);
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100404}
405
Jitse Niesen837db082011-05-09 10:17:41 +0100406template<typename ArrayType> void array_complex(const ArrayType& m)
407{
408 typedef typename ArrayType::Index Index;
Deanna Hood41b717d2015-03-18 03:11:03 +1000409 typedef typename ArrayType::Scalar Scalar;
410 typedef typename NumTraits<Scalar>::Real RealScalar;
Jitse Niesen837db082011-05-09 10:17:41 +0100411
412 Index rows = m.rows();
413 Index cols = m.cols();
414
415 ArrayType m1 = ArrayType::Random(rows, cols),
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200416 m2(rows, cols),
417 m4 = m1;
418
419 m4.real() = (m4.real().abs()==RealScalar(0)).select(RealScalar(1),m4.real());
420 m4.imag() = (m4.imag().abs()==RealScalar(0)).select(RealScalar(1),m4.imag());
Jitse Niesen837db082011-05-09 10:17:41 +0100421
Deanna Hood41b717d2015-03-18 03:11:03 +1000422 Array<RealScalar, -1, -1> m3(rows, cols);
423
Jitse Niesen837db082011-05-09 10:17:41 +0100424 for (Index i = 0; i < m.rows(); ++i)
425 for (Index j = 0; j < m.cols(); ++j)
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100426 m2(i,j) = sqrt(m1(i,j));
Jitse Niesen837db082011-05-09 10:17:41 +0100427
Deanna Hood41b717d2015-03-18 03:11:03 +1000428 // these tests are mostly to check possible compilation issues with free-functions.
Deanna Hoodf89fcef2015-03-11 13:13:30 +1000429 VERIFY_IS_APPROX(m1.sin(), sin(m1));
430 VERIFY_IS_APPROX(m1.cos(), cos(m1));
431 VERIFY_IS_APPROX(m1.tan(), tan(m1));
432 VERIFY_IS_APPROX(m1.sinh(), sinh(m1));
433 VERIFY_IS_APPROX(m1.cosh(), cosh(m1));
434 VERIFY_IS_APPROX(m1.tanh(), tanh(m1));
Deanna Hood41b717d2015-03-18 03:11:03 +1000435 VERIFY_IS_APPROX(m1.arg(), arg(m1));
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200436 VERIFY((m1.isNaN() == (Eigen::isnan)(m1)).all());
437 VERIFY((m1.isInf() == (Eigen::isinf)(m1)).all());
438 VERIFY((m1.isFinite() == (Eigen::isfinite)(m1)).all());
Deanna Hood41b717d2015-03-18 03:11:03 +1000439 VERIFY_IS_APPROX(m1.inverse(), inverse(m1));
440 VERIFY_IS_APPROX(m1.log(), log(m1));
441 VERIFY_IS_APPROX(m1.log10(), log10(m1));
442 VERIFY_IS_APPROX(m1.abs(), abs(m1));
443 VERIFY_IS_APPROX(m1.abs2(), abs2(m1));
444 VERIFY_IS_APPROX(m1.sqrt(), sqrt(m1));
445 VERIFY_IS_APPROX(m1.square(), square(m1));
446 VERIFY_IS_APPROX(m1.cube(), cube(m1));
447 VERIFY_IS_APPROX(cos(m1+RealScalar(3)*m2), cos((m1+RealScalar(3)*m2).eval()));
Gael Guennebaud3f32f5e2015-11-27 16:27:53 +0100448 VERIFY_IS_APPROX(m1.sign(), sign(m1));
Deanna Hood41b717d2015-03-18 03:11:03 +1000449
450
451 VERIFY_IS_APPROX(m1.exp() * m2.exp(), exp(m1+m2));
452 VERIFY_IS_APPROX(m1.exp(), exp(m1));
453 VERIFY_IS_APPROX(m1.exp() / m2.exp(),(m1-m2).exp());
454
455 VERIFY_IS_APPROX(sinh(m1), 0.5*(exp(m1)-exp(-m1)));
456 VERIFY_IS_APPROX(cosh(m1), 0.5*(exp(m1)+exp(-m1)));
457 VERIFY_IS_APPROX(tanh(m1), (0.5*(exp(m1)-exp(-m1)))/(0.5*(exp(m1)+exp(-m1))));
458
459 for (Index i = 0; i < m.rows(); ++i)
460 for (Index j = 0; j < m.cols(); ++j)
461 m3(i,j) = std::atan2(imag(m1(i,j)), real(m1(i,j)));
462 VERIFY_IS_APPROX(arg(m1), m3);
463
464 std::complex<RealScalar> zero(0.0,0.0);
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200465 VERIFY((Eigen::isnan)(m1*zero/zero).all());
Gael Guennebaud4a985e72015-11-20 14:52:08 +0100466#if EIGEN_COMP_MSVC
467 // msvc complex division is not robust
468 VERIFY((Eigen::isinf)(m4/RealScalar(0)).all());
469#else
Gael Guennebaud40f326e2015-06-17 15:33:09 +0200470#if EIGEN_COMP_CLANG
Gael Guennebaud4a985e72015-11-20 14:52:08 +0100471 // clang's complex division is notoriously broken too
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200472 if((numext::isinf)(m4(0,0)/RealScalar(0))) {
Gael Guennebaud40f326e2015-06-17 15:33:09 +0200473#endif
Gael Guennebaud4a985e72015-11-20 14:52:08 +0100474 VERIFY((Eigen::isinf)(m4/zero).all());
Gael Guennebaud40f326e2015-06-17 15:33:09 +0200475#if EIGEN_COMP_CLANG
476 }
477 else
478 {
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200479 VERIFY((Eigen::isinf)(m4.real()/zero.real()).all());
Gael Guennebaud40f326e2015-06-17 15:33:09 +0200480 }
481#endif
Gael Guennebaud4a985e72015-11-20 14:52:08 +0100482#endif // MSVC
483
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200484 VERIFY(((Eigen::isfinite)(m1) && (!(Eigen::isfinite)(m1*zero/zero)) && (!(Eigen::isfinite)(m1/zero))).all());
Deanna Hood41b717d2015-03-18 03:11:03 +1000485
486 VERIFY_IS_APPROX(inverse(inverse(m1)),m1);
487 VERIFY_IS_APPROX(conj(m1.conjugate()), m1);
488 VERIFY_IS_APPROX(abs(m1), sqrt(square(real(m1))+square(imag(m1))));
489 VERIFY_IS_APPROX(abs(m1), sqrt(abs2(m1)));
490 VERIFY_IS_APPROX(log10(m1), log(m1)/log(10));
491
Gael Guennebaud3f32f5e2015-11-27 16:27:53 +0100492 VERIFY_IS_APPROX( m1.sign(), -(-m1).sign() );
493 VERIFY_IS_APPROX( m1.sign() * m1.abs(), m1);
494
Deanna Hood41b717d2015-03-18 03:11:03 +1000495 // scalar by array division
Eugene Brevdof7362772015-12-24 21:15:38 -0800496 Scalar s1 = internal::random<Scalar>();
Benoit Steiner2d74ef92016-05-17 07:20:11 -0700497 const RealScalar tiny = std::sqrt(std::numeric_limits<RealScalar>::epsilon());
Deanna Hood41b717d2015-03-18 03:11:03 +1000498 s1 += Scalar(tiny);
499 m1 += ArrayType::Constant(rows,cols,Scalar(tiny));
500 VERIFY_IS_APPROX(s1/m1, s1 * m1.inverse());
501
502 // check inplace transpose
503 m2 = m1;
504 m2.transposeInPlace();
505 VERIFY_IS_APPROX(m2, m1.transpose());
506 m2.transposeInPlace();
507 VERIFY_IS_APPROX(m2, m1);
Deanna Hood31fdd672015-03-11 06:39:23 +1000508
Jitse Niesen837db082011-05-09 10:17:41 +0100509}
510
Abraham Bachrach039408c2012-01-11 11:00:30 -0500511template<typename ArrayType> void min_max(const ArrayType& m)
512{
513 typedef typename ArrayType::Index Index;
514 typedef typename ArrayType::Scalar Scalar;
515
516 Index rows = m.rows();
517 Index cols = m.cols();
518
519 ArrayType m1 = ArrayType::Random(rows, cols);
520
521 // min/max with array
522 Scalar maxM1 = m1.maxCoeff();
523 Scalar minM1 = m1.minCoeff();
524
525 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, minM1), (m1.min)(ArrayType::Constant(rows,cols, minM1)));
526 VERIFY_IS_APPROX(m1, (m1.min)(ArrayType::Constant(rows,cols, maxM1)));
527
528 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, maxM1), (m1.max)(ArrayType::Constant(rows,cols, maxM1)));
529 VERIFY_IS_APPROX(m1, (m1.max)(ArrayType::Constant(rows,cols, minM1)));
530
531 // min/max with scalar input
532 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, minM1), (m1.min)( minM1));
533 VERIFY_IS_APPROX(m1, (m1.min)( maxM1));
534
535 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, maxM1), (m1.max)( maxM1));
536 VERIFY_IS_APPROX(m1, (m1.max)( minM1));
537
538}
539
Gael Guennebaudccb408e2016-05-19 18:34:41 +0200540template<typename X, typename Y>
541void verify_component_wise(const X& x, const Y& y)
542{
543 for(Index i=0; i<x.size(); ++i)
544 {
545 if((numext::isfinite)(y(i)))
546 VERIFY_IS_APPROX( x(i), y(i) );
547 else if((numext::isnan)(y(i)))
548 VERIFY((numext::isnan)(x(i)));
549 else
550 VERIFY_IS_EQUAL( x(i), y(i) );
551 }
552}
553
554// check special functions (comparing against numpy implementation)
555template<typename ArrayType> void array_special_functions()
556{
557 using std::abs;
558 using std::sqrt;
559 typedef typename ArrayType::Scalar Scalar;
560 typedef typename NumTraits<Scalar>::Real RealScalar;
561
562 Scalar plusinf = std::numeric_limits<Scalar>::infinity();
563 Scalar nan = std::numeric_limits<Scalar>::quiet_NaN();
564
565 // Check the zeta function against scipy.special.zeta
566 {
567 ArrayType x(7), q(7), res(7), ref(7);
568 x << 1.5, 4, 10.5, 10000.5, 3, 1, 0.9;
569 q << 2, 1.5, 3, 1.0001, -2.5, 1.2345, 1.2345;
570 ref << 1.61237534869, 0.234848505667, 1.03086757337e-5, 0.367879440865, 0.054102025820864097, plusinf, nan;
571 CALL_SUBTEST( verify_component_wise(ref, ref); );
572 CALL_SUBTEST( res = x.zeta(q); verify_component_wise(res, ref); );
573 CALL_SUBTEST( res = zeta(x,q); verify_component_wise(res, ref); );
574 }
575
576 // digamma
577 {
578 ArrayType x(7), res(7), ref(7);
579 x << 1, 1.5, 4, -10.5, 10000.5, 0, -1;
580 ref << -0.5772156649015329, 0.03648997397857645, 1.2561176684318, 2.398239129535781, 9.210340372392849, plusinf, plusinf;
581 CALL_SUBTEST( verify_component_wise(ref, ref); );
582
583 CALL_SUBTEST( res = x.digamma(); verify_component_wise(res, ref); );
584 CALL_SUBTEST( res = digamma(x); verify_component_wise(res, ref); );
585 }
586
587
Gael Guennebaud13950562016-05-20 14:58:19 +0200588#if EIGEN_HAS_C99_MATH
Gael Guennebaudccb408e2016-05-19 18:34:41 +0200589 {
590 ArrayType n(11), x(11), res(11), ref(11);
591 n << 1, 1, 1, 1.5, 17, 31, 28, 8, 42, 147, 170;
592 x << 2, 3, 25.5, 1.5, 4.7, 11.8, 17.7, 30.2, 15.8, 54.1, 64;
593 ref << 0.644934066848, 0.394934066848, 0.0399946696496, nan, 293.334565435, 0.445487887616, -2.47810300902e-07, -8.29668781082e-09, -0.434562276666, 0.567742190178, -0.0108615497927;
594 CALL_SUBTEST( verify_component_wise(ref, ref); );
595
Eugene Brevdo39baff82016-06-02 17:04:19 -0700596 if(sizeof(RealScalar)>=8) { // double
597 // Reason for commented line: http://eigen.tuxfamily.org/bz/show_bug.cgi?id=1232
598 // CALL_SUBTEST( res = x.polygamma(n); verify_component_wise(res, ref); );
Gael Guennebaudccb408e2016-05-19 18:34:41 +0200599 CALL_SUBTEST( res = polygamma(n,x); verify_component_wise(res, ref); );
600 }
601 else {
Eugene Brevdo39baff82016-06-02 17:04:19 -0700602 // CALL_SUBTEST( res = x.polygamma(n); verify_component_wise(res.head(8), ref.head(8)); );
Gael Guennebaudccb408e2016-05-19 18:34:41 +0200603 CALL_SUBTEST( res = polygamma(n,x); verify_component_wise(res.head(8), ref.head(8)); );
604 }
Gael Guennebaudccb408e2016-05-19 18:34:41 +0200605 }
Benoit Steiner7aa5bc92016-05-23 14:39:51 -0700606#endif
Eugene Brevdo39baff82016-06-02 17:04:19 -0700607
608#if EIGEN_HAS_C99_MATH
609 {
610 // Inputs and ground truth generated with scipy via:
611 // a = np.logspace(-3, 3, 5) - 1e-3
612 // b = np.logspace(-3, 3, 5) - 1e-3
613 // x = np.linspace(-0.1, 1.1, 5)
614 // (full_a, full_b, full_x) = np.vectorize(lambda a, b, x: (a, b, x))(*np.ix_(a, b, x))
615 // full_a = full_a.flatten().tolist() # same for full_b, full_x
616 // v = scipy.special.betainc(full_a, full_b, full_x).flatten().tolist()
617 //
618 // Note in Eigen, we call betainc with arguments in the order (x, a, b).
619 ArrayType a(125);
620 ArrayType b(125);
621 ArrayType x(125);
622 ArrayType v(125);
623 ArrayType res(125);
624
625 a << 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0,
626 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0,
627 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
628 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
629 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
630 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
631 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
632 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
633 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
634 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
635 0.03062277660168379, 0.999, 0.999, 0.999, 0.999, 0.999, 0.999, 0.999,
636 0.999, 0.999, 0.999, 0.999, 0.999, 0.999, 0.999, 0.999, 0.999, 0.999,
637 0.999, 0.999, 0.999, 0.999, 0.999, 0.999, 0.999, 0.999,
638 31.62177660168379, 31.62177660168379, 31.62177660168379,
639 31.62177660168379, 31.62177660168379, 31.62177660168379,
640 31.62177660168379, 31.62177660168379, 31.62177660168379,
641 31.62177660168379, 31.62177660168379, 31.62177660168379,
642 31.62177660168379, 31.62177660168379, 31.62177660168379,
643 31.62177660168379, 31.62177660168379, 31.62177660168379,
644 31.62177660168379, 31.62177660168379, 31.62177660168379,
645 31.62177660168379, 31.62177660168379, 31.62177660168379,
646 31.62177660168379, 999.999, 999.999, 999.999, 999.999, 999.999, 999.999,
647 999.999, 999.999, 999.999, 999.999, 999.999, 999.999, 999.999, 999.999,
648 999.999, 999.999, 999.999, 999.999, 999.999, 999.999, 999.999, 999.999,
649 999.999, 999.999, 999.999;
650
651 b << 0.0, 0.0, 0.0, 0.0, 0.0, 0.03062277660168379, 0.03062277660168379,
652 0.03062277660168379, 0.03062277660168379, 0.03062277660168379, 0.999,
653 0.999, 0.999, 0.999, 0.999, 31.62177660168379, 31.62177660168379,
654 31.62177660168379, 31.62177660168379, 31.62177660168379, 999.999,
655 999.999, 999.999, 999.999, 999.999, 0.0, 0.0, 0.0, 0.0, 0.0,
656 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
657 0.03062277660168379, 0.03062277660168379, 0.999, 0.999, 0.999, 0.999,
658 0.999, 31.62177660168379, 31.62177660168379, 31.62177660168379,
659 31.62177660168379, 31.62177660168379, 999.999, 999.999, 999.999,
660 999.999, 999.999, 0.0, 0.0, 0.0, 0.0, 0.0, 0.03062277660168379,
661 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
662 0.03062277660168379, 0.999, 0.999, 0.999, 0.999, 0.999,
663 31.62177660168379, 31.62177660168379, 31.62177660168379,
664 31.62177660168379, 31.62177660168379, 999.999, 999.999, 999.999,
665 999.999, 999.999, 0.0, 0.0, 0.0, 0.0, 0.0, 0.03062277660168379,
666 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
667 0.03062277660168379, 0.999, 0.999, 0.999, 0.999, 0.999,
668 31.62177660168379, 31.62177660168379, 31.62177660168379,
669 31.62177660168379, 31.62177660168379, 999.999, 999.999, 999.999,
670 999.999, 999.999, 0.0, 0.0, 0.0, 0.0, 0.0, 0.03062277660168379,
671 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
672 0.03062277660168379, 0.999, 0.999, 0.999, 0.999, 0.999,
673 31.62177660168379, 31.62177660168379, 31.62177660168379,
674 31.62177660168379, 31.62177660168379, 999.999, 999.999, 999.999,
675 999.999, 999.999;
676
677 x << -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5,
678 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2,
679 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1,
680 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1,
681 -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8,
682 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5,
683 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2,
684 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1,
685 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5,
686 0.8, 1.1;
687
688 v << nan, nan, nan, nan, nan, nan, nan, nan, nan, nan, nan, nan, nan, nan,
689 nan, nan, nan, nan, nan, nan, nan, nan, nan, nan, nan, nan, nan, nan,
690 nan, nan, nan, 0.47972119876364683, 0.5, 0.5202788012363533, nan, nan,
691 0.9518683957740043, 0.9789663010413743, 0.9931729188073435, nan, nan,
692 0.999995949033062, 0.9999999999993698, 0.9999999999999999, nan, nan,
693 0.9999999999999999, 0.9999999999999999, 0.9999999999999999, nan, nan,
694 nan, nan, nan, nan, nan, 0.006827081192655869, 0.0210336989586256,
695 0.04813160422599567, nan, nan, 0.20014344256217678, 0.5000000000000001,
696 0.7998565574378232, nan, nan, 0.9991401428435834, 0.999999999698403,
697 0.9999999999999999, nan, nan, 0.9999999999999999, 0.9999999999999999,
698 0.9999999999999999, nan, nan, nan, nan, nan, nan, nan,
699 1.0646600232370887e-25, 6.301722877826246e-13, 4.050966937974938e-06,
700 nan, nan, 7.864342668429763e-23, 3.015969667594166e-10,
701 0.0008598571564165444, nan, nan, 6.031987710123844e-08,
702 0.5000000000000007, 0.9999999396801229, nan, nan, 0.9999999999999999,
703 0.9999999999999999, 0.9999999999999999, nan, nan, nan, nan, nan, nan,
704 nan, 0.0, 7.029920380986636e-306, 2.2450728208591345e-101, nan, nan,
705 0.0, 9.275871147869727e-302, 1.2232913026152827e-97, nan, nan, 0.0,
706 3.0891393081932924e-252, 2.9303043666183996e-60, nan, nan,
707 2.248913486879199e-196, 0.5000000000004947, 0.9999999999999999, nan;
708
709 CALL_SUBTEST(res = betainc(a, b, x);
710 verify_component_wise(res, v););
711 }
Eugene Brevdoc53687d2016-06-05 11:10:30 -0700712
713 // Test various properties of betainc
714 {
715 ArrayType m1 = ArrayType::Random(32);
716 ArrayType m2 = ArrayType::Random(32);
717 ArrayType m3 = ArrayType::Random(32);
718 ArrayType one = ArrayType::Constant(32, Scalar(1.0));
719 const Scalar eps = std::numeric_limits<Scalar>::epsilon();
720 ArrayType a = (m1 * 4.0).exp();
721 ArrayType b = (m2 * 4.0).exp();
722 ArrayType x = m3.abs();
723
724 // betainc(a, 1, x) == x**a
725 CALL_SUBTEST(
726 ArrayType test = betainc(a, one, x);
727 ArrayType expected = x.pow(a);
728 verify_component_wise(test, expected););
729
730 // betainc(1, b, x) == 1 - (1 - x)**b
731 CALL_SUBTEST(
732 ArrayType test = betainc(one, b, x);
733 ArrayType expected = one - (one - x).pow(b);
734 verify_component_wise(test, expected););
735
736 // betainc(a, b, x) == 1 - betainc(b, a, 1-x)
737 CALL_SUBTEST(
738 ArrayType test = betainc(a, b, x) + betainc(b, a, one - x);
739 ArrayType expected = one;
740 verify_component_wise(test, expected););
741
742 // betainc(a+1, b, x) = betainc(a, b, x) - x**a * (1 - x)**b / (a * beta(a, b))
743 CALL_SUBTEST(
744 ArrayType num = x.pow(a) * (one - x).pow(b);
745 ArrayType denom = a * (a.lgamma() + b.lgamma() - (a + b).lgamma()).exp();
746 // Add eps to rhs and lhs so that component-wise test doesn't result in
747 // nans when both outputs are zeros.
748 ArrayType expected = betainc(a, b, x) - num / denom + eps;
749 ArrayType test = betainc(a + one, b, x) + eps;
750 if (sizeof(Scalar) >= 8) { // double
751 verify_component_wise(test, expected);
752 } else {
753 // Reason for limited test: http://eigen.tuxfamily.org/bz/show_bug.cgi?id=1232
754 verify_component_wise(test.head(8), expected.head(8));
755 });
756
757 // betainc(a, b+1, x) = betainc(a, b, x) + x**a * (1 - x)**b / (b * beta(a, b))
758 CALL_SUBTEST(
759 // Add eps to rhs and lhs so that component-wise test doesn't result in
760 // nans when both outputs are zeros.
761 ArrayType num = x.pow(a) * (one - x).pow(b);
762 ArrayType denom = b * (a.lgamma() + b.lgamma() - (a + b).lgamma()).exp();
763 ArrayType expected = betainc(a, b, x) + num / denom + eps;
764 ArrayType test = betainc(a, b + one, x) + eps;
765 verify_component_wise(test, expected););
766 }
Eugene Brevdo39baff82016-06-02 17:04:19 -0700767#endif
Gael Guennebaudccb408e2016-05-19 18:34:41 +0200768}
769
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000770void test_array()
771{
Gael Guennebaud8d6bd562016-05-20 14:45:33 +0200772#ifndef EIGEN_HAS_C99_MATH
773 std::cerr << "WARNING: testing of special math functions disabled" << std::endl;
774#endif
775
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000776 for(int i = 0; i < g_repeat; i++) {
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100777 CALL_SUBTEST_1( array(Array<float, 1, 1>()) );
Hauke Heibelb31e1242010-12-16 19:07:23 +0100778 CALL_SUBTEST_2( array(Array22f()) );
779 CALL_SUBTEST_3( array(Array44d()) );
Gael Guennebauda8f66fe2011-07-12 14:41:00 +0200780 CALL_SUBTEST_4( array(ArrayXXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
781 CALL_SUBTEST_5( array(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
782 CALL_SUBTEST_6( array(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000783 }
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100784 for(int i = 0; i < g_repeat; i++) {
785 CALL_SUBTEST_1( comparisons(Array<float, 1, 1>()) );
Hauke Heibelb31e1242010-12-16 19:07:23 +0100786 CALL_SUBTEST_2( comparisons(Array22f()) );
787 CALL_SUBTEST_3( comparisons(Array44d()) );
Gael Guennebauda8f66fe2011-07-12 14:41:00 +0200788 CALL_SUBTEST_5( comparisons(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
789 CALL_SUBTEST_6( comparisons(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100790 }
791 for(int i = 0; i < g_repeat; i++) {
Abraham Bachrach039408c2012-01-11 11:00:30 -0500792 CALL_SUBTEST_1( min_max(Array<float, 1, 1>()) );
793 CALL_SUBTEST_2( min_max(Array22f()) );
794 CALL_SUBTEST_3( min_max(Array44d()) );
795 CALL_SUBTEST_5( min_max(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
796 CALL_SUBTEST_6( min_max(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
797 }
798 for(int i = 0; i < g_repeat; i++) {
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100799 CALL_SUBTEST_1( array_real(Array<float, 1, 1>()) );
Hauke Heibelb31e1242010-12-16 19:07:23 +0100800 CALL_SUBTEST_2( array_real(Array22f()) );
801 CALL_SUBTEST_3( array_real(Array44d()) );
Gael Guennebauda8f66fe2011-07-12 14:41:00 +0200802 CALL_SUBTEST_5( array_real(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100803 }
Jitse Niesen837db082011-05-09 10:17:41 +0100804 for(int i = 0; i < g_repeat; i++) {
Gael Guennebauda8f66fe2011-07-12 14:41:00 +0200805 CALL_SUBTEST_4( array_complex(ArrayXXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
Jitse Niesen837db082011-05-09 10:17:41 +0100806 }
Benoit Jacob5d63d2c2010-04-28 22:42:34 -0400807
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100808 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<int>::type, int >::value));
809 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<float>::type, float >::value));
810 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Array2i>::type, ArrayBase<Array2i> >::value));
Gael Guennebaud64fcfd32016-06-14 11:26:57 +0200811 typedef CwiseUnaryOp<internal::scalar_abs_op<double>, ArrayXd > Xpr;
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100812 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Xpr>::type,
813 ArrayBase<Xpr>
814 >::value));
Gael Guennebaudccb408e2016-05-19 18:34:41 +0200815
816 CALL_SUBTEST_7(array_special_functions<ArrayXf>());
817 CALL_SUBTEST_7(array_special_functions<ArrayXd>());
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000818}