aboutsummaryrefslogtreecommitdiffhomepage
path: root/unsupported/Eigen/src/SpecialFunctions/SpecialFunctionsArrayAPI.h
blob: 30cdf4751791a0e1efbb36e033fe68f20e07bbfa (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2016 Gael Guennebaud <gael.guennebaud@inria.fr>
//
// This Source Code Form is subject to the terms of the Mozilla
// Public License v. 2.0. If a copy of the MPL was not distributed
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.


#ifndef EIGEN_SPECIALFUNCTIONS_ARRAYAPI_H
#define EIGEN_SPECIALFUNCTIONS_ARRAYAPI_H

namespace Eigen {

/** \cpp11 \returns an expression of the coefficient-wise igamma(\a a, \a x) to the given arrays.
  *
  * This function computes the coefficient-wise incomplete gamma function.
  *
  * \note This function supports only float and double scalar types in c++11 mode. To support other scalar types,
  * or float/double in non c++11 mode, the user has to provide implementations of igammac(T,T) for any scalar
  * type T to be supported.
  *
  * \sa Eigen::igammac(), Eigen::lgamma()
  */
template<typename Derived,typename ExponentDerived>
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)
{
  return Eigen::CwiseBinaryOp<Eigen::internal::scalar_igamma_op<typename Derived::Scalar>, const Derived, const ExponentDerived>(
    a.derived(),
    x.derived()
  );
}

/** \cpp11 \returns an expression of the coefficient-wise igamma_der_a(\a a, \a x) to the given arrays.
  *
  * This function computes the coefficient-wise derivative of the incomplete
  * gamma function with respect to the parameter a.
  *
  * \note This function supports only float and double scalar types in c++11
  * mode. To support other scalar types,
  * or float/double in non c++11 mode, the user has to provide implementations
  * of igamma_der_a(T,T) for any scalar
  * type T to be supported.
  *
  * \sa Eigen::igamma(), Eigen::lgamma()
  */
template <typename Derived, typename ExponentDerived>
inline const Eigen::CwiseBinaryOp<Eigen::internal::scalar_igamma_der_a_op<typename Derived::Scalar>, const Derived, const ExponentDerived>
igamma_der_a(const Eigen::ArrayBase<Derived>& a, const Eigen::ArrayBase<ExponentDerived>& x) {
  return Eigen::CwiseBinaryOp<Eigen::internal::scalar_igamma_der_a_op<typename Derived::Scalar>, const Derived, const ExponentDerived>(
    a.derived(),
    x.derived());
}

/** \cpp11 \returns an expression of the coefficient-wise gamma_sample_der_alpha(\a alpha, \a sample) to the given arrays.
  *
  * This function computes the coefficient-wise derivative of the sample
  * of a Gamma(alpha, 1) random variable with respect to the parameter alpha.
  *
  * \note This function supports only float and double scalar types in c++11
  * mode. To support other scalar types,
  * or float/double in non c++11 mode, the user has to provide implementations
  * of gamma_sample_der_alpha(T,T) for any scalar
  * type T to be supported.
  *
  * \sa Eigen::igamma(), Eigen::lgamma()
  */
template <typename AlphaDerived, typename SampleDerived>
inline const Eigen::CwiseBinaryOp<Eigen::internal::scalar_gamma_sample_der_alpha_op<typename AlphaDerived::Scalar>, const AlphaDerived, const SampleDerived>
gamma_sample_der_alpha(const Eigen::ArrayBase<AlphaDerived>& alpha, const Eigen::ArrayBase<SampleDerived>& sample) {
  return Eigen::CwiseBinaryOp<Eigen::internal::scalar_gamma_sample_der_alpha_op<typename AlphaDerived::Scalar>, const AlphaDerived, const SampleDerived>(
      alpha.derived(),
      sample.derived());
}

/** \cpp11 \returns an expression of the coefficient-wise igammac(\a a, \a x) to the given arrays.
  *
  * This function computes the coefficient-wise complementary incomplete gamma function.
  *
  * \note This function supports only float and double scalar types in c++11 mode. To support other scalar types,
  * or float/double in non c++11 mode, the user has to provide implementations of igammac(T,T) for any scalar
  * type T to be supported.
  *
  * \sa Eigen::igamma(), Eigen::lgamma()
  */
template<typename Derived,typename ExponentDerived>
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)
{
  return Eigen::CwiseBinaryOp<Eigen::internal::scalar_igammac_op<typename Derived::Scalar>, const Derived, const ExponentDerived>(
    a.derived(),
    x.derived()
  );
}

