blob: 41d2bf61c1fa3c44658db01148cbd8d6bc97ff2e [file] [log] [blame]
Austin Schuhc55b0172022-02-20 17:52:35 -08001// 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_BESSELFUNCTIONS_ARRAYAPI_H
12#define EIGEN_BESSELFUNCTIONS_ARRAYAPI_H
13
14namespace Eigen {
15
16/** \returns an expression of the coefficient-wise i0(\a x) to the given
17 * arrays.
18 *
19 * It returns the modified Bessel function of the first kind of order zero.
20 *
21 * \param x is the argument
22 *
23 * \note This function supports only float and double scalar types. To support
24 * other scalar types, the user has to provide implementations of i0(T) for
25 * any scalar type T to be supported.
26 *
27 * \sa ArrayBase::bessel_i0()
28 */
29template <typename Derived>
30EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
31 Eigen::internal::scalar_bessel_i0_op<typename Derived::Scalar>, const Derived>
32bessel_i0(const Eigen::ArrayBase<Derived>& x) {
33 return Eigen::CwiseUnaryOp<
34 Eigen::internal::scalar_bessel_i0_op<typename Derived::Scalar>,
35 const Derived>(x.derived());
36}
37
38/** \returns an expression of the coefficient-wise i0e(\a x) to the given
39 * arrays.
40 *
41 * It returns the exponentially scaled modified Bessel
42 * function of the first kind of order zero.
43 *
44 * \param x is the argument
45 *
46 * \note This function supports only float and double scalar types. To support
47 * other scalar types, the user has to provide implementations of i0e(T) for
48 * any scalar type T to be supported.
49 *
50 * \sa ArrayBase::bessel_i0e()
51 */
52template <typename Derived>
53EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
54 Eigen::internal::scalar_bessel_i0e_op<typename Derived::Scalar>, const Derived>
55bessel_i0e(const Eigen::ArrayBase<Derived>& x) {
56 return Eigen::CwiseUnaryOp<
57 Eigen::internal::scalar_bessel_i0e_op<typename Derived::Scalar>,
58 const Derived>(x.derived());
59}
60
61/** \returns an expression of the coefficient-wise i1(\a x) to the given
62 * arrays.
63 *
64 * It returns the modified Bessel function of the first kind of order one.
65 *
66 * \param x is the argument
67 *
68 * \note This function supports only float and double scalar types. To support
69 * other scalar types, the user has to provide implementations of i1(T) for
70 * any scalar type T to be supported.
71 *
72 * \sa ArrayBase::bessel_i1()
73 */
74template <typename Derived>
75EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
76 Eigen::internal::scalar_bessel_i1_op<typename Derived::Scalar>, const Derived>
77bessel_i1(const Eigen::ArrayBase<Derived>& x) {
78 return Eigen::CwiseUnaryOp<
79 Eigen::internal::scalar_bessel_i1_op<typename Derived::Scalar>,
80 const Derived>(x.derived());
81}
82
83/** \returns an expression of the coefficient-wise i1e(\a x) to the given
84 * arrays.
85 *
86 * It returns the exponentially scaled modified Bessel
87 * function of the first kind of order one.
88 *
89 * \param x is the argument
90 *
91 * \note This function supports only float and double scalar types. To support
92 * other scalar types, the user has to provide implementations of i1e(T) for
93 * any scalar type T to be supported.
94 *
95 * \sa ArrayBase::bessel_i1e()
96 */
97template <typename Derived>
98EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
99 Eigen::internal::scalar_bessel_i1e_op<typename Derived::Scalar>, const Derived>
100bessel_i1e(const Eigen::ArrayBase<Derived>& x) {
101 return Eigen::CwiseUnaryOp<
102 Eigen::internal::scalar_bessel_i1e_op<typename Derived::Scalar>,
103 const Derived>(x.derived());
104}
105
106/** \returns an expression of the coefficient-wise k0(\a x) to the given
107 * arrays.
108 *
109 * It returns the modified Bessel function of the second kind of order zero.
110 *
111 * \param x is the argument
112 *
113 * \note This function supports only float and double scalar types. To support
114 * other scalar types, the user has to provide implementations of k0(T) for
115 * any scalar type T to be supported.
116 *
117 * \sa ArrayBase::bessel_k0()
118 */
119template <typename Derived>
120EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
121 Eigen::internal::scalar_bessel_k0_op<typename Derived::Scalar>, const Derived>
122bessel_k0(const Eigen::ArrayBase<Derived>& x) {
123 return Eigen::CwiseUnaryOp<
124 Eigen::internal::scalar_bessel_k0_op<typename Derived::Scalar>,
125 const Derived>(x.derived());
126}
127
128/** \returns an expression of the coefficient-wise k0e(\a x) to the given
129 * arrays.
130 *
131 * It returns the exponentially scaled modified Bessel
132 * function of the second kind of order zero.
133 *
134 * \param x is the argument
135 *
136 * \note This function supports only float and double scalar types. To support
137 * other scalar types, the user has to provide implementations of k0e(T) for
138 * any scalar type T to be supported.
139 *
140 * \sa ArrayBase::bessel_k0e()
141 */
142template <typename Derived>
143EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
144 Eigen::internal::scalar_bessel_k0e_op<typename Derived::Scalar>, const Derived>
145bessel_k0e(const Eigen::ArrayBase<Derived>& x) {
146 return Eigen::CwiseUnaryOp<
147 Eigen::internal::scalar_bessel_k0e_op<typename Derived::Scalar>,
148 const Derived>(x.derived());
149}
150
151/** \returns an expression of the coefficient-wise k1(\a x) to the given
152 * arrays.
153 *
154 * It returns the modified Bessel function of the second kind of order one.
155 *
156 * \param x is the argument
157 *
158 * \note This function supports only float and double scalar types. To support
159 * other scalar types, the user has to provide implementations of k1(T) for
160 * any scalar type T to be supported.
161 *
162 * \sa ArrayBase::bessel_k1()
163 */
164template <typename Derived>
165EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
166 Eigen::internal::scalar_bessel_k1_op<typename Derived::Scalar>, const Derived>
167bessel_k1(const Eigen::ArrayBase<Derived>& x) {
168 return Eigen::CwiseUnaryOp<
169 Eigen::internal::scalar_bessel_k1_op<typename Derived::Scalar>,
170 const Derived>(x.derived());
171}
172
173/** \returns an expression of the coefficient-wise k1e(\a x) to the given
174 * arrays.
175 *
176 * It returns the exponentially scaled modified Bessel
177 * function of the second kind of order one.
178 *
179 * \param x is the argument
180 *
181 * \note This function supports only float and double scalar types. To support
182 * other scalar types, the user has to provide implementations of k1e(T) for
183 * any scalar type T to be supported.
184 *
185 * \sa ArrayBase::bessel_k1e()
186 */
187template <typename Derived>
188EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
189 Eigen::internal::scalar_bessel_k1e_op<typename Derived::Scalar>, const Derived>
190bessel_k1e(const Eigen::ArrayBase<Derived>& x) {
191 return Eigen::CwiseUnaryOp<
192 Eigen::internal::scalar_bessel_k1e_op<typename Derived::Scalar>,
193 const Derived>(x.derived());
194}
195
196/** \returns an expression of the coefficient-wise j0(\a x) to the given
197 * arrays.
198 *
199 * It returns the Bessel function of the first kind of order zero.
200 *
201 * \param x is the argument
202 *
203 * \note This function supports only float and double scalar types. To support
204 * other scalar types, the user has to provide implementations of j0(T) for
205 * any scalar type T to be supported.
206 *
207 * \sa ArrayBase::bessel_j0()
208 */
209template <typename Derived>
210EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
211 Eigen::internal::scalar_bessel_j0_op<typename Derived::Scalar>, const Derived>
212bessel_j0(const Eigen::ArrayBase<Derived>& x) {
213 return Eigen::CwiseUnaryOp<
214 Eigen::internal::scalar_bessel_j0_op<typename Derived::Scalar>,
215 const Derived>(x.derived());
216}
217
218/** \returns an expression of the coefficient-wise y0(\a x) to the given
219 * arrays.
220 *
221 * It returns the Bessel function of the second kind of order zero.
222 *
223 * \param x is the argument
224 *
225 * \note This function supports only float and double scalar types. To support
226 * other scalar types, the user has to provide implementations of y0(T) for
227 * any scalar type T to be supported.
228 *
229 * \sa ArrayBase::bessel_y0()
230 */
231template <typename Derived>
232EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
233 Eigen::internal::scalar_bessel_y0_op<typename Derived::Scalar>, const Derived>
234bessel_y0(const Eigen::ArrayBase<Derived>& x) {
235 return Eigen::CwiseUnaryOp<
236 Eigen::internal::scalar_bessel_y0_op<typename Derived::Scalar>,
237 const Derived>(x.derived());
238}
239
240/** \returns an expression of the coefficient-wise j1(\a x) to the given
241 * arrays.
242 *
243 * It returns the modified Bessel function of the first kind of order one.
244 *
245 * \param x is the argument
246 *
247 * \note This function supports only float and double scalar types. To support
248 * other scalar types, the user has to provide implementations of j1(T) for
249 * any scalar type T to be supported.
250 *
251 * \sa ArrayBase::bessel_j1()
252 */
253template <typename Derived>
254EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
255 Eigen::internal::scalar_bessel_j1_op<typename Derived::Scalar>, const Derived>
256bessel_j1(const Eigen::ArrayBase<Derived>& x) {
257 return Eigen::CwiseUnaryOp<
258 Eigen::internal::scalar_bessel_j1_op<typename Derived::Scalar>,
259 const Derived>(x.derived());
260}
261
262/** \returns an expression of the coefficient-wise y1(\a x) to the given
263 * arrays.
264 *
265 * It returns the Bessel function of the second kind of order one.
266 *
267 * \param x is the argument
268 *
269 * \note This function supports only float and double scalar types. To support
270 * other scalar types, the user has to provide implementations of y1(T) for
271 * any scalar type T to be supported.
272 *
273 * \sa ArrayBase::bessel_y1()
274 */
275template <typename Derived>
276EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
277 Eigen::internal::scalar_bessel_y1_op<typename Derived::Scalar>, const Derived>
278bessel_y1(const Eigen::ArrayBase<Derived>& x) {
279 return Eigen::CwiseUnaryOp<
280 Eigen::internal::scalar_bessel_y1_op<typename Derived::Scalar>,
281 const Derived>(x.derived());
282}
283
284} // end namespace Eigen
285
286#endif // EIGEN_BESSELFUNCTIONS_ARRAYAPI_H