00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025 #include "main.h"
00026
00027 template<typename ArrayType> void array(const ArrayType& m)
00028 {
00029 typedef typename ArrayType::Index Index;
00030 typedef typename ArrayType::Scalar Scalar;
00031 typedef typename NumTraits<Scalar>::Real RealScalar;
00032 typedef Array<Scalar, ArrayType::RowsAtCompileTime, 1> ColVectorType;
00033 typedef Array<Scalar, 1, ArrayType::ColsAtCompileTime> RowVectorType;
00034
00035 Index rows = m.rows();
00036 Index cols = m.cols();
00037
00038 ArrayType m1 = ArrayType::Random(rows, cols),
00039 m2 = ArrayType::Random(rows, cols),
00040 m3(rows, cols);
00041
00042 ColVectorType cv1 = ColVectorType::Random(rows);
00043 RowVectorType rv1 = RowVectorType::Random(cols);
00044
00045 Scalar s1 = internal::random<Scalar>(),
00046 s2 = internal::random<Scalar>();
00047
00048
00049 VERIFY_IS_APPROX(m1 + s1, s1 + m1);
00050 VERIFY_IS_APPROX(m1 + s1, ArrayType::Constant(rows,cols,s1) + m1);
00051 VERIFY_IS_APPROX(s1 - m1, (-m1)+s1 );
00052 VERIFY_IS_APPROX(m1 - s1, m1 - ArrayType::Constant(rows,cols,s1));
00053 VERIFY_IS_APPROX(s1 - m1, ArrayType::Constant(rows,cols,s1) - m1);
00054 VERIFY_IS_APPROX((m1*Scalar(2)) - s2, (m1+m1) - ArrayType::Constant(rows,cols,s2) );
00055 m3 = m1;
00056 m3 += s2;
00057 VERIFY_IS_APPROX(m3, m1 + s2);
00058 m3 = m1;
00059 m3 -= s1;
00060 VERIFY_IS_APPROX(m3, m1 - s1);
00061
00062
00063 m3 = m1;
00064 ArrayType::Map(m1.data(), m1.rows(), m1.cols()) -= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
00065 VERIFY_IS_APPROX(m1, m3 - m2);
00066
00067 m3 = m1;
00068 ArrayType::Map(m1.data(), m1.rows(), m1.cols()) += ArrayType::Map(m2.data(), m2.rows(), m2.cols());
00069 VERIFY_IS_APPROX(m1, m3 + m2);
00070
00071 m3 = m1;
00072 ArrayType::Map(m1.data(), m1.rows(), m1.cols()) *= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
00073 VERIFY_IS_APPROX(m1, m3 * m2);
00074
00075 m3 = m1;
00076 m2 = ArrayType::Random(rows,cols);
00077 m2 = (m2==0).select(1,m2);
00078 ArrayType::Map(m1.data(), m1.rows(), m1.cols()) /= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
00079 VERIFY_IS_APPROX(m1, m3 / m2);
00080
00081
00082 VERIFY_IS_APPROX(m1.colwise().sum().sum(), m1.sum());
00083 VERIFY_IS_APPROX(m1.rowwise().sum().sum(), m1.sum());
00084 if (!internal::isApprox(m1.sum(), (m1+m2).sum(), test_precision<Scalar>()))
00085 VERIFY_IS_NOT_APPROX(((m1+m2).rowwise().sum()).sum(), m1.sum());
00086 VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(internal::scalar_sum_op<Scalar>()));
00087
00088
00089 m3 = m1;
00090 VERIFY_IS_APPROX(m3.colwise() += cv1, m1.colwise() + cv1);
00091 m3 = m1;
00092 VERIFY_IS_APPROX(m3.colwise() -= cv1, m1.colwise() - cv1);
00093 m3 = m1;
00094 VERIFY_IS_APPROX(m3.rowwise() += rv1, m1.rowwise() + rv1);
00095 m3 = m1;
00096 VERIFY_IS_APPROX(m3.rowwise() -= rv1, m1.rowwise() - rv1);
00097 }
00098
00099 template<typename ArrayType> void comparisons(const ArrayType& m)
00100 {
00101 typedef typename ArrayType::Index Index;
00102 typedef typename ArrayType::Scalar Scalar;
00103 typedef typename NumTraits<Scalar>::Real RealScalar;
00104 typedef Array<Scalar, ArrayType::RowsAtCompileTime, 1> VectorType;
00105
00106 Index rows = m.rows();
00107 Index cols = m.cols();
00108
00109 Index r = internal::random<Index>(0, rows-1),
00110 c = internal::random<Index>(0, cols-1);
00111
00112 ArrayType m1 = ArrayType::Random(rows, cols),
00113 m2 = ArrayType::Random(rows, cols),
00114 m3(rows, cols);
00115
00116 VERIFY(((m1 + Scalar(1)) > m1).all());
00117 VERIFY(((m1 - Scalar(1)) < m1).all());
00118 if (rows*cols>1)
00119 {
00120 m3 = m1;
00121 m3(r,c) += 1;
00122 VERIFY(! (m1 < m3).all() );
00123 VERIFY(! (m1 > m3).all() );
00124 }
00125
00126
00127 VERIFY( (m1 != (m1(r,c)+1) ).any() );
00128 VERIFY( (m1 > (m1(r,c)-1) ).any() );
00129 VERIFY( (m1 < (m1(r,c)+1) ).any() );
00130 VERIFY( (m1 == m1(r,c) ).any() );
00131
00132
00133 VERIFY_IS_APPROX( (m1<m2).select(m1,m2), m1.cwiseMin(m2) );
00134 VERIFY_IS_APPROX( (m1>m2).select(m1,m2), m1.cwiseMax(m2) );
00135 Scalar mid = (m1.cwiseAbs().minCoeff() + m1.cwiseAbs().maxCoeff())/Scalar(2);
00136 for (int j=0; j<cols; ++j)
00137 for (int i=0; i<rows; ++i)
00138 m3(i,j) = internal::abs(m1(i,j))<mid ? 0 : m1(i,j);
00139 VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid))
00140 .select(ArrayType::Zero(rows,cols),m1), m3);
00141
00142 VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid))
00143 .select(0,m1), m3);
00144 VERIFY_IS_APPROX( (m1.abs()>=ArrayType::Constant(rows,cols,mid))
00145 .select(m1,0), m3);
00146
00147 VERIFY_IS_APPROX( (m1.abs()<mid).select(0,m1), m3);
00148
00149
00150 VERIFY(((m1.abs()+1)>RealScalar(0.1)).count() == rows*cols);
00151
00152 typedef Array<typename ArrayType::Index, Dynamic, 1> ArrayOfIndices;
00153
00154
00155 VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).colwise().count(), ArrayOfIndices::Constant(cols,rows).transpose());
00156 VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).rowwise().count(), ArrayOfIndices::Constant(rows, cols));
00157 }
00158
00159 template<typename ArrayType> void array_real(const ArrayType& m)
00160 {
00161 typedef typename ArrayType::Index Index;
00162 typedef typename ArrayType::Scalar Scalar;
00163 typedef typename NumTraits<Scalar>::Real RealScalar;
00164
00165 Index rows = m.rows();
00166 Index cols = m.cols();
00167
00168 ArrayType m1 = ArrayType::Random(rows, cols),
00169 m2 = ArrayType::Random(rows, cols),
00170 m3(rows, cols);
00171
00172
00173 VERIFY_IS_APPROX(m1.sin(), std::sin(m1));
00174 VERIFY_IS_APPROX(m1.sin(), internal::sin(m1));
00175 VERIFY_IS_APPROX(m1.cos(), std::cos(m1));
00176 VERIFY_IS_APPROX(m1.cos(), internal::cos(m1));
00177 VERIFY_IS_APPROX(m1.asin(), std::asin(m1));
00178 VERIFY_IS_APPROX(m1.asin(), internal::asin(m1));
00179 VERIFY_IS_APPROX(m1.acos(), std::acos(m1));
00180 VERIFY_IS_APPROX(m1.acos(), internal::acos(m1));
00181 VERIFY_IS_APPROX(m1.tan(), std::tan(m1));
00182 VERIFY_IS_APPROX(m1.tan(), internal::tan(m1));
00183
00184 VERIFY_IS_APPROX(internal::cos(m1+RealScalar(3)*m2), internal::cos((m1+RealScalar(3)*m2).eval()));
00185 VERIFY_IS_APPROX(std::cos(m1+RealScalar(3)*m2), std::cos((m1+RealScalar(3)*m2).eval()));
00186
00187 VERIFY_IS_APPROX(m1.abs().sqrt(), std::sqrt(std::abs(m1)));
00188 VERIFY_IS_APPROX(m1.abs().sqrt(), internal::sqrt(internal::abs(m1)));
00189 VERIFY_IS_APPROX(m1.abs(), internal::sqrt(internal::abs2(m1)));
00190
00191 VERIFY_IS_APPROX(internal::abs2(internal::real(m1)) + internal::abs2(internal::imag(m1)), internal::abs2(m1));
00192 VERIFY_IS_APPROX(internal::abs2(std::real(m1)) + internal::abs2(std::imag(m1)), internal::abs2(m1));
00193 if(!NumTraits<Scalar>::IsComplex)
00194 VERIFY_IS_APPROX(internal::real(m1), m1);
00195
00196 VERIFY_IS_APPROX(m1.abs().log(), std::log(std::abs(m1)));
00197 VERIFY_IS_APPROX(m1.abs().log(), internal::log(internal::abs(m1)));
00198
00199 VERIFY_IS_APPROX(m1.exp(), std::exp(m1));
00200 VERIFY_IS_APPROX(m1.exp() * m2.exp(), std::exp(m1+m2));
00201 VERIFY_IS_APPROX(m1.exp(), internal::exp(m1));
00202 VERIFY_IS_APPROX(m1.exp() / m2.exp(), std::exp(m1-m2));
00203
00204 VERIFY_IS_APPROX(m1.pow(2), m1.square());
00205 VERIFY_IS_APPROX(std::pow(m1,2), m1.square());
00206 m3 = m1.abs();
00207 VERIFY_IS_APPROX(m3.pow(RealScalar(0.5)), m3.sqrt());
00208 VERIFY_IS_APPROX(std::pow(m3,RealScalar(0.5)), m3.sqrt());
00209 }
00210
00211 void test_array()
00212 {
00213 for(int i = 0; i < g_repeat; i++) {
00214 CALL_SUBTEST_1( array(Array<float, 1, 1>()) );
00215 CALL_SUBTEST_2( array(Array22f()) );
00216 CALL_SUBTEST_3( array(Array44d()) );
00217 CALL_SUBTEST_4( array(ArrayXXcf(3, 3)) );
00218 CALL_SUBTEST_5( array(ArrayXXf(8, 12)) );
00219 CALL_SUBTEST_6( array(ArrayXXi(8, 12)) );
00220 }
00221 for(int i = 0; i < g_repeat; i++) {
00222 CALL_SUBTEST_1( comparisons(Array<float, 1, 1>()) );
00223 CALL_SUBTEST_2( comparisons(Array22f()) );
00224 CALL_SUBTEST_3( comparisons(Array44d()) );
00225 CALL_SUBTEST_5( comparisons(ArrayXXf(8, 12)) );
00226 CALL_SUBTEST_6( comparisons(ArrayXXi(8, 12)) );
00227 }
00228 for(int i = 0; i < g_repeat; i++) {
00229 CALL_SUBTEST_1( array_real(Array<float, 1, 1>()) );
00230 CALL_SUBTEST_2( array_real(Array22f()) );
00231 CALL_SUBTEST_3( array_real(Array44d()) );
00232 CALL_SUBTEST_5( array_real(ArrayXXf(8, 12)) );
00233 }
00234
00235 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<int>::type, int >::value));
00236 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<float>::type, float >::value));
00237 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Array2i>::type, ArrayBase<Array2i> >::value));
00238 typedef CwiseUnaryOp<internal::scalar_sum_op<double>, ArrayXd > Xpr;
00239 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Xpr>::type,
00240 ArrayBase<Xpr>
00241 >::value));
00242 }