/** \cpp11 \returns an expression of the coefficient-wise polygamma(\a n, \a x) to the given arrays.
  *
  * It returns the \a n -th derivative of the digamma(psi) evaluated at \c x.
  *
  * \note This function supports only float and double scalar types in c++11 mode. To support other scalar types,
  * or float/double in non c++11 mode, the user has to provide implementations of polygamma(T,T) for any scalar
  * type T to be supported.
  *
  * \sa Eigen::digamma()
  */
// * \warning Be careful with the order of the parameters: x.polygamma(n) is equivalent to polygamma(n,x)
// * \sa ArrayBase::polygamma()
template<typename DerivedN,typename DerivedX>
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)
{
  return Eigen::CwiseBinaryOp<Eigen::internal::scalar_polygamma_op<typename DerivedX::Scalar>, const DerivedN, const DerivedX>(
    n.derived(),
    x.derived()
  );
}

/** \cpp11 \returns an expression of the coefficient-wise betainc(\a x, \a a, \a b) to the given arrays.
  *
  * This function computes the regularized incomplete beta function (integral).
  *
  * \note This function supports only float and double scalar types in c++11 mode. To support other scalar types,
  * or float/double in non c++11 mode, the user has to provide implementations of betainc(T,T,T) for any scalar
  * type T to be supported.
  *
  * \sa Eigen::betainc(), Eigen::lgamma()
  */
template<typename ArgADerived, typename ArgBDerived, typename ArgXDerived>
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)
{
  return Eigen::CwiseTernaryOp<Eigen::internal::scalar_betainc_op<typename ArgXDerived::Scalar>, const ArgADerived, const ArgBDerived, const ArgXDerived>(
    a.derived(),
    b.derived(),
    x.derived()
  );
}


/** \returns an expression of the coefficient-wise zeta(\a x, \a q) to the given arrays.
  *
  * It returns the Riemann zeta function of two arguments \a x and \a q:
  *
  * \param x is the exposent, it must be > 1
  * \param q is the shift, it must be > 0
  *
  * \note This function supports only float and double scalar types. To support other scalar types, the user has
  * to provide implementations of zeta(T,T) for any scalar type T to be supported.
  *
  * \sa ArrayBase::zeta()
  */
template<typename DerivedX,typename DerivedQ>
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)
{
  return Eigen::CwiseBinaryOp<Eigen::internal::scalar_zeta_op<typename DerivedX::Scalar>, const DerivedX, const DerivedQ>(
    x.derived(),
    q.derived()
  );
}

/** \returns an expression of the coefficient-wise i0e(\a x) to the given
 * arrays.
  *
  * It returns the exponentially scaled modified Bessel
  * function of order zero.
  *
  * \param x is the argument
  *
  * \note This function supports only float and double scalar types. To support
  * other scalar types, the user has to provide implementations of i0e(T) for
  * any scalar type T to be supported.
  *
  * \sa ArrayBase::i0e()
  */
template <typename Derived>
inline const Eigen::CwiseUnaryOp<
    Eigen::internal::scalar_i0e_op<typename Derived::Scalar>, const Derived>
i0e(const Eigen::ArrayBase<Derived>& x) {
  return Eigen::CwiseUnaryOp<
      Eigen::internal::scalar_i0e_op<typename Derived::Scalar>,
      const Derived>(x.derived());
}

/** \returns an expression of the coefficient-wise i1e(\a x) to the given
 * arrays.
  *
  * It returns the exponentially scaled modified Bessel
  * function of order one.
  *
  * \param x is the argument
  *
  * \note This function supports only float and double scalar types. To support
  * other scalar types, the user has to provide implementations of i1e(T) for
  * any scalar type T to be supported.
  *
  * \sa ArrayBase::i1e()
  */
template <typename Derived>
inline const Eigen::CwiseUnaryOp<
    Eigen::internal::scalar_i1e_op<typename Derived::Scalar>, const Derived>
i1e(const Eigen::ArrayBase<Derived>& x) {
  return Eigen::CwiseUnaryOp<
      Eigen::internal::scalar_i1e_op<typename Derived::Scalar>,
      const Derived>(x.derived());
}

} // end namespace Eigen

#endif // EIGEN_SPECIALFUNCTIONS_ARRAYAPI_H