BornAgain  1.19.0
Simulate and fit neutron and x-ray scattering at grazing incidence
Math::Bessel Namespace Reference

Real and complex Bessel functions. More...

Functions

double I0 (double x)
 Modified Bessel function of the first kind and order 0. More...
 
complex_t J0 (const complex_t z)
 Complex Bessel function of the first kind and order 0. More...
 
double J0 (double x)
 Bessel function of the first kind and order 0. More...
 
complex_t J1 (const complex_t z)
 Complex Bessel function of the first kind and order 1. More...
 
double J1 (double x)
 Bessel function of the first kind and order 1. More...
 
complex_t J1c (const complex_t z)
 Complex Bessel function J1(x)/x. More...
 
double J1c (double x)
 Bessel function J1(x)/x. More...
 

Detailed Description

Real and complex Bessel functions.

Function Documentation

◆ I0()

double Math::Bessel::I0 ( double  x)

Modified Bessel function of the first kind and order 0.

Definition at line 177 of file Bessel.cpp.

178 {
179  return gsl_sf_bessel_I0(x);
180 }

◆ J0() [1/2]

complex_t Math::Bessel::J0 ( const complex_t  z)

Complex Bessel function of the first kind and order 0.

Definition at line 182 of file Bessel.cpp.

183 {
184  if (std::imag(z) == 0)
185  return gsl_sf_bessel_J0(std::real(z));
186  return J0_PowSer(z);
187 }

◆ J0() [2/2]

double Math::Bessel::J0 ( double  x)

Bessel function of the first kind and order 0.

Definition at line 162 of file Bessel.cpp.

163 {
164  return gsl_sf_bessel_J0(x);
165 }

Referenced by FTDistribution2DCone::evaluate(), and InterferenceFunctionHardDisk::integrand().

◆ J1() [1/2]

complex_t Math::Bessel::J1 ( const complex_t  z)

Complex Bessel function of the first kind and order 1.

Definition at line 189 of file Bessel.cpp.

190 {
191  if (std::imag(z) == 0)
192  return gsl_sf_bessel_J1(std::real(z));
193  return J1_PowSer(z);
194 }

◆ J1() [2/2]

double Math::Bessel::J1 ( double  x)

Bessel function of the first kind and order 1.

Definition at line 167 of file Bessel.cpp.

168 {
169  return gsl_sf_bessel_J1(x);
170 }

◆ J1c() [1/2]

complex_t Math::Bessel::J1c ( const complex_t  z)

Complex Bessel function J1(x)/x.

Definition at line 196 of file Bessel.cpp.

197 {
198  if (std::imag(z) == 0) {
199  double xv = std::real(z);
200  return xv == 0 ? 0.5 : gsl_sf_bessel_J1(xv) / xv;
201  }
202  return z == 0. ? 0.5 : J1_PowSer(z) / z;
203 }

◆ J1c() [2/2]