gcc/libstdc++-v3/include/tr1/cmath

508 lines
16 KiB
C++

// TR1 cmath -*- C++ -*-
// Copyright (C) 2006, 2007, 2008, 2009 Free Software Foundation, Inc.
//
// This file is part of the GNU ISO C++ Library. This library is free
// software; you can redistribute it and/or modify it under the
// terms of the GNU General Public License as published by the
// Free Software Foundation; either version 3, or (at your option)
// any later version.
// This library is distributed in the hope that it will be useful,
// but WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
// GNU General Public License for more details.
// Under Section 7 of GPL version 3, you are granted additional
// permissions described in the GCC Runtime Library Exception, version
// 3.1, as published by the Free Software Foundation.
// You should have received a copy of the GNU General Public License and
// a copy of the GCC Runtime Library Exception along with this program;
// see the files COPYING3 and COPYING.RUNTIME respectively. If not, see
// <http://www.gnu.org/licenses/>.
/** @file tr1/cmath
* This is a TR1 C++ Library header.
*/
#ifndef _GLIBCXX_TR1_CMATH
#define _GLIBCXX_TR1_CMATH 1
#pragma GCC system_header
#if defined(_GLIBCXX_INCLUDE_AS_CXX0X)
# error TR1 header cannot be included from C++0x header
#endif
#include <cmath>
#if defined(_GLIBCXX_INCLUDE_AS_TR1)
# include <tr1_impl/cmath>
#else
# define _GLIBCXX_INCLUDE_AS_TR1
# define _GLIBCXX_BEGIN_NAMESPACE_TR1 namespace tr1 {
# define _GLIBCXX_END_NAMESPACE_TR1 }
# define _GLIBCXX_TR1 tr1::
# include <tr1_impl/cmath>
# undef _GLIBCXX_TR1
# undef _GLIBCXX_END_NAMESPACE_TR1
# undef _GLIBCXX_BEGIN_NAMESPACE_TR1
# undef _GLIBCXX_INCLUDE_AS_TR1
#endif
namespace std
{
namespace tr1
{
// DR 550. What should the return type of pow(float,int) be?
// NB: C++0x and TR1 != C++03.
inline double
pow(double __x, double __y)
{ return std::pow(__x, __y); }
inline float
pow(float __x, float __y)
{ return std::pow(__x, __y); }
inline long double
pow(long double __x, long double __y)
{ return std::pow(__x, __y); }
template<typename _Tp, typename _Up>
inline typename __gnu_cxx::__promote_2<_Tp, _Up>::__type
pow(_Tp __x, _Up __y)
{
typedef typename __gnu_cxx::__promote_2<_Tp, _Up>::__type __type;
return std::pow(__type(__x), __type(__y));
}
}
}
#include <bits/stl_algobase.h>
#include <limits>
#include <tr1/type_traits>
#include <tr1/gamma.tcc>
#include <tr1/bessel_function.tcc>
#include <tr1/beta_function.tcc>
#include <tr1/ell_integral.tcc>
#include <tr1/exp_integral.tcc>
#include <tr1/hypergeometric.tcc>
#include <tr1/legendre_function.tcc>
#include <tr1/modified_bessel_func.tcc>
#include <tr1/poly_hermite.tcc>
#include <tr1/poly_laguerre.tcc>
#include <tr1/riemann_zeta.tcc>
namespace std
{
namespace tr1
{
/**
* @defgroup tr1_math_spec_func Mathematical Special Functions
* @ingroup numerics
*
* A collection of advanced mathematical special functions.
* @{
*/
inline float
assoc_laguerref(unsigned int __n, unsigned int __m, float __x)
{ return __detail::__assoc_laguerre<float>(__n, __m, __x); }
inline long double
assoc_laguerrel(unsigned int __n, unsigned int __m, long double __x)
{
return __detail::__assoc_laguerre<long double>(__n, __m, __x);
}
/// 5.2.1.1 Associated Laguerre polynomials.
template<typename _Tp>
inline typename __gnu_cxx::__promote<_Tp>::__type
assoc_laguerre(unsigned int __n, unsigned int __m, _Tp __x)
{
typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
return __detail::__assoc_laguerre<__type>(__n, __m, __x);
}
inline float
assoc_legendref(unsigned int __l, unsigned int __m, float __x)
{ return __detail::__assoc_legendre_p<float>(__l, __m, __x); }
inline long double
assoc_legendrel(unsigned int __l, unsigned int __m, long double __x)
{ return __detail::__assoc_legendre_p<long double>(__l, __m, __x); }
/// 5.2.1.2 Associated Legendre functions.
template<typename _Tp>
inline typename __gnu_cxx::__promote<_Tp>::__type
assoc_legendre(unsigned int __l, unsigned int __m, _Tp __x)
{
typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
return __detail::__assoc_legendre_p<__type>(__l, __m, __x);
}
inline float
betaf(float __x, float __y)
{ return __detail::__beta<float>(__x, __y); }
inline long double
betal(long double __x, long double __y)
{ return __detail::__beta<long double>(__x, __y); }
/// 5.2.1.3 Beta functions.
template<typename _Tpx, typename _Tpy>
inline typename __gnu_cxx::__promote_2<_Tpx, _Tpy>::__type
beta(_Tpx __x, _Tpy __y)
{
typedef typename __gnu_cxx::__promote_2<_Tpx, _Tpy>::__type __type;
return __detail::__beta<__type>(__x, __y);
}
inline float
comp_ellint_1f(float __k)
{ return __detail::__comp_ellint_1<float>(__k); }
inline long double
comp_ellint_1l(long double __k)
{ return __detail::__comp_ellint_1<long double>(__k); }
/// 5.2.1.4 Complete elliptic integrals of the first kind.
template<typename _Tp>
inline typename __gnu_cxx::__promote<_Tp>::__type
comp_ellint_1(_Tp __k)
{
typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
return __detail::__comp_ellint_1<__type>(__k);
}
inline float
comp_ellint_2f(float __k)
{ return __detail::__comp_ellint_2<float>(__k); }
inline long double
comp_ellint_2l(long double __k)
{ return __detail::__comp_ellint_2<long double>(__k); }
/// 5.2.1.5 Complete elliptic integrals of the second kind.
template<typename _Tp>
inline typename __gnu_cxx::__promote<_Tp>::__type
comp_ellint_2(_Tp __k)
{
typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
return __detail::__comp_ellint_2<__type>(__k);
}
inline float
comp_ellint_3f(float __k, float __nu)
{ return __detail::__comp_ellint_3<float>(__k, __nu); }
inline long double
comp_ellint_3l(long double __k, long double __nu)
{ return __detail::__comp_ellint_3<long double>(__k, __nu); }
/// 5.2.1.6 Complete elliptic integrals of the third kind.
template<typename _Tp, typename _Tpn>
inline typename __gnu_cxx::__promote_2<_Tp, _Tpn>::__type
comp_ellint_3(_Tp __k, _Tpn __nu)
{
typedef typename __gnu_cxx::__promote_2<_Tp, _Tpn>::__type __type;
return __detail::__comp_ellint_3<__type>(__k, __nu);
}
inline float
conf_hypergf(float __a, float __c, float __x)
{ return __detail::__conf_hyperg<float>(__a, __c, __x); }
inline long double
conf_hypergl(long double __a, long double __c, long double __x)
{ return __detail::__conf_hyperg<long double>(__a, __c, __x); }
/// 5.2.1.7 Confluent hypergeometric functions.
template<typename _Tpa, typename _Tpc, typename _Tp>
inline typename __gnu_cxx::__promote_3<_Tpa, _Tpc, _Tp>::__type
conf_hyperg(_Tpa __a, _Tpc __c, _Tp __x)
{
typedef typename __gnu_cxx::__promote_3<_Tpa, _Tpc, _Tp>::__type __type;
return __detail::__conf_hyperg<__type>(__a, __c, __x);
}
inline float
cyl_bessel_if(float __nu, float __x)
{ return __detail::__cyl_bessel_i<float>(__nu, __x); }
inline long double
cyl_bessel_il(long double __nu, long double __x)
{ return __detail::__cyl_bessel_i<long double>(__nu, __x); }
/// 5.2.1.8 Regular modified cylindrical Bessel functions.
template<typename _Tpnu, typename _Tp>
inline typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type
cyl_bessel_i(_Tpnu __nu, _Tp __x)
{
typedef typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type __type;
return __detail::__cyl_bessel_i<__type>(__nu, __x);
}
inline float
cyl_bessel_jf(float __nu, float __x)
{ return __detail::__cyl_bessel_j<float>(__nu, __x); }
inline long double
cyl_bessel_jl(long double __nu, long double __x)
{ return __detail::__cyl_bessel_j<long double>(__nu, __x); }
/// 5.2.1.9 Cylindrical Bessel functions (of the first kind).
template<typename _Tpnu, typename _Tp>
inline typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type
cyl_bessel_j(_Tpnu __nu, _Tp __x)
{
typedef typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type __type;
return __detail::__cyl_bessel_j<__type>(__nu, __x);
}
inline float
cyl_bessel_kf(float __nu, float __x)
{ return __detail::__cyl_bessel_k<float>(__nu, __x); }
inline long double
cyl_bessel_kl(long double __nu, long double __x)
{ return __detail::__cyl_bessel_k<long double>(__nu, __x); }
/// 5.2.1.10 Irregular modified cylindrical Bessel functions.
template<typename _Tpnu, typename _Tp>
inline typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type
cyl_bessel_k(_Tpnu __nu, _Tp __x)
{
typedef typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type __type;
return __detail::__cyl_bessel_k<__type>(__nu, __x);
}
inline float
cyl_neumannf(float __nu, float __x)
{ return __detail::__cyl_neumann_n<float>(__nu, __x); }
inline long double
cyl_neumannl(long double __nu, long double __x)
{ return __detail::__cyl_neumann_n<long double>(__nu, __x); }
/// 5.2.1.11 Cylindrical Neumann functions.
template<typename _Tpnu, typename _Tp>
inline typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type
cyl_neumann(_Tpnu __nu, _Tp __x)
{
typedef typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type __type;
return __detail::__cyl_neumann_n<__type>(__nu, __x);
}
inline float
ellint_1f(float __k, float __phi)
{ return __detail::__ellint_1<float>(__k, __phi); }
inline long double
ellint_1l(long double __k, long double __phi)
{ return __detail::__ellint_1<long double>(__k, __phi); }
/// 5.2.1.12 Incomplete elliptic integrals of the first kind.
template<typename _Tp, typename _Tpp>
inline typename __gnu_cxx::__promote_2<_Tp, _Tpp>::__type
ellint_1(_Tp __k, _Tpp __phi)
{
typedef typename __gnu_cxx::__promote_2<_Tp, _Tpp>::__type __type;
return __detail::__ellint_1<__type>(__k, __phi);
}
inline float
ellint_2f(float __k, float __phi)
{ return __detail::__ellint_2<float>(__k, __phi); }
inline long double
ellint_2l(long double __k, long double __phi)
{ return __detail::__ellint_2<long double>(__k, __phi); }
/// 5.2.1.13 Incomplete elliptic integrals of the second kind.
template<typename _Tp, typename _Tpp>
inline typename __gnu_cxx::__promote_2<_Tp, _Tpp>::__type
ellint_2(_Tp __k, _Tpp __phi)
{
typedef typename __gnu_cxx::__promote_2<_Tp, _Tpp>::__type __type;
return __detail::__ellint_2<__type>(__k, __phi);
}
inline float
ellint_3f(float __k, float __nu, float __phi)
{ return __detail::__ellint_3<float>(__k, __nu, __phi); }
inline long double
ellint_3l(long double __k, long double __nu, long double __phi)
{ return __detail::__ellint_3<long double>(__k, __nu, __phi); }
/// 5.2.1.14 Incomplete elliptic integrals of the third kind.
template<typename _Tp, typename _Tpn, typename _Tpp>
inline typename __gnu_cxx::__promote_3<_Tp, _Tpn, _Tpp>::__type
ellint_3(_Tp __k, _Tpn __nu, _Tpp __phi)
{
typedef typename __gnu_cxx::__promote_3<_Tp, _Tpn, _Tpp>::__type __type;
return __detail::__ellint_3<__type>(__k, __nu, __phi);
}
inline float
expintf(float __x)
{ return __detail::__expint<float>(__x); }
inline long double
expintl(long double __x)
{ return __detail::__expint<long double>(__x); }
/// 5.2.1.15 Exponential integrals.
template<typename _Tp>
inline typename __gnu_cxx::__promote<_Tp>::__type
expint(_Tp __x)
{
typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
return __detail::__expint<__type>(__x);
}
inline float
hermitef(unsigned int __n, float __x)
{ return __detail::__poly_hermite<float>(__n, __x); }
inline long double
hermitel(unsigned int __n, long double __x)
{ return __detail::__poly_hermite<long double>(__n, __x); }
/// 5.2.1.16 Hermite polynomials.
template<typename _Tp>
inline typename __gnu_cxx::__promote<_Tp>::__type
hermite(unsigned int __n, _Tp __x)
{
typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
return __detail::__poly_hermite<__type>(__n, __x);
}
inline float
hypergf(float __a, float __b, float __c, float __x)
{ return __detail::__hyperg<float>(__a, __b, __c, __x); }
inline long double
hypergl(long double __a, long double __b, long double __c, long double __x)
{ return __detail::__hyperg<long double>(__a, __b, __c, __x); }
/// 5.2.1.17 Hypergeometric functions.
template<typename _Tpa, typename _Tpb, typename _Tpc, typename _Tp>
inline typename __gnu_cxx::__promote_4<_Tpa, _Tpb, _Tpc, _Tp>::__type
hyperg(_Tpa __a, _Tpb __b, _Tpc __c, _Tp __x)
{
typedef typename __gnu_cxx::__promote_4<_Tpa, _Tpb, _Tpc, _Tp>::__type __type;
return __detail::__hyperg<__type>(__a, __b, __c, __x);
}
inline float
laguerref(unsigned int __n, float __x)
{ return __detail::__laguerre<float>(__n, __x); }
inline long double
laguerrel(unsigned int __n, long double __x)
{ return __detail::__laguerre<long double>(__n, __x); }
/// 5.2.1.18 Laguerre polynomials.
template<typename _Tp>
inline typename __gnu_cxx::__promote<_Tp>::__type
laguerre(unsigned int __n, _Tp __x)
{
typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
return __detail::__laguerre<__type>(__n, __x);
}
inline float
legendref(unsigned int __n, float __x)
{ return __detail::__poly_legendre_p<float>(__n, __x); }
inline long double
legendrel(unsigned int __n, long double __x)
{ return __detail::__poly_legendre_p<long double>(__n, __x); }
/// 5.2.1.19 Legendre polynomials.
template<typename _Tp>
inline typename __gnu_cxx::__promote<_Tp>::__type
legendre(unsigned int __n, _Tp __x)
{
typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
return __detail::__poly_legendre_p<__type>(__n, __x);
}
inline float
riemann_zetaf(float __x)
{ return __detail::__riemann_zeta<float>(__x); }
inline long double
riemann_zetal(long double __x)
{ return __detail::__riemann_zeta<long double>(__x); }
/// 5.2.1.20 Riemann zeta function.
template<typename _Tp>
inline typename __gnu_cxx::__promote<_Tp>::__type
riemann_zeta(_Tp __x)
{
typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
return __detail::__riemann_zeta<__type>(__x);
}
inline float
sph_besself(unsigned int __n, float __x)
{ return __detail::__sph_bessel<float>(__n, __x); }
inline long double
sph_bessell(unsigned int __n, long double __x)
{ return __detail::__sph_bessel<long double>(__n, __x); }
/// 5.2.1.21 Spherical Bessel functions.
template<typename _Tp>
inline typename __gnu_cxx::__promote<_Tp>::__type
sph_bessel(unsigned int __n, _Tp __x)
{
typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
return __detail::__sph_bessel<__type>(__n, __x);
}
inline float
sph_legendref(unsigned int __l, unsigned int __m, float __theta)
{ return __detail::__sph_legendre<float>(__l, __m, __theta); }
inline long double
sph_legendrel(unsigned int __l, unsigned int __m, long double __theta)
{ return __detail::__sph_legendre<long double>(__l, __m, __theta); }
/// 5.2.1.22 Spherical associated Legendre functions.
template<typename _Tp>
inline typename __gnu_cxx::__promote<_Tp>::__type
sph_legendre(unsigned int __l, unsigned int __m, _Tp __theta)
{
typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
return __detail::__sph_legendre<__type>(__l, __m, __theta);
}
inline float
sph_neumannf(unsigned int __n, float __x)
{ return __detail::__sph_neumann<float>(__n, __x); }
inline long double
sph_neumannl(unsigned int __n, long double __x)
{ return __detail::__sph_neumann<long double>(__n, __x); }
/// 5.2.1.23 Spherical Neumann functions.
template<typename _Tp>
inline typename __gnu_cxx::__promote<_Tp>::__type
sph_neumann(unsigned int __n, _Tp __x)
{
typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
return __detail::__sph_neumann<__type>(__n, __x);
}
/* @} */ // tr1_math_spec_func
}
}
#endif // _GLIBCXX_TR1_CMATH