• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 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 
11 #ifndef EIGEN_SPECIALFUNCTIONS_ARRAYAPI_H
12 #define EIGEN_SPECIALFUNCTIONS_ARRAYAPI_H
13 
14 namespace Eigen {
15 
16 /** \cpp11 \returns an expression of the coefficient-wise igamma(\a a, \a x) to the given arrays.
17   *
18   * This function computes the coefficient-wise incomplete gamma function.
19   *
20   * \note This function supports only float and double scalar types in c++11 mode. To support other scalar types,
21   * or float/double in non c++11 mode, the user has to provide implementations of igammac(T,T) for any scalar
22   * type T to be supported.
23   *
24   * \sa Eigen::igammac(), Eigen::lgamma()
25   */
26 template<typename Derived,typename ExponentDerived>
27 inline const Eigen::CwiseBinaryOp<Eigen::internal::scalar_igamma_op<typename Derived::Scalar>, const Derived, const ExponentDerived>
igamma(const Eigen::ArrayBase<Derived> & a,const Eigen::ArrayBase<ExponentDerived> & x)28 igamma(const Eigen::ArrayBase<Derived>& a, const Eigen::ArrayBase<ExponentDerived>& x)
29 {
30   return Eigen::CwiseBinaryOp<Eigen::internal::scalar_igamma_op<typename Derived::Scalar>, const Derived, const ExponentDerived>(
31     a.derived(),
32     x.derived()
33   );
34 }
35 
36 /** \cpp11 \returns an expression of the coefficient-wise igammac(\a a, \a x) to the given arrays.
37   *
38   * This function computes the coefficient-wise complementary incomplete gamma function.
39   *
40   * \note This function supports only float and double scalar types in c++11 mode. To support other scalar types,
41   * or float/double in non c++11 mode, the user has to provide implementations of igammac(T,T) for any scalar
42   * type T to be supported.
43   *
44   * \sa Eigen::igamma(), Eigen::lgamma()
45   */
46 template<typename Derived,typename ExponentDerived>
47 inline const Eigen::CwiseBinaryOp<Eigen::internal::scalar_igammac_op<typename Derived::Scalar>, const Derived, const ExponentDerived>
igammac(const Eigen::ArrayBase<Derived> & a,const Eigen::ArrayBase<ExponentDerived> & x)48 igammac(const Eigen::ArrayBase<Derived>& a, const Eigen::ArrayBase<ExponentDerived>& x)
49 {
50   return Eigen::CwiseBinaryOp<Eigen::internal::scalar_igammac_op<typename Derived::Scalar>, const Derived, const ExponentDerived>(
51     a.derived(),
52     x.derived()
53   );
54 }
55 
56 /** \cpp11 \returns an expression of the coefficient-wise polygamma(\a n, \a x) to the given arrays.
57   *
58   * It returns the \a n -th derivative of the digamma(psi) evaluated at \c x.
59   *
60   * \note This function supports only float and double scalar types in c++11 mode. To support other scalar types,
61   * or float/double in non c++11 mode, the user has to provide implementations of polygamma(T,T) for any scalar
62   * type T to be supported.
63   *
64   * \sa Eigen::digamma()
65   */
66 // * \warning Be careful with the order of the parameters: x.polygamma(n) is equivalent to polygamma(n,x)
67 // * \sa ArrayBase::polygamma()
68 template<typename DerivedN,typename DerivedX>
69 inline const Eigen::CwiseBinaryOp<Eigen::internal::scalar_polygamma_op<typename DerivedX::Scalar>, const DerivedN, const DerivedX>
polygamma(const Eigen::ArrayBase<DerivedN> & n,const Eigen::ArrayBase<DerivedX> & x)70 polygamma(const Eigen::ArrayBase<DerivedN>& n, const Eigen::ArrayBase<DerivedX>& x)
71 {
72   return Eigen::CwiseBinaryOp<Eigen::internal::scalar_polygamma_op<typename DerivedX::Scalar>, const DerivedN, const DerivedX>(
73     n.derived(),
74     x.derived()
75   );
76 }
77 
78 /** \cpp11 \returns an expression of the coefficient-wise betainc(\a x, \a a, \a b) to the given arrays.
79   *
80   * This function computes the regularized incomplete beta function (integral).
81   *
82   * \note This function supports only float and double scalar types in c++11 mode. To support other scalar types,
83   * or float/double in non c++11 mode, the user has to provide implementations of betainc(T,T,T) for any scalar
84   * type T to be supported.
85   *
86   * \sa Eigen::betainc(), Eigen::lgamma()
87   */
88 template<typename ArgADerived, typename ArgBDerived, typename ArgXDerived>
89 inline const Eigen::CwiseTernaryOp<Eigen::internal::scalar_betainc_op<typename ArgXDerived::Scalar>, const ArgADerived, const ArgBDerived, const ArgXDerived>
betainc(const Eigen::ArrayBase<ArgADerived> & a,const Eigen::ArrayBase<ArgBDerived> & b,const Eigen::ArrayBase<ArgXDerived> & x)90 betainc(const Eigen::ArrayBase<ArgADerived>& a, const Eigen::ArrayBase<ArgBDerived>& b, const Eigen::ArrayBase<ArgXDerived>& x)
91 {
92   return Eigen::CwiseTernaryOp<Eigen::internal::scalar_betainc_op<typename ArgXDerived::Scalar>, const ArgADerived, const ArgBDerived, const ArgXDerived>(
93     a.derived(),
94     b.derived(),
95     x.derived()
96   );
97 }
98 
99 
100 /** \returns an expression of the coefficient-wise zeta(\a x, \a q) to the given arrays.
101   *
102   * It returns the Riemann zeta function of two arguments \a x and \a q:
103   *
104   * \param x is the exposent, it must be > 1
105   * \param q is the shift, it must be > 0
106   *
107   * \note This function supports only float and double scalar types. To support other scalar types, the user has
108   * to provide implementations of zeta(T,T) for any scalar type T to be supported.
109   *
110   * \sa ArrayBase::zeta()
111   */
112 template<typename DerivedX,typename DerivedQ>
113 inline const Eigen::CwiseBinaryOp<Eigen::internal::scalar_zeta_op<typename DerivedX::Scalar>, const DerivedX, const DerivedQ>
zeta(const Eigen::ArrayBase<DerivedX> & x,const Eigen::ArrayBase<DerivedQ> & q)114 zeta(const Eigen::ArrayBase<DerivedX>& x, const Eigen::ArrayBase<DerivedQ>& q)
115 {
116   return Eigen::CwiseBinaryOp<Eigen::internal::scalar_zeta_op<typename DerivedX::Scalar>, const DerivedX, const DerivedQ>(
117     x.derived(),
118     q.derived()
119   );
120 }
121 
122 } // end namespace Eigen
123 
124 #endif // EIGEN_SPECIALFUNCTIONS_ARRAYAPI_H
125