File size: 10,015 Bytes
a42735d |
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 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 |
// 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_BESSELFUNCTIONS_ARRAYAPI_H
#define EIGEN_BESSELFUNCTIONS_ARRAYAPI_H
namespace Eigen {
/** \returns an expression of the coefficient-wise i0(\a x) to the given
* arrays.
*
* It returns the modified Bessel function of the first kind 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 i0(T) for
* any scalar type T to be supported.
*
* \sa ArrayBase::bessel_i0()
*/
template <typename Derived>
EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_i0_op<typename Derived::Scalar>, const Derived>
bessel_i0(const Eigen::ArrayBase<Derived>& x) {
return Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_i0_op<typename Derived::Scalar>,
const Derived>(x.derived());
}
/** \returns an expression of the coefficient-wise i0e(\a x) to the given
* arrays.
*
* It returns the exponentially scaled modified Bessel
* function of the first kind 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::bessel_i0e()
*/
template <typename Derived>
EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_i0e_op<typename Derived::Scalar>, const Derived>
bessel_i0e(const Eigen::ArrayBase<Derived>& x) {
return Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_i0e_op<typename Derived::Scalar>,
const Derived>(x.derived());
}
/** \returns an expression of the coefficient-wise i1(\a x) to the given
* arrays.
*
* It returns the modified Bessel function of the first kind 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 i1(T) for
* any scalar type T to be supported.
*
* \sa ArrayBase::bessel_i1()
*/
template <typename Derived>
EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_i1_op<typename Derived::Scalar>, const Derived>
bessel_i1(const Eigen::ArrayBase<Derived>& x) {
return Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_i1_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 the first kind 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::bessel_i1e()
*/
template <typename Derived>
EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_i1e_op<typename Derived::Scalar>, const Derived>
bessel_i1e(const Eigen::ArrayBase<Derived>& x) {
return Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_i1e_op<typename Derived::Scalar>,
const Derived>(x.derived());
}
/** \returns an expression of the coefficient-wise k0(\a x) to the given
* arrays.
*
* It returns the modified Bessel function of the second kind 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 k0(T) for
* any scalar type T to be supported.
*
* \sa ArrayBase::bessel_k0()
*/
template <typename Derived>
EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_k0_op<typename Derived::Scalar>, const Derived>
bessel_k0(const Eigen::ArrayBase<Derived>& x) {
return Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_k0_op<typename Derived::Scalar>,
const Derived>(x.derived());
}
/** \returns an expression of the coefficient-wise k0e(\a x) to the given
* arrays.
*
* It returns the exponentially scaled modified Bessel
* function of the second kind 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 k0e(T) for
* any scalar type T to be supported.
*
* \sa ArrayBase::bessel_k0e()
*/
template <typename Derived>
EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_k0e_op<typename Derived::Scalar>, const Derived>
bessel_k0e(const Eigen::ArrayBase<Derived>& x) {
return Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_k0e_op<typename Derived::Scalar>,
const Derived>(x.derived());
}
/** \returns an expression of the coefficient-wise k1(\a x) to the given
* arrays.
*
* It returns the modified Bessel function of the second kind 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 k1(T) for
* any scalar type T to be supported.
*
* \sa ArrayBase::bessel_k1()
*/
template <typename Derived>
EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_k1_op<typename Derived::Scalar>, const Derived>
bessel_k1(const Eigen::ArrayBase<Derived>& x) {
return Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_k1_op<typename Derived::Scalar>,
const Derived>(x.derived());
}
/** \returns an expression of the coefficient-wise k1e(\a x) to the given
* arrays.
*
* It returns the exponentially scaled modified Bessel
* function of the second kind 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 k1e(T) for
* any scalar type T to be supported.
*
* \sa ArrayBase::bessel_k1e()
*/
template <typename Derived>
EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_k1e_op<typename Derived::Scalar>, const Derived>
bessel_k1e(const Eigen::ArrayBase<Derived>& x) {
return Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_k1e_op<typename Derived::Scalar>,
const Derived>(x.derived());
}
/** \returns an expression of the coefficient-wise j0(\a x) to the given
* arrays.
*
* It returns the Bessel function of the first kind 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 j0(T) for
* any scalar type T to be supported.
*
* \sa ArrayBase::bessel_j0()
*/
template <typename Derived>
EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_j0_op<typename Derived::Scalar>, const Derived>
bessel_j0(const Eigen::ArrayBase<Derived>& x) {
return Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_j0_op<typename Derived::Scalar>,
const Derived>(x.derived());
}
/** \returns an expression of the coefficient-wise y0(\a x) to the given
* arrays.
*
* It returns the Bessel function of the second kind 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 y0(T) for
* any scalar type T to be supported.
*
* \sa ArrayBase::bessel_y0()
*/
template <typename Derived>
EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_y0_op<typename Derived::Scalar>, const Derived>
bessel_y0(const Eigen::ArrayBase<Derived>& x) {
return Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_y0_op<typename Derived::Scalar>,
const Derived>(x.derived());
}
/** \returns an expression of the coefficient-wise j1(\a x) to the given
* arrays.
*
* It returns the modified Bessel function of the first kind 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 j1(T) for
* any scalar type T to be supported.
*
* \sa ArrayBase::bessel_j1()
*/
template <typename Derived>
EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_j1_op<typename Derived::Scalar>, const Derived>
bessel_j1(const Eigen::ArrayBase<Derived>& x) {
return Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_j1_op<typename Derived::Scalar>,
const Derived>(x.derived());
}
/** \returns an expression of the coefficient-wise y1(\a x) to the given
* arrays.
*
* It returns the Bessel function of the second kind 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 y1(T) for
* any scalar type T to be supported.
*
* \sa ArrayBase::bessel_y1()
*/
template <typename Derived>
EIGEN_STRONG_INLINE const Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_y1_op<typename Derived::Scalar>, const Derived>
bessel_y1(const Eigen::ArrayBase<Derived>& x) {
return Eigen::CwiseUnaryOp<
Eigen::internal::scalar_bessel_y1_op<typename Derived::Scalar>,
const Derived>(x.derived());
}
} // end namespace Eigen
#endif // EIGEN_BESSELFUNCTIONS_ARRAYAPI_H
|