This bugzilla service is closed. All entries have been migrated to https://gitlab.com/libeigen/eigen
View | Details | Raw Unified | Return to bug 564
Collapse All | Expand All

(-)a/Eigen/src/Core/BooleanRedux.h (+20 lines)
Lines 124-134 inline bool DenseBase<Derived>::any() co Link Here
124
  * \sa all(), any()
124
  * \sa all(), any()
125
  */
125
  */
126
template<typename Derived>
126
template<typename Derived>
127
inline typename DenseBase<Derived>::Index DenseBase<Derived>::count() const
127
inline typename DenseBase<Derived>::Index DenseBase<Derived>::count() const
128
{
128
{
129
  return derived().template cast<bool>().template cast<Index>().sum();
129
  return derived().template cast<bool>().template cast<Index>().sum();
130
}
130
}
131
131
132
/** \returns true is \c *this contains at least one Not A Number (NaN).
133
  *
134
  * \sa isFinite()
135
  */
136
template<typename Derived>
137
inline bool DenseBase<Derived>::hasNaN() const
138
{
139
  return !((derived().array()==derived().array()).all());
140
}
141
142
/** \returns true if \c *this contains only finite numbers, i.e., no NaN and no +/-INF values.
143
  *
144
  * \sa hasNaN()
145
  */
146
template<typename Derived>
147
inline bool DenseBase<Derived>::isFinite() const
148
{
149
  return !((derived()-derived()).hasNaN());
150
}
151
    
132
} // end namespace Eigen
152
} // end namespace Eigen
133
153
134
#endif // EIGEN_ALLANDANY_H
154
#endif // EIGEN_ALLANDANY_H
(-)a/Eigen/src/Core/DenseBase.h (-2 / +3 lines)
Lines 331-346 template<typename Derived> class DenseBa Link Here
331
    template<typename OtherDerived>
331
    template<typename OtherDerived>
332
    bool isMuchSmallerThan(const DenseBase<OtherDerived>& other,
332
    bool isMuchSmallerThan(const DenseBase<OtherDerived>& other,
333
                           const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
333
                           const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
334
334
335
    bool isApproxToConstant(const Scalar& value, const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
335
    bool isApproxToConstant(const Scalar& value, const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
336
    bool isConstant(const Scalar& value, const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
336
    bool isConstant(const Scalar& value, const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
337
    bool isZero(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
337
    bool isZero(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
338
    bool isOnes(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
338
    bool isOnes(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
339
    
340
    inline bool hasNaN() const;
341
    inline bool isFinite() const;
339
342
340
    inline Derived& operator*=(const Scalar& other);
343
    inline Derived& operator*=(const Scalar& other);
341
    inline Derived& operator/=(const Scalar& other);
344
    inline Derived& operator/=(const Scalar& other);
342
345
343
    typedef typename internal::add_const_on_value_type<typename internal::eval<Derived>::type>::type EvalReturnType;
346
    typedef typename internal::add_const_on_value_type<typename internal::eval<Derived>::type>::type EvalReturnType;
344
    /** \returns the matrix or vector obtained by evaluating this expression.
347
    /** \returns the matrix or vector obtained by evaluating this expression.
345
      *
348
      *
346
      * Notice that in the case of a plain matrix or vector (not an expression) this function just returns
349
      * Notice that in the case of a plain matrix or vector (not an expression) this function just returns
Lines 410-427 template<typename Derived> class DenseBa Link Here
410
    /** \returns the unique coefficient of a 1x1 expression */
413
    /** \returns the unique coefficient of a 1x1 expression */
411
    CoeffReturnType value() const
414
    CoeffReturnType value() const
412
    {
415
    {
413
      EIGEN_STATIC_ASSERT_SIZE_1x1(Derived)
416
      EIGEN_STATIC_ASSERT_SIZE_1x1(Derived)
414
      eigen_assert(this->rows() == 1 && this->cols() == 1);
417
      eigen_assert(this->rows() == 1 && this->cols() == 1);
415
      return derived().coeff(0,0);
418
      return derived().coeff(0,0);
416
    }
419
    }
417
420
418
/////////// Array module ///////////
419
420
    bool all(void) const;
421
    bool all(void) const;
421
    bool any(void) const;
422
    bool any(void) const;
422
    Index count() const;
423
    Index count() const;
423
424
424
    typedef VectorwiseOp<Derived, Horizontal> RowwiseReturnType;
425
    typedef VectorwiseOp<Derived, Horizontal> RowwiseReturnType;
425
    typedef const VectorwiseOp<const Derived, Horizontal> ConstRowwiseReturnType;
426
    typedef const VectorwiseOp<const Derived, Horizontal> ConstRowwiseReturnType;
426
    typedef VectorwiseOp<Derived, Vertical> ColwiseReturnType;
427
    typedef VectorwiseOp<Derived, Vertical> ColwiseReturnType;
427
    typedef const VectorwiseOp<const Derived, Vertical> ConstColwiseReturnType;
428
    typedef const VectorwiseOp<const Derived, Vertical> ConstColwiseReturnType;
(-)a/test/CMakeLists.txt (+1 lines)
Lines 217-232 ei_add_test(sparse_permutations) Link Here
217
ei_add_test(nullary)
217
ei_add_test(nullary)
218
ei_add_test(nesting_ops "${CMAKE_CXX_FLAGS_DEBUG}")
218
ei_add_test(nesting_ops "${CMAKE_CXX_FLAGS_DEBUG}")
219
ei_add_test(zerosized)
219
ei_add_test(zerosized)
220
ei_add_test(dontalign)
220
ei_add_test(dontalign)
221
ei_add_test(evaluators)
221
ei_add_test(evaluators)
222
ei_add_test(sizeoverflow)
222
ei_add_test(sizeoverflow)
223
ei_add_test(prec_inverse_4x4)
223
ei_add_test(prec_inverse_4x4)
224
ei_add_test(vectorwiseop)
224
ei_add_test(vectorwiseop)
225
ei_add_test(special_numbers)
225
226
226
ei_add_test(simplicial_cholesky)
227
ei_add_test(simplicial_cholesky)
227
ei_add_test(conjugate_gradient)
228
ei_add_test(conjugate_gradient)
228
ei_add_test(bicgstab)
229
ei_add_test(bicgstab)
229
ei_add_test(sparselu)
230
ei_add_test(sparselu)
230
ei_add_test(sparseqr)
231
ei_add_test(sparseqr)
231
232
232
# ei_add_test(denseLM)
233
# ei_add_test(denseLM)
(-)a/test/special_numbers.cpp (+59 lines)
Line 0 Link Here
1
// This file is part of Eigen, a lightweight C++ template library
2
// for linear algebra.
3
//
4
// Copyright (C) 2013 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
#include "main.h"
11
12
template<typename Scalar> void special_numbers()
13
{
14
  typedef typename NumTraits<Scalar>::Real RealScalar;
15
  typedef Matrix<Scalar, Dynamic,Dynamic> MatType;
16
  int rows = internal::random<int>(1,300);
17
  int cols = internal::random<int>(1,300);
18
  
19
  Scalar nan = Scalar(0)/Scalar(0);
20
  Scalar inf = Scalar(1)/Scalar(0);
21
  Scalar s1 = internal::random<Scalar>();
22
  
23
  MatType m1    = MatType::Random(rows,cols),
24
          mnan  = MatType::Random(rows,cols),
25
          minf  = MatType::Random(rows,cols),
26
          mboth = MatType::Random(rows,cols);
27
          
28
  int n = internal::random<int>(1,10);
29
  for(int k=0; k<n; ++k)
30
  {
31
    mnan(internal::random<int>(0,rows-1), internal::random<int>(0,cols-1)) = nan;
32
    minf(internal::random<int>(0,rows-1), internal::random<int>(0,cols-1)) = inf;
33
  }
34
  mboth = mnan + minf;
35
  
36
  VERIFY(!m1.hasNaN());
37
  VERIFY(m1.isFinite());
38
  
39
  VERIFY(mnan.hasNaN());
40
  VERIFY((s1*mnan).hasNaN());
41
  VERIFY(!minf.hasNaN());
42
  VERIFY(!(2*minf).hasNaN());
43
  VERIFY(mboth.hasNaN());
44
  VERIFY(mboth.array().hasNaN());
45
  
46
  VERIFY(!mnan.isFinite());
47
  VERIFY(!minf.isFinite());
48
  VERIFY(!(minf-mboth).isFinite());
49
  VERIFY(!mboth.isFinite());
50
  VERIFY(!mboth.array().isFinite());
51
}
52
53
void test_special_numbers()
54
{
55
  for(int i = 0; i < 10*g_repeat; i++) {
56
    CALL_SUBTEST_1( special_numbers<float>() );
57
    CALL_SUBTEST_1( special_numbers<double>() );
58
  }
59
}

Return to bug 564