특수 함수

이 단원에서는 GSL이 제공하는 특수 함수 구현체들을 기술합니다. 라이브러리에서 에어리 함수, 베셀 함수, 클라우센 함수, 쿨룽 파동 함수, 커플링 계수, 도슨 함수, 디바이 함수, 다이로그, 타원 적분, 자코비 타원적분, 오류 함수, 지수적 적분, 페르미-디렉 함수, 감마 함수, 구겐바우어 함수, 에르미트 다항식과 함수, 초기하 함수, 라게르 함수, 르장그르 함수, 구면 조화 함수, 프사이(디감마)함수, 싱크로트론 함수, 전송 함수, 삼각함수, 제타 함수 등의 구현체들을 제공합니다.

각각의 구현체는 수치적 오류와 함수의 값을 함께 계산합니다.

기술된 함수 구현체들은 각각 다른 헤더파일에 정의되어 있습니다. 예를 들어서 gsl_sf_airy.hgsl_sf_bessel.h 로 각각 나뉘어져 있습니다. 모두 포함시켜서 계산하려면, gsl_sf.h 을 포함시키면 됩니다.

함수의 사용

함수 구현체들은 두 가지 형태로 호출할 수 있습니다. 하나는 함수의 값을 불러오는 호출입니다. 이는 주어진 인자에 대한 함수의 수치값을 반환합니다. 다른 하나는 오류 관리를 위한 호출 형태입니다. 이 호출은 함수 값과 함께 오류 값을 반환합니다. 이 두 가지 호출 형태는 모두 동일한 코드에 접근합니다. 즉, 오류 관리 호출에서 얻은 계산 값과 곧바로 계산 값을 받는 호출에서 얻은 계산 값은 정확히 동일합니다.

참고

적어도 동일한 과정을 거친 계산입니다(*).

함수 값을 불러오는 호출은 말 그대로 인자에 대해 값을 계산한 수치 값만을 반환합니다. 이는 바로 수학 수식 표현에 사용할 수 있습니다. 예를 들어서 다음의 함수 호출은 베셀 함수 \(J_0 (x)\) 의 값을 계산합니다.

double y = gsl_sf_bessel_J0(x);

이 방법은 오류 값를 확인하거나 오차 값을 계산할 수 없습니다. 해당 정보를 얻기 위해서는 오류 관리 호출 형태로 대체해서 해당 값들을 저장할 변수를 인자에 하나 더 추가해야 합니다.

gsl_sf_result result;
int status = gsl_sf_bessel_J0_e (x, &result);

이러한 오류 관리 함수는 _e 접미사를 가지고 있습니다. 반환 값은 오류의 조건을 나타냅니다. 대표적으로 오버플로우, 언더플로우 아니면 오차를 나타냅니다. 만약, 오류가 없다면 오류 관리 함수는 GSL_SUCCESS 값을 반환합니다.

gsl_sf_result 구조체

특수 함수들의 구현체들은 결과 값을 계산하는 과정에서 오차를 함께 계산합니다. 이를 위해, 함수의 오차와 결과값을 함께 다룰 수 있는 구조체를 함께 제공합니다. 이 구조체는 헤더파일 gsl_sf_result.h 정의되어 있습니다.

이 구조체는 다음과 같이 함수 값과 오차 값으로 이루어져 있습니다.

type gsl_sf_result
typedef struct
{
  double val;
  double err;
} gsl_sf_result;

val 는 함수의 값을 저장하고, err 는 함수 값의 절대 오차를 나타냅니다.

몇몇 경우에, 오버플로우나 언더플로우가 함수에 의해 감지될 수 있습니다. 그러한 경우에 크기를 나타내는 지수 값을 오차/값 쌍과 함께 반환할 수도 있습니다. 이러한 방법은 내장된 변수의 동적 크기를 초과해서 결과를 나타낼 수 있는 장점이 있습니다. 이를 위해 다음과 같이 함수 값, 오차, 그리고 지수 값을 포함하는 구조체를 사용합니다. 실제 값은 result * 10^(e10) 을 계산해야 합니다.

type gsl_sf_result_e10
typedef struct
{
  double val;
  double err;
  int    e10;
} gsl_sf_result_e10;

모드

특수 함수 구현체들의 제공 목적은 배정밀도를 가지는 특수 함수 값의 제공을 위함입니다. 그러나 몇몇 특수 함수들의 경우, 배정밀도 수준의 정밀도를 얻기 위해서 매우 높은 항이 필요한 경우도 있습니다. 이러한 경우에 gsl_mode_t 자료형으로 정의된 mode 인자를 이용해 요구 정밀도를 낮추어 연산 효율을 높일 수 있습니다. 다음의 정밀도 값들이 mode 인자에 들어갈 수 있습니다.

type gsl_mode_t
GSL_PREC_DOUBLE

배정밀도, 근사적으로 \(2 \cdot 10^{-16}\) 의 상대 정확도를 가짐.

GSL_PREC_SINGLE

단정밀도, 근사적으로 \(1 \cdot 10^{-7}\) 의 상대 정확도를 가짐.

GSL_PREC_APPROX

근사 값, 근사적으로 \(5 \cdot 10^{-4}\) 의 상대 정확도를 가짐.

근사 상태는 가장 빠른 계산 시간을 보장하지만 가장 낮은 정확도를 가집니다.

에어리 함수와 도함수(Airy function & Derivative)

에어리 함수(Airy function) \(Ai(x)\)\(Bi(x)\) 는 적분표현으로 다음과 같이 정의됩니다.

\[\begin{split}Ai(x) & = {1\over\pi} \int_0^\infty \cos(t^3/3 + xt ) \,dt \\ Bi(x) & = {1\over\pi} \int_0^\infty (e^{-t^3/3 + xt} + \sin(t^3/3 + xt)) \,dt\end{split}\]

더 자세한 정보는 Abramowitz & Stengu, Section 10.4를 참고할 수 있습니다. 에어리 함수들은 헤더 파일 gsl_sf_airy.h 에 정의되어 있습니다.

에어리 함수

double gsl_sf_airy_Ai(double x, gsl_mode_t mode)
int gsl_sf_airy_Ai_e(double x, gsl_mode_t mode, gsl_sf_result *result)

에어리 함수 \(Ai(x)\) 정밀도에 따라 계산합니다.

double gsl_sf_airy_Bi(double x, gsl_mode_t mode)
int gsl_sf_airy_Bi_e(double x, gsl_mode_t mode, gsl_sf_result *result)

에어리 함수 \(Bi(x)\) 정밀도에 따라 계산합니다.

double gsl_sf_airy_Ai_scaled(double x, gsl_mode_t mode)
int gsl_sf_airy_Ai_scaled_e(double x, gsl_mode_t mode, gsl_sf_result *result)

조정된 에어리 함수 \(S_A(x) Ai (x)\) 의 값을 계산합니다. \(x>0\) 의 값은 \(exp{(\frac{2}{3} x^\frac{3}{2})}\) 이고, \(x<0\) 일 때는 \(1\) 입니다.

double gsl_sf_airy_Bi_scaled(double x, gsl_mode_t mode)
int gsl_sf_airy_Bi_scaled_e(double x, gsl_mode_t mode, gsl_sf_result *result)

조정된 에어리 함수 \(S_A(x) Bi (x)\) 의 값을 계산합니다. \(x>0\) 일 때, 조정 계수 \(S_A(x)\) 의 값은 \(exp{(\frac{2}{3} x^\frac{3}{2})}\) 이고, \(x<0\) 일 때는 \(1\) 입니다.

에어리 함수의 도함수

double gsl_sf_airy_Ai_deriv(double x, gsl_mode_t mode)
int gsl_sf_airy_Ai_deriv_e(double x, gsl_mode_t mode, gsl_sf_result *result)

에어리 함수의 도함수 \(Ai'(x)\) 의 값을 주어진 mode 정밀도에 따라 계산합니다.

double gsl_sf_airy_Bi_deriv(double x, gsl_mode_t mode)
int gsl_sf_airy_Bi_deriv_e(double x, gsl_mode_t mode, gsl_sf_result *result)

에어리 함수의 도함수 \(Bi'(x)\) 의 값을 주어진 mode 정밀도에 따라 계산합니다.

double gsl_sf_airy_Ai_deriv_scaled(double x, gsl_mode_t mode)
int gsl_sf_airy_Ai_deriv_scaled_e(double x, gsl_mode_t mode, gsl_sf_result *result)

조정된 에어리 도함수 \(S_A(x) Ai' (x)\) 의 값은 \(exp{(\frac{2}{3} x^\frac{3}{2})}\) 이고, \(x<0\) 일 때는 \(1\) 입니다.

double gsl_sf_airy_Bi_deriv_scaled(double x, gsl_mode_t mode)
int gsl_sf_airy_Bi_deriv_scaled_e(double x, gsl_mode_t mode, gsl_sf_result *result)

조정된 에어리 도함수 \(S_A(x) Bi' (x)\) 의 값을 계산합니다. \(x>0\) 의 값은 \(exp{(\frac{2}{3} x^\frac{3}{2})}\) 이고, \(x<0\) 일 때는 \(1\) 입니다.

에어리 함수의 근

double gsl_sf_airy_zero_Ai(unsigned int s)
int gsl_sf_airy_zero_Ai_e(unsigned int s, gsl_sf_result *result)

에어리 함수 \(Ai(x)\) 번째 근을 찾아 반환합니다.

double gsl_sf_airy_zero_Bi(unsigned int s)
int gsl_sf_airy_zero_Bi_e(unsigned int s, gsl_sf_result *result)

에어리 함수 \(Bi(x)\) 번째 근을 찾아 반환합니다.

에어리 도함수의 근

double gsl_sf_airy_zero_Ai_deriv(unsigned int s)
int gsl_sf_airy_zero_Ai_deriv_e(unsigned int s, gsl_sf_result *result)

에어리 도함수 \(Ai'(x)\) 번째 근을 찾아 반환합니다.

double gsl_sf_airy_zero_Bi_deriv(unsigned int s)
int gsl_sf_airy_zero_Bi_deriv_e(unsigned int s, gsl_sf_result *result)

에어리 도함수 \(Bi'(x)\) 번째 근을 찾아 반환합니다.

베셀 함수(Bessel function)

이 단원은 원통형 베셀 함수 \(J_n(x), Y_n(x)\), 수정 베셀 함수 \(I_n(x), K_n(x)\), 구면 베셀 함수 \(j_l(x), y_l(x)\) 그리고 수정 구면 베셀 함수 \(i_l(x), k_l(x)\) 에 정의되어 있습니다.

1종 베셀 함수

동의어: Regular Bessel function

double gsl_sf_bessel_J0(double x)
int gsl_sf_bessel_J0_e(double x, gsl_sf_result *result)

\(0\) 차수의 1종 베셀 함수 \(J_0(x)\) 의 값을 계산합니다.

double gsl_sf_bessel_J1(double x)
int gsl_sf_bessel_J1_e(double x, gsl_sf_result *result)

\(1\) 차수의 1종 베셀 함수 \(J_1(x)\) 의 값을 계산합니다.

double gsl_sf_bessel_Jn(int n, double x)
int gsl_sf_bessel_Jn_e(int n, double x, gsl_sf_result *result)

주어진 \(n\) 차수의 1종 베셀 함수 \(J_n(x)\) 의 값을 계산합니다.

int gsl_sf_bessel_Jn_array(int nmin, int nmax, double x, double result_array[])

nmin 에서 nmax 까지, 1종 베셀 함수 \(J_n(x)\) 의 값을 계산합니다. 계산 결과값은 result_array 배열에 저장됩니다. 재귀식을 이용해 계산 효율을 높여, 실제 값과는 조금 다를 수 있습니다.

2종 베셀 함수

동의어: 노이먼 함수(Neumann funtion), Irregular Bessel function

double gsl_sf_bessel_Y0(double x)
int gsl_sf_bessel_Y0_e(double x, gsl_sf_result *result)

\(x>0\) 에 대해, \(0\) 차수의 2종 베셀 함수 \(Y_0(x)\) 값을 계산합니다.

double gsl_sf_bessel_Y1(double x)
int gsl_sf_bessel_Y1_e(double x, gsl_sf_result *result)

\(x>0\) 에 대해, \(1\) 차수의 2종 베셀 함수 \(Y_1(x)\) 값을 계산합니다.

double gsl_sf_bessel_Yn(int n, double x)
int gsl_sf_bessel_Yn_e(int n, double x, gsl_sf_result *result)

\(x>0\) 에 대해, 주어진 \(n\) 차수의 2종 베셀 함수 \(Y_1(x)\) 값을 계산합니다.

int gsl_sf_bessel_Yn_array(int nmin, int nmax, double x, double result_array[])

nmin 에서 nmax 까지, 2종 베셀 함수 \(Y_n(x)\) 의 값을 계산합니다. 계산 결과값은 result_array 배열에 저장됩니다. 이 함수의 정의역은 \(x>0\) 입니다. 재귀식을 이용해 계산 효율을 높여, 실제 값과는 조금 다를 수 있습니다.

1종 변형 베셀 함수(Regular Modified Cylindrical Bessel Functions)

double gsl_sf_bessel_I0(double x)
int gsl_sf_bessel_I0_e(double x, gsl_sf_result *result)

\(0\) 차수의 1종 변형 베셀 함수 \(I_0(x)\) 의 값을 계산합니다.

double gsl_sf_bessel_I1(double x)
int gsl_sf_bessel_I1_e(double x, gsl_sf_result *result)

\(1\) 차수의 1종 변형 베셀 함수 \(1_0(x)\) 의 값을 계산합니다.

double gsl_sf_bessel_In(int n, double x)
int gsl_sf_bessel_In_e(int n, double x, gsl_sf_result *result)

주어진 \(n\) 차수의 1종 변형 베셀 함수 \(I_n(x)\) 의 값을 계산합니다.

int gsl_sf_bessel_In_array(int nmin, int nmax, double x, double result_array[])

nmin 에서 nmax 까지, 1종 변형 베셀 함수 \(I_n(x)\) 의 값을 계산합니다. 계산 결과값은 result_array 배열에 저장됩니다. nmin 는 반드시 양수이거나 \(0\) 이어야 합니다. 재귀식을 이용해 계산 효율을 높여, 실제 값과는 조금 다를 수 있습니다.

double gsl_sf_bessel_I0_scaled(double x)
int gsl_sf_bessel_I0_scaled_e(double x, gsl_sf_result *result)

조정 계수가 곱해진 \(0\) 차수의 1종 변형 베셀 함수 \(\text{exp}(-|x|) I_0(x)\) 를 계산합니다.

double gsl_sf_bessel_I1_scaled(double x)
int gsl_sf_bessel_I1_scaled_e(double x, gsl_sf_result *result)

조정 계수가 곱해진 \(1\) 차수의 1종 변형 베셀 함수 \(\text{exp}(-|x|) I_1(x)\) 를 계산합니다.

double gsl_sf_bessel_In_scaled(int n, double x)
int gsl_sf_bessel_In_scaled_e(int n, double x, gsl_sf_result *result)

조정 계수가 곱해진, \(n\) 차수의 1종 변형 베셀 함수 \(\text{exp}(-|x|) I_n(x)\) 를 계산합니다.

int gsl_sf_bessel_In_scaled_array(int nmin, int nmax, double x, double result_array[])

nmin 에서 nmax 까지, 조정 계수가 곱해진, 1종 변형 베셀 함수 \(\text{exp}(-|x|)I_n(x)\) 의 값을 계산합니다. 계산 결과값은 result_array 배열에 저장됩니다. nmin 반드시 양수이거나 \(0\) 이어야 합니다. 재귀식을 이용해 계산 효율을 높여, 실제 값과는 조금 다를 수 있습니다.

2종 변형 베셀 함수

double gsl_sf_bessel_K0(double x)
int gsl_sf_bessel_K0_e(double x, gsl_sf_result *result)

\(x>0\) 에 대해, \(0\) 차수의 2종 변형 베셀 함수 \(K_0(x)\) 값을 계산합니다.

double gsl_sf_bessel_K1(double x)
int gsl_sf_bessel_K1_e(double x, gsl_sf_result *result)

\(x>0\) 차수의 2종 변형 베셀 함수 \(K_1(x)\) 값을 계산합니다.

double gsl_sf_bessel_Kn(int n, double x)
int gsl_sf_bessel_Kn_e(int n, double x, gsl_sf_result *result)

\(x>0\) 차수의 2종 변형 베셀 함수 \(K_n(x)\) 값을 계산합니다.

int gsl_sf_bessel_Kn_array(int nmin, int nmax, double x, double result_array[])

nmin 에서 nmax 까지, 2종 변형 베셀 함수 \(K_n(x)\) 의 값을 계산합니다. 계산 결과값은 result_array 배열에 저장됩니다. nmin 반드시 양수이거나 \(0\) 이어야 합니다. 함수의 정의역은 \(x>0\) 입니다. 재귀식을 이용해 계산 효율을 높여, 실제 값과는 조금 다를 수 있습니다.

double gsl_sf_bessel_K0_scaled(double x)
int gsl_sf_bessel_K0_scaled_e(double x, gsl_sf_result *result)

\(x>0\) 에 대해, 조정 계수가 곱해진 \(0\) 차수의 2종 변형 베셀 함수 \(\text{exp}(x) K_0(x)\) 를 계산합니다.

double gsl_sf_bessel_K1_scaled(double x)
int gsl_sf_bessel_K1_scaled_e(double x, gsl_sf_result *result)

\(x>0\) 에 대해, 조정 계수가 곱해진 \(1\) 차수의 2종 변형 베셀 함수 \(\text{exp}(x) K_1(x)\) 를 계산합니다.

double gsl_sf_bessel_Kn_scaled(int n, double x)
int gsl_sf_bessel_Kn_scaled_e(int n, double x, gsl_sf_result *result)

\(x>0\) 차수의 2종 변형 베셀 함수 \(\text{exp}(x) K_n(x)\) 를 계산합니다.

int gsl_sf_bessel_Kn_scaled_array(int nmin, int nmax, double x, double result_array[])

nmin 에서 nmax 까지, 조정 계수가 곱해진 2종 변형 베셀 함수 \(\text{exp}(x) K_n(x)\) 의 값을 계산합니다. 계산 결과값은 result_array 배열에 저장됩니다. nmin 는 반드시 양수이거나 \(0\) 이어야 합니다. 함수의 정의역은 \(x>0\) 입니다. 재귀식을 이용해 계산 효율을 높여, 실제 값과는 조금 다를 수 있습니다.

1종 구면 베셀 함수

double gsl_sf_bessel_j0(double x)
int gsl_sf_bessel_j0_e(double x, gsl_sf_result *result)

\(0\) 차수의 1종 구면 베셀 함수 \(j_0 (x) = \sin (x) /x\) 의 값을 계산합니다.

double gsl_sf_bessel_j1(double x)
int gsl_sf_bessel_j1_e(double x, gsl_sf_result *result)

\(1\) 차수의 1종 구면 베셀 함수 \(j_1 (x) = (\sin (x) /x - \cos(x)) /x\) 의 값을 계산합니다.

double gsl_sf_bessel_j2(double x)
int gsl_sf_bessel_j2_e(double x, gsl_sf_result *result)

\(2\) 차수의 1종 구면 베셀 함수 \(j_2 (x) = ((3/x^2 -1)\sin(x) -3 \cos(x)/x) /x\) 의 값을 계산합니다.

double gsl_sf_bessel_jl(int l, double x)
int gsl_sf_bessel_jl_e(int l, double x, gsl_sf_result *result)

\(l\) 차수의 1종 구면 베셀 함수 \(j_l (x)\) 의 값을 계산합니다. \(x,l\)\(l \geq 0, x \geq 0\) 이어야 합니다.

int gsl_sf_bessel_jl_array(int lmax, double x, double result_array[])

\(lmax \geq 0, x \geq 0\) 에 대해, 1종 구면 베셀 함수 \(j_l(x)\) 의 값을 \(l=0\) 에서 \(l=lmax\) 까지 계산합니다. 계산 결과값은 result_array 배열에 저장됩니다. 재귀식을 이용해 계산 효율을 높여, 실제 값과는 조금 다를 수 있습니다.

int gsl_sf_bessel_jl_steed_array(int lmax, double x, double *result_array)

Steed 방법을 이용해 1종 구면 베셀 함수 \(j_l(x)\) 의 값을 \(l=0\) 에서 \(l=lmax\) 까지 계산합니다. \(lmax, x\)\(lmax \geq 0, x \geq 0\) 이어야 합니다. 계산 결과값은 result_array 배열에 저장됩니다. Steed/Barnett 알고리즘은 Comp. Phys. Comm. 21, 297(1981)에 기술되어 있습니다. Steed 방법은 다른 함수의 재귀적 방법보다 더 안정적이지만, 그 대신 더 느립니다.

2종 구면 베셀 함수

double gsl_sf_bessel_y0(double x)
int gsl_sf_bessel_y0_e(double x, gsl_sf_result *result)

\(0\) 차수의 2종 구면 베셀 함수 \(y_0 (x) = -\cos (x) /x\) 의 값을 계산합니다.

double gsl_sf_bessel_y1(double x)
int gsl_sf_bessel_y1_e(double x, gsl_sf_result *result)

\(1\) 차수의 2종 구면 베셀 함수 \(y_1 (x) = -(\cos (x) /x + \sin (x))/x\) 의 값을 계산합니다.

double gsl_sf_bessel_y2(double x)
int gsl_sf_bessel_y2_e(double x, gsl_sf_result *result)

\(2\) 차수의 2종 구면 베셀 함수 \(y_2 (x) = (-3/x^3 + 1/x)\cos(x) - (3/x^2)\sin(x)\) 의 값을 계산합니다.

double gsl_sf_bessel_yl(int l, double x)
int gsl_sf_bessel_yl_e(int l, double x, gsl_sf_result *result)

\(l \geq 0\) 에 대해, \(l\) 차수의 2종 구면 베셀 함수 \(y_l (x)\) 의 값을 계산합니다.

int gsl_sf_bessel_yl_array(int lmax, double x, double result_array[])

\(lmax \geq 0\) 에 대해, 2종 구면 베셀 함수 \(y_l(x)\) 의 값을 \(l=0\) 에서 \(l=lmax\) 까지 계산합니다. 계산 결과값은 result_array 배열에 저장됩니다. 재귀식을 이용해 계산 효율을 높여, 실제 값과는 조금 다를 수 있습니다.

1종 변형 구면 베셀 함수

1종 변형 구면 베셀 함수 \(i_l(x)\) 는 분수 차수의 1종 수정 베셀 함수와 다음과 같은 관계를 가집니다.

\[i_l(x) = \sqrt{\pi/(2x)}I_{l + 1/2}(x)\]
double gsl_sf_bessel_i0_scaled(double x)
int gsl_sf_bessel_i0_scaled_e(double x, gsl_sf_result *result)

조정 계수가 곱해진, \(0\) 차수의 1종 변형 구면 베셀 함수 \(\text{exp}(-|x|) i_0 (x)\) 를 계산합니다.

double gsl_sf_bessel_i1_scaled(double x)
int gsl_sf_bessel_i1_scaled_e(double x, gsl_sf_result *result)

조정 계수가 곱해진, \(1\) 차수의 1종 변형 구면 베셀 함수 \(\text{exp}(-|x|) i_1 (x)\) 를 계산합니다.

double gsl_sf_bessel_i2_scaled(double x)
int gsl_sf_bessel_i2_scaled_e(double x, gsl_sf_result *result)

조정 계수가 곱해진, \(2\) 차수의 1종 변형 구면 베셀 함수 \(\text{exp}(-|x|) i_2 (x)\) 를 계산합니다.

double gsl_sf_bessel_il_scaled(int l, double x)
int gsl_sf_bessel_il_scaled_e(int l, double x, gsl_sf_result *result)

조정 계수가 곱해진, \(l\) 차수의 1종 변형 구면 베셀 함수 \(\text{exp}(-|x|) i_2 (x)\) 를 계산합니다.

int gsl_sf_bessel_il_scaled_array(int lmax, double x, double result_array[])

\(lmax \geq 0, x \geq 0\) 에 대해, 조정 계수가 곱해진 1종 변형 구면 베셀 함수 \(\text{exp}(-|x|) i_l(x)\) 의 값을 \(l=0\) 에서 \(l=lmax\) 까지 계산합니다. 계산 결과값은 result_array 배열에 저장됩니다. 재귀식을 이용해 계산 효율을 높여, 실제 값과는 조금 다를 수 있습니다.

2종 변형 구면 베셀 함수

2종 변형 구면 베셀 함수 \(k_l(x)\) 는 분수 차수 2종 구면 베셀 함수와 다음과 같은 관계를 가집니다.

\[k_l(x) = \sqrt{\pi / (2x) K_{l+1/2}(x)}\]
double gsl_sf_bessel_k0_scaled(double x)
int gsl_sf_bessel_k0_scaled_e(double x, gsl_sf_result *result)

\(x>0\) 에 대해, 조정 계수가 곱해진 \(0\) 차수의 2종 변형 구면 베셀 함수 \(\text{exp}(x)k_0(x)\) 의 값을 계산합니다.

double gsl_sf_bessel_k1_scaled(double x)
int gsl_sf_bessel_k1_scaled_e(double x, gsl_sf_result *result)

\(x>0\) 에 대해, 조정 계수가 곱해진 \(1\) 차수의 2종 변형 구면 베셀 함수 \(\text{exp}(x)k_1(x)\) 의 값을 계산합니다.

double gsl_sf_bessel_k2_scaled(double x)
int gsl_sf_bessel_k2_scaled_e(double x, gsl_sf_result *result)

\(x>0\) 에 대해, 조정 계수가 곱해진 \(2\) 차수의 2종 변형 구면 베셀 함수 \(\text{exp}(x)k_2(x)\) 의 값을 계산합니다.

double gsl_sf_bessel_kl_scaled(int l, double x)
int gsl_sf_bessel_kl_scaled_e(int l, double x, gsl_sf_result *result)

\(x>0\) 에 대해, 조정 계수가 곱해진 \(l\) 차수의 2종 변형 구면 베셀 함수 \(\text{exp}(x)k_l(x)\) 의 값을 계산합니다.

int gsl_sf_bessel_kl_scaled_array(int lmax, double x, double result_array[])

\(lmax \geq 0, x \geq 0\) 에 대해, 조정 계수가 곱해진 1종 변형 구면 베셀 함수 \(\text{exp}(x) k_l(x)\) 의 값을 \(l=0\) 에서 \(l=lmax\) 까지 계산합니다. 계산 결과값은 result_array 배열에 저장됩니다. 재귀식을 이용해 계산 효율을 높여, 실제 값과는 조금 다를 수 있습니다.

1종 베셀 함수-분수 차수

double gsl_sf_bessel_Jnu(double nu, double x)
int gsl_sf_bessel_Jnu_e(double nu, double x, gsl_sf_result *result)

분수 차수 \(\nu\) 에 대해, 1종 베셀 함수 \(J_\nu (x)\) 의 값을 계산합니다.

int gsl_sf_bessel_sequence_Jnu_e(double nu, gsl_mode_t mode, size_t size, double v[])

분수 차수 \(\nu\) 의 1종 배셀함수 \(J_\nu (x)\) 의 값을 주어진 \(x\) 값 배열에 대해 계산합니다. \(size\) 길이의 배열 \(v\)\(x\) 값들을 담고있습니다. 함수는 이 배열이 양수가 순차적으로 배열되어 있다 가정합니다. \(v\) 배열을 수정해 \(J_\nu (x_i)\) 의 값을 덮어 씌웁니다.

2종 베셀 함수-분수 차수

double gsl_sf_bessel_Ynu(double nu, double x)
int gsl_sf_bessel_Ynu_e(double nu, double x, gsl_sf_result *result)

분수 차수 \(\nu\) 에 대해, 2종 베셀 함수 \(Y_\nu (x)\) 의 값을 계산합니다.

1종 변형 베셀 함수-분수 차수

double gsl_sf_bessel_Inu(double nu, double x)
int gsl_sf_bessel_Inu_e(double nu, double x, gsl_sf_result *result)

\(x>0, \nu>0\) 에 대해, 분수 차수 \(\nu\) 의 1종 변형 베셀 함수 \(I_\nu(x)\) 를 계산합니다.

double gsl_sf_bessel_Inu_scaled(double nu, double x)
int gsl_sf_bessel_Inu_scaled_e(double nu, double x, gsl_sf_result *result)

\(x>0, \nu>0\) 에 대해, 조정 계수가 곱해진 분수 차수 \(\nu\) 의 2종 변형 베셀 함수 \(\text{exp}(-|x|)I_\nu (x)\) 의 값을 계산합니다.

2종 변형 베셀 함수-분수 차수

double gsl_sf_bessel_Knu(double nu, double x)
int gsl_sf_bessel_Knu_e(double nu, double x, gsl_sf_result *result)

\(x>0, \nu>0\) 에 대해, 분수 차수 \(\nu\) 의 2종 변형 베셀 함수 \(K_\nu (x)\) 의 값을 계산합니다.

double gsl_sf_bessel_lnKnu(double nu, double x)
int gsl_sf_bessel_lnKnu_e(double nu, double x, gsl_sf_result *result)

\(x>0, \nu>0\) 에 대해, 로그가 씌워진, 분수 차수 \(\nu\) 의 2종 변형 베셀 함수 \(\ln(K_\nu (x))\) 의 값을 계산합니다.

double gsl_sf_bessel_Knu_scaled(double nu, double x)
int gsl_sf_bessel_Knu_scaled_e(double nu, double x, gsl_sf_result *result)

\(x>0, \nu>0\) 에 대해, 조정 계수가 곱해진 분수 차수 \(\nu\) 의 2종 변형 베셀 함수 \(\text{exp}(+|x|) K_\nu (x)\) 의 값을 계산합니다.

1종 베셀 함수의 근

double gsl_sf_bessel_zero_J0(unsigned int s)
int gsl_sf_bessel_zero_J0_e(unsigned int s, gsl_sf_result *result)

\(0\) 차수의 베셀 함수 \(J_0(x)\) 번째, 양수 근의 위치를 찾습니다.

double gsl_sf_bessel_zero_J1(unsigned int s)
int gsl_sf_bessel_zero_J1_e(unsigned int s, gsl_sf_result *result)

\(1\) 차수의 베셀 함수 \(J_1(x)\) 번째, 양수 근의 위치를 찾습니다.

double gsl_sf_bessel_zero_Jnu(double nu, unsigned int s)
int gsl_sf_bessel_zero_Jnu_e(double nu, unsigned int s, gsl_sf_result *result)

베셀 함수 \(J_\nu (x)\) 번째, 양수 근의 위치를 찾습니다. 현재 구현체는 nu 가 음수일 때를 지원하지 않습니다.

클라우센 함수 (Clausen Functions)

클라우센 함수는 다음과 같은 적분으로 정의됩니다.

\[Cl_2(x) = - \int_0^x \log(2 \sin (\frac{t}{2})) \, dt\]

이 함수는 다이로그 함수와 다음과 같은 관계를 가집니다.

\[Cl_2 (\theta) = \Im[Li_2 (e^{i \theta})]\]

클라우센 함수들은 헤더 파일 gsl_sf_clausen.h 에 정의되어 있습니다.

double gsl_sf_clausen(double x)
int gsl_sf_clausen_e(double x, gsl_sf_result *result)

이 함수들은 클라우센 적분 \(Cl_2(x)\) 의 값을 계산합니다.

쿨롱 함수 (Coulomb Functions)

쿨롱 함수의 초안은 헤더 파일 gsl_sf_coulomb.h 에 정의되어 있습니다. 구속 상태와 산란해를 모두 포함합니다.

정규화 된 수소의 구속 상태

double gsl_sf_hydrogenicR_1(double Z, double r)
int gsl_sf_hydrogenicR_1_e(double Z, double r, gsl_sf_result *result)

가장 낮은 차수의 정규화 된 수소의 구속 방사형 상태 함수 \(R_1 = 2 Z \sqrt{Z} \exp (-Zr)\) 를 계산합니다.

double gsl_sf_hydrogenicR(int n, int l, double Z, double r)
int gsl_sf_hydrogenicR_e(int n, int l, double Z, double r, gsl_sf_result *result)

\(n\) 차수의 정규화 된 수소의 구속 방사형 상태 함수

\[R_n := \frac{2 Z^{3/2}}{n^2} (\frac{2Zr}{n})^l \sqrt{\frac{(n-l-1)!}{n+l}} \exp(-Zr /n) L_{n-l-1}^{2l+1} (2Zr/n)\]

를 계산합니다. \(L_b^a(x)\) 는 일반화 된 르장드르 다항식입니다. 이 정규화 방법은 수소의 파동 함수 \(\psi\)\(\psi(n,l,r) = R_n Y_{lm}\) 로 주어지도록 정해졌습니다.

쿨롱 파동 함수

쿨롱 파동함수 \(F_L(\eta,x), G_L (\eta,x)\) 들은 Abramowitz & Stegun Chapter 14에 기술되어 있습니다. 다양하고 넓은 범위의 값들을 가지기 때문에, 오버플로우를 이용할 수 있도록 구현되었습니다. 만약, 오버플로우가 발생한다면, GSL_EOVERFLW 값을 반환하고 exponent(s) 는 수정 가능한 인자 exp_Fexp_G 를 반환합니다. 완전해는 다음과 같이 구할 수 있습니다.

\[\begin{split}F_L(\eta,x) = f_c[k_L] * \exp(exp_F) \\ G_L(\eta,x) = g_c[k_L] * \exp(exp_G)\end{split}\]
\[\begin{split}{F'}_L(\eta,x)= f_{cp}[k_L] * \exp(exp_F)\\ {G'}_L(\eta,x)= g_{cp}[k_L] * \exp(exp_G)\end{split}\]
int gsl_sf_coulomb_wave_FG_e(double eta, double x, double L_F, int k, gsl_sf_result *F, gsl_sf_result *Fp, gsl_sf_result *G, gsl_sf_result *Gp, double *exp_F, double *exp_G)

이 함수는 쿨롱 파동 함수 \(F_L(\eta,x), G_{L-k}(\eta,x)\) 와 그 도함수 \(F'_L(\eta,x), G'_{L-k}(\eta,x)\) 를 인자 x`에 대해 계산합니다. 계수들은 :math:`L 에 대해 다음의 제약조건을 가집니다.

\[L-k > -1/2, x>0, k \in \mathbb{Z}\]

유의할 점은 \(L\) 이 반드시 정수로 있을 필요는 없다는 점입니다. 결과 값들은 인자 FG 에 함수 값이 저장되고, 도함수의 값은 FpGp 에 저장됩니다. 오버플로우가 발생하면, GSL_EOVERFLW 가 반환되고 조정된 지수값이 수정 가능한 인자 exp_Fexp_G 에 저장됩니다.

int gsl_sf_coulomb_wave_F_array(double L_min, int kmax, double eta, double x, double fc_array[], double *F_exponent)

이 함수는 \(L = L_{min} \dots L_{min} + k_{max}\) 에 대해, 함수 \(F_L(\eta,x)\) 의 값을 계산합니다. 계산 결과값은 fc_array 배열에 저장됩니다. 오버플로우가 발생하면 지수값이 F_exponenet 에 저장됩니다.

int gsl_sf_coulomb_wave_FG_array(double L_min, int kmax, double eta, double x, double fc_array[], double gc_array[], double *F_exponent, double *G_exponent)

이 함수는 \(L = L_{min} \dots L_{min} + k_{max}\) 에 대해, 함수 \(F_L(\eta,x), G_L(\eta,x)\) 의 값을 계산합니다. 계산 결과값은 각각 fc_arraygc_array 배열에 저장됩니다. 오버플로우가 발생하면 지수값이 F_exponenetG_exponent 에 저장됩니다.

int gsl_sf_coulomb_wave_FGp_array(double L_min, int kmax, double eta, double x, double fc_array[], double fcp_array[], double gc_array[], double gcp_array[], double *F_exponent, double *G_exponent)

이 함수는 \(L = L_{min} \dots L_{min} + k_{max}\) 에 대해, 함수 \(F_L(\eta,x), G_L(\eta,x)\) 와 그 도함수 \(F'_L(\eta,x), G'_L(\eta,x)\) 의 값을 계산합니다. 계산 결과값은 각각 fc_array , gc_array , fcp_array 그리고 gcp_array 배열에 저장됩니다. 오버플로우가 발생하면 지수값이 F_exponenetG_exponent 에 저장됩니다.

int gsl_sf_coulomb_wave_sphF_array(double L_min, int kmax, double eta, double x, double fc_array[], double F_exponent[])

이 함수는 \(L = L_{min} \dots L_{min} + k_{max}\) 에 대해, 인자로 나누어진 쿨롱 함수 \(F_L(\eta,x)/x\) 값을 계산합니다. 계산 결과값은 fc_array 배열에 저장됩니다. 오버플로우가 발생하면 지수값이 F_exponenet 에 저장됩니다. \(\eta \rightarrow 0\) 이 함수는 구면 베셀 함수로 수렴합니다.

쿨롱 파동함수의 정규화 계수

쿨롱 파동 함수의 정규화 상수들은 Abramowitz 14.1.7에 정의되어 있습니다.

int gsl_sf_coulomb_CL_e(double L, double eta, gsl_sf_result *result)

\(L>-1\) 에 대해, 쿨롱 파동 함수의 정규화 계수 \(C_L (\eta)\) 를 계산합니다.

int gsl_sf_coulomb_CL_array(double Lmin, int kmax, double eta, double cl[])

\(L = L_{min} \dots L_{min} + k_{max}, L_{min} > -1\) 에 대해, 쿨롱 파동 함수의 정규화 계수 \(C_L(\eta)\) 를 계산합니다.

상호작용 계수 (Coupling Coefficients)

위그너 \(3-j\) , \(6-j\) , 그리고 \(9-j\) 기호들은, 결합된 각 운동량 벡터들의 상호 작용 계수들을 나타냅니다. 표준 상호작용 계수 함수들의 인자가 정수거나 반-정수(half-integer)이기 때문에, 이를 계산하는 함수들의 인자들도 이와 같습니다. 이 정수들은 실제 스핀 값들의 두 배를 의미합니다. \(3-j\) 계수들에 대한 더 자세한 정보는 Abramowitz & Stegun, Section 27.9를 참고할 수 있습니다. 이 단원의 함수들은 헤더 파일 gsl_sf_coupling.h 에 기술되어 있습니다.

3-j 기호

double gsl_sf_coupling_3j(int two_ja, int two_jb, int two_jc, int two_ma, int two_mb, int two_mc)
int gsl_sf_coupling_3j_e(int two_ja, int two_jb, int two_jc, int two_ma, int two_mb, int two_mc, gsl_sf_result *result)

위그너 \(3-j\) 계수

\[\begin{split}\begin{pmatrix} ja & jb & jc \\\\ ma & mb & mc \end{pmatrix}\end{split}\]

를 계산합니다. 인자들은 반-정수 값을 가집니다. 예를 들어, \(ja=\) two_ja /2 이고, \(ma=\) two_ma /2 입니다.

6-j 기호(6-j Symbols)

double gsl_sf_coupling_6j(int two_ja, int two_jb, int two_jc, int two_jd, int two_je, int two_jf)
int gsl_sf_coupling_6j_e(int two_ja, int two_jb, int two_jc, int two_jd, int two_je, int two_jf, gsl_sf_result *result)

위그너 \(6-j\) 계수

\[\begin{split}\begin{Bmatrix} ja & jb & jc \\\\ jd & je & jf \end{Bmatrix}\end{split}\]

를 계산합니다. 인자들은 반-정수 값을 가집니다. 예를 들어, \(ja=\) two_ja /2 이고 \(ma=\) two_ma /2 입니다.

9-j 기호(9-j Symbols)

double gsl_sf_coupling_9j(int two_ja, int two_jb, int two_jc, int two_jd, int two_je, int two_jf, int two_jg, int two_jh, int two_ji)
int gsl_sf_coupling_9j_e(int two_ja, int two_jb, int two_jc, int two_jd, int two_je, int two_jf, int two_jg, int two_jh, int two_ji, gsl_sf_result *result)

위그너 \(9-j\) 계수

\[\begin{split}\begin{Bmatrix} ja & jb & jc \\\\ jd & je & jf \\\\ jg & jh & ji \end{Bmatrix}\end{split}\]

를 계산합니다. 인자들은 반-정수 값을 값을 가집니다. 예를 들어, \(ja=\) two_ja /2 이고 \(ma=\) two_ma /2 입니다.

도슨 함수 (Dawson Function)

도슨 적분은 다음의 적분을 의미합니다.

\[D_+(x) = e^{-x^2} \int_0^x e^{t^2} \, dt\]

도슨 적분표는 Abramowitz & Stegun의 표7.5에서 찾을 수 있습니다. 이 도슨 함수는 헤더 파일 gsl_sf_dawson.h 에 정의되어 있습니다.

double gsl_sf_dawson(double x)
int gsl_sf_dawson_e(double x, gsl_sf_result *result)

이 함수들은 주어진 값 \(x\) 대해, 도슨 적분 값을 계산합니다.

디바이 함수 (Debye Functions)

디바이 함수 \(D_n(x)\) 는 다음과 같이 정의됩니다.

\[D_n(x) = \frac{n}{x^n} \int_0^x \frac{t^n}{e^t -1} \, dt\]

더 자세한 정보는 Abramowitz & Stegun, Section 27.1을 참고할 수 있습니다. 디바이 함수는 헤더 파일 gsl_sf_debye.h 에 정의되어 있습니다.

double gsl_sf_debye_1(double x)
int gsl_sf_debye_1_e(double x, gsl_sf_result *result)

\(1\) 차 디바이 함수 \(D_1(x)\) 의 값을 계산합니다.

double gsl_sf_debye_2(double x)
int gsl_sf_debye_2_e(double x, gsl_sf_result *result)

\(2\) 차 디바이 함수 \(D_2(x)\) 의 값을 계산합니다.

double gsl_sf_debye_3(double x)
int gsl_sf_debye_3_e(double x, gsl_sf_result *result)

\(3\) 차 디바이 함수 \(D_3(x)\) 의 값을 계산합니다.

double gsl_sf_debye_4(double x)
int gsl_sf_debye_4_e(double x, gsl_sf_result *result)

\(4\) 차 디바이 함수 \(D_4(x)\) 의 값을 계산합니다.

double gsl_sf_debye_5(double x)
int gsl_sf_debye_5_e(double x, gsl_sf_result *result)

\(5\) 차 디바이 함수 \(D_5(x)\) 의 값을 계산합니다.

double gsl_sf_debye_6(double x)
int gsl_sf_debye_6_e(double x, gsl_sf_result *result)

\(6\) 차 디바이 함수 \(D_6(x)\) 의 값을 계산합니다.

다이로그 함수 (Dilogarithm)

다이 로그는 다음과 같이 정의됩니다.

\[Li_2 (z) = - \int_0^z \frac{\log(1-s)}{s} \, ds\]

이 함수들은 헤더 파일 gsl_sf_dilog.h 에 정의되어 있습니다.

실수 인자 (Real Argument)

double gsl_sf_dilog(double x)
int gsl_sf_dilog_e(double x, gsl_sf_result *result)

이 함수들은 실수 값에 대한 다이 로그값을 계산합니다. 르윈의 표기법으로 \(Li_2(x)\) 로 표기되며, 이는 실수 값 \(x\) 의 다이 로그 값을 나타냅니다. 이 값은 다음과 같이 적분형 표현으로 정의됩니다.

\[Li_2(x) = - \mathfrak{R} \int_0^x \log(1-s) /s \, ds\]

참고

  • \(\mathfrak{I}(Li_2(x))\)\(x \leq 1\) 에서 \(0\) , \(x>1\) 에서 \(- \pi \log(x)\) 입니다.

  • Abramowitz & Stegun에서는 스펜스 적분을 \(Li_2(x)\) 대신 \(S(x) = Li_2(1-x)\) 로 정의합니다.

복소수 인자 (Complex Argument)

int gsl_sf_complex_dilog_e(double r, double theta, gsl_sf_result *result_re, gsl_sf_result *result_im)

복소수 인자 \(z = r \exp(i\theta)\) 에 대해, 복소수 값을 가지는 완전한 다이 로그 값을 계산합니다. 실수, 허수 값은 각각 result_reresult_im 에 반환됩니다.

기초 연산 (Elementary Operations)

다음 함수들은 곱셈 과정에서 오차의 전파를 같이 계산합니다. 이 함수들은 헤더 파일 gsl_sf_elementary.h 에 정의되어 있습니다.

double gsl_sf_multiply(double x, double y)
int gsl_sf_multiply_e(double x, double y, gsl_sf_result *result)

\(x\)\(y\) 에 저장합니다.

int gsl_sf_multiply_err_e(double x, double dx, double y, double dy, gsl_sf_result *result)

\(x\)\(y\) 의 곱셈을 절대 오차 \(dx\)\(dy\) 와 함께 연산합니다. result 에는 \(xy \pm xy\sqrt{(dx/x)^2 + (dy/y)^2}\) 가 저장됩니다.

타원 적분 (Elliptic Integrals)

이 단원에서 기술된 함수들은 헤더 파일 gsl_sf_ellint.h 에 정의되어 있습니다. 타원 적분에 관한 더 자세한 정보들은 Abramowitz & Stegun, Chpater 17를 참고해 볼 수 있습니다.

르장드르 형태 정의 (Definition of Legendre Forms)

타원 적분의 르장드르 형태 \(F(\phi, k), E(\phi, k)\) 그리고 \(\Pi(\phi, k, n)\) 는 다음과 같이 정의됩니다.

\[\begin{split}F(\phi,k) = \int_0^{\phi} dt \frac{1}{\sqrt(1- k^2 \sin^2 (t))}\\ E(\phi,k) = \int_0^{\phi} dt {\sqrt(1- k^2 \sin^2 (t))}\\ \Pi(\phi, k, n) = \int_0^{\phi} dt \frac{1}{(1+n \sin^2(t)) \sqrt{1-k^2 \sin^2 (t)}}\end{split}\]

르장드르 형태의 완전 타원 적분은 \(K(k) = F(\pi/2, k)\)\(E(k) = E(\pi/2, k)\) 로 표기할 수 있습니다.

이 단원에서 사용하는 표기법은 Carlson, “Numerische Mathematik” 33 (1979) 1 에 기반해 있습니다. 이 표기는 Abramowitz & Stegun의 표기와 조금 다른데, 함수의 표기에서 인자가 \(m=k^2\) 이고, \(n\)\(-n\) 으로 바뀌어 있습니다.

칼슨 형태 정의 (Definition of Carlson Forms)

킬슨 대칭 형태 함수 \(RC(x,y), RD(x,y,z), RF(x,y,z)\)\(RJ(x,y,z,p)\) 는 다음과 같이 정의됩니다.

\[\begin{split}RC(x,y) &= 1/2 \int_0^\infty dt (t+x)^{-1/2} (t+y)^{-1} \\ RD(x,y,z) &= 3/2 \int_0^\infty dt (t+x)^{-1/2} (t+y)^{-1/2} (t+z)^{-3/2} \\ RF(x,y,z) &= 1/2 \int_0^\infty dt (t+x)^{-1/2} (t+y)^{-1/2} (t+z)^{-1/2} \\ RJ(x,y,z,p) &= 3/2 \int_0^\infty dt (t+x)^{-1/2} (t+y)^{-1/2} (t+z)^{-1/2} (t+p)^{-1}\end{split}\]

르장드르 형태-완전 타원 적분 (Legendre Form of Complete Elliptic Integrals)

double gsl_sf_ellint_Kcomp(double k, gsl_mode_t mode)
int gsl_sf_ellint_Kcomp_e(double k, gsl_mode_t mode, gsl_sf_result *result)

완전 타원 적분의 르장드르 형태 \(K(k)\)mode 변수의 값에 따라 정확도를 결정해 계산합니다. Abramowitz & Stegun에서 이 함수는 \(m=k^2\) 이라는 점에 유의해야 합니다.

double gsl_sf_ellint_Ecomp(double k, gsl_mode_t mode)
int gsl_sf_ellint_Ecomp_e(double k, gsl_mode_t mode, gsl_sf_result *result)

완전 타원 적분의 르장드르 형태 \(E(k)\)mode 변수의 값에 따라 정확도를 결정해 계산합니다. Abramowitz & Stegun에서 이 함수는 \(m=k^2\) 이라는 점에 유의해야 합니다.

double gsl_sf_ellint_Pcomp(double k, double n, gsl_mode_t mode)
int gsl_sf_ellint_Pcomp_e(double k, double n, gsl_mode_t mode, gsl_sf_result *result)

완전 타원 적분의 르장드르 형태 \(\Pi(k, n)\)mode 변수의 값에 따라 정확도를 결정해 계산합니다. Abramowitz & Stegun에서 이 함수는 \(m=k^2\) 이고 \(\sin^2(\alpha) = k^2\) 이고, \(n\) 의 부호를 \(-n\) 으로 바뀌었다는 점에 유의해야 합니다.

르장드르 형태-불완전 타원 적분 (Legendre Form of Incomplete Elliptic Integrals)

double gsl_sf_ellint_F(double phi, double k, gsl_mode_t mode)
int gsl_sf_ellint_F_e(double phi, double k, gsl_mode_t mode, gsl_sf_result *result)

타원 적분의 르장드르 형태 \(K(\phi, k)\)mode 변수의 값에 따라 정확도를 결정해 계산합니다. Abramowitz & Stegun에서 이 함수는 \(m=k^2\) 이라는 점에 유의해야 합니다.

double gsl_sf_ellint_E(double phi, double k, gsl_mode_t mode)
int gsl_sf_ellint_E_e(double phi, double k, gsl_mode_t mode, gsl_sf_result *result)

완전 타원 적분의 르장드르 형태 \(E(\phi, k)\)mode 변수의 값에 따라 정확도를 결정해 계산합니다. Abramowitz & Stegun에서 이 함수는 \(m=k^2\) 이라는 점에 유의해야 합니다.

double gsl_sf_ellint_P(double phi, double k, double n, gsl_mode_t mode)
int gsl_sf_ellint_P_e(double phi, double k, double n, gsl_mode_t mode, gsl_sf_result *result)

완전 타원 적분의 르장드르 형태 \(\Pi(\phi, k, n)\)mode 변수의 값에 따라 정확도를 결정해 계산합니다. Abramowitz & Stegun에서 이 함수는 \(m=k^2\) 이고 \(\sin^2(\alpha) = k^2\) 이고, \(n\) 의 부호를 \(-n\) 으로 바뀌었다는 점에 유의해야 합니다.

double gsl_sf_ellint_D(double phi, double k, gsl_mode_t mode)
int gsl_sf_ellint_D_e(double phi, double k, gsl_mode_t mode, gsl_sf_result *result)

타원 적분 \(D(\phi, k)\) 을 계산합니다. 이 함수는 칼슨 형태 타원 함수 \(RD(x,y,z)\) 와 다음의 관계로 정의되어 있습니다.

\[D(\phi, k) = \frac{1}{3}(\sin \phi)^3 RD(1-\sin^2(\phi), 1-k^2\sin^2(\phi),1)\]

칼슨 형태 (Carlson Forms)

double gsl_sf_ellint_RC(double x, double y, gsl_mode_t mode)
int gsl_sf_ellint_RC_e(double x, double y, gsl_mode_t mode, gsl_sf_result *result)

타원적분 \(RC(x, y)\)mode 변수의 값에 따라 정확도를 결정해 계산합니다.

double gsl_sf_ellint_RD(double x, double y, double z, gsl_mode_t mode)
int gsl_sf_ellint_RD_e(double x, double y, double z, gsl_mode_t mode, gsl_sf_result *result)

타원적분 \(RD(x ,y, z)\)mode 변수의 값에 따라 정확도를 결정해 계산합니다.

double gsl_sf_ellint_RF(double x, double y, double z, gsl_mode_t mode)
int gsl_sf_ellint_RF_e(double x, double y, double z, gsl_mode_t mode, gsl_sf_result *result)

타원적분 \(RF(x, y, z)\)mode 변수의 값에 따라 정확도를 결정해 계산합니다.

double gsl_sf_ellint_RJ(double x, double y, double z, double p, gsl_mode_t mode)
int gsl_sf_ellint_RJ_e(double x, double y, double z, double p, gsl_mode_t mode, gsl_sf_result *result)

타원적분 \(RJ(x, y, z, p)\)mode 변수의 값에 따라 정확도를 결정해 계산합니다.

자코비 타원 적분 (Elliptic Functions (Jacobi))

자코비 타원 함수들은 Abramowitz & Stegun, 16 단원의 정의를 따릅니다. 이 함수들은 gsl_sf_elljac.h 에 정의되어 있습니다.

int gsl_sf_elljac_e(double u, double m, double *sn, double *cn, double *dn)

자코비 타원 함수 \(sn(u\|m), cn(u\|m) , dn(u\|m)\) 을 내림차순 란덴 변환을 이용해 계산합니다.

오차 함수 (Error Functions)

오차 함수는 Abramowitz & Stegun, Chapter 7.에 기술되어 있습니다. 이 단원에서 기술된 함수들은 헤더파일 gs_sf_erf.h 에 기술되어 있습니다.

오차 함수(Error Functionn)

double gsl_sf_erf(double x)
int gsl_sf_erf_e(double x, gsl_sf_result *result)

오차함수 \(\text{erf}(x)\) 의 값을 계산합니다. \(\text{erf}(x)\) 는 다음과 같이 정의됩니다.

\[\text{erf}(x) = \frac{2}{\sqrt{\pi}} \int_0^x dt \exp(-t^2)\]

상보 오차 함수(Complementary Error Functionn)

double gsl_sf_erfc(double x)
int gsl_sf_erfc_e(double x, gsl_sf_result *result)

상보 오차 함수 \(\text{erfc}(x) = 1- \text{erf}(x)\) 의 값을 계산합니다. \(\text{erfc}(x)\) 는 다음과 같이 정의됩니다.

\[\text{erfc}(x) = \frac{2}{\sqrt{\pi}} \int_x^\infty dt \exp(-t^2)\]

로그 상보 오차 함수(Log Complementary Error Functionn)

double gsl_sf_log_erfc(double x)
int gsl_sf_log_erfc_e(double x, gsl_sf_result *result)

상보 오차 함수의 로그값 \(\log(\text{erfc}(x))\) 를 계산합니다.

확률 함수 (Probability functions)

표준/가우스 분포의 확률 함수들은 Abramowitz & Stegun, Section 26.2에 기술되어 있습니다.

double gsl_sf_erf_Z(double x)
int gsl_sf_erf_Z_e(double x, gsl_sf_result *result)

가우스 확률 밀도 함수 \(Z(x) = \frac{1}{\sqrt{2\pi}} \exp(- \frac{x^2}{2})\) 값을 계산합니다.

double gsl_sf_erf_Q(double x)
int gsl_sf_erf_Q_e(double x, gsl_sf_result *result)

가우스 확률 함수의 \(Q(x) = (1/\sqrt{2\pi}) \int_x^\infty dt \exp(-t^2/2)\) 값을 계산합니다.

표준 분포의 하자드 함수(Hazard function) 는 Mills’ 비의 역으로도 알려져 있습니다. 이는 다음과 같이 정의됩니다.

\[h(x) = \frac{Z(x)}{Q(x)} = \sqrt{\frac{2}{\pi}} \frac{\exp(- x^2/2)}{\text{erfc}(x/\sqrt{2})}\]

\(x\)\(-\infty\) 에 가까우질 수록 급격히 감소하며, \(x\)\(+\infty\) 에 가까워질 수록 \(h(x) \approx\) 로 점근합니다.

double gsl_sf_hazard(double x)
int gsl_sf_hazard_e(double x, gsl_sf_result *result)

표준 분포의 하자드 함수를 계산합니다.

지수 함수 (Exponential Functions)

이 단원에서 기술하는 함수들은 헤더 파일 gsl_sf_exp.h 에 정의되어 있습니다.

지수 함수 (Exponential Function)

double gsl_sf_exp(double x)
int gsl_sf_exp_e(double x, gsl_sf_result *result)

지수 함수 \(\exp(x)\) 의 값을 계산합니다. GSL semantics와 오차 검사를 함께 진행합니다.

int gsl_sf_exp_e10_e(double x, gsl_sf_result_e10 *result)

지수 함수 \(\exp(x)\) 의 값을 계산하는 데, 반환값의 자료형으로 gsl_sf_result_e10 을 사용해 확장된 크기의 반환값을 계산합니다. 이 함수는 \(double\) 자료형의 범주를 초과한 \(\exp(x)\) 값을 구할 때, 사용할 수 있습니다.

double gsl_sf_exp_mult(double x, double y)
int gsl_sf_exp_mult_e(double x, double y, gsl_sf_result *result)

주어진 실수 \(x\) 지수 함수 값에 계수 \(y\) 곱한 값 \(y \exp(x)\) 를 계산합니다.

int gsl_sf_exp_mult_e10_e(const double x, const double y, gsl_sf_result_e10 *result)

확장된 범위를 가지는 gsl_sf_result_e10 자료형을 반환 값에 사용해 \(y \exp(x)\) 값을 계산합니다.

상대 지수 함수 (Relative Exponential Functions)

double gsl_sf_expm1(double x)
int gsl_sf_expm1_e(double x, gsl_sf_result *result)

\(\exp(x)-1\) 를 계산합니다. 이 함수에 사용된 알고리즘은 작은 \(x\) 에서만 정확합니다.

double gsl_sf_exprel(double x)
int gsl_sf_exprel_e(double x, gsl_sf_result *result)

작은 \(x\) 값에서 정확한 알고리즘을 이용해 \((\exp(x)-1) /x\) 값을 계산합니다. 작은 \(x\) 값에 대해, 알고리즘은 다음과 같은 확장을 이용합니다.

\[(\exp(x)-1) /x = 1 + x/2 + x^2/(2 \cdot 3) + x^3/(2 \cdot 3 \cdot 4) + \dots\]
double gsl_sf_exprel_2(double x)
int gsl_sf_exprel_2_e(double x, gsl_sf_result *result)

작은 \(x\) 값에서 정확한 알고리즘을 이용해 \(2(\exp(x)-1 -x) /x^2\) 값을 계산합니다. 작은 \(x\) 값에 대해, 알고리즘은 다음과 같은 확장을 이용합니다.

\[2(\exp(x)-1 -x) /x^2 = 1 + x/3 + x^2/(3 \cdot 4) + x^3/(3 \cdot 4 \cdot 5) + \dots\]
double gsl_sf_exprel_n(int n, double x)
int gsl_sf_exprel_n_e(int n, double x, gsl_sf_result *result)

\(N-\) 상대 지수 함수 값을 계산합니다. 이 함수는 gsl_sf_exprel()gsl_sf_exprel_2()\(n\) 차 일반화 함수입니다. 다음과 같이 주어집니다.

\[\begin{split}\begin{flalign} \text{exprel}_N (x) &= N!/x^N (\exp(x) - \sum_{k=0}^{N-1} x^k/k!) \\\\ &= 1+x/(N+1) + x^2/((N+1)(N+2)) + \dots \\\\ &= {}_1F_1(1,1+N,x) \end{flalign}\end{split}\]

오차 평가가 있는 지수 함수 (Exponentiation With Error Estimate)

int gsl_sf_exp_err_e(double x, double dx, gsl_sf_result *result)

주어진 \(x\) 의 지수함수 값을 절대 오차 \(dx\) 와 함께 반환합니다.

int gsl_sf_exp_err_e10_e(double x, double dx, gsl_sf_result_e10 *result)

주어진 \(x\) 의 지수함수 값을 절대 오차 \(dx\) 와 함께 반환합니다. 이때, 반환 값의 자료형을 gsl_sf_result_e10 을 사용해 확장된 범위의 결과를 계산할 수 있습니다.

int gsl_sf_exp_mult_err_e(double x, double dx, double y, double dy, gsl_sf_result *result)

\(x\)\(y\) 대해, \(y \exp(x)\) 값을 절대 오차 \(dx\)\(dy\) 함께 계산합니다.

int gsl_sf_exp_mult_err_e10_e(double x, double dx, double y, double dy, gsl_sf_result_e10 *result)

\(x\)\(y\) 대해, \(y \exp(x)\) 값을 절대 오차 \(dx\)\(dy\) 함께 계산합니다. 이때, 반환 값의 자료형을 gsl_sf_result_e10 을 사용해 확장된 범위의 결과를 계산할 수 있습니다.

지수 적분 함수 (Exponential Integrals)

지수 적분의 자세한 정보는 Abramowitz & Stegun, Chapter 5. 에서 찾아볼 수 있습니다. 이 함수들은 헤더 파일 gsl_sf_expint.h 정의되어 있습니다.

지수 적분 (Exponential Integral)

double gsl_sf_expint_E1(double x)
int gsl_sf_expint_E1_e(double x, gsl_sf_result *result)

지수 적분 함수 \(E_1(x)\) 를 계산합니다.

\[E_1(x) := \mathfrak{R} \int_1^\infty \frac{\exp(-xt)}{t} dt\]
double gsl_sf_expint_E2(double x)
int gsl_sf_expint_E2_e(double x, gsl_sf_result *result)

\(2\) 차 지수 적분 함수 \(E_2(x)\) 를 계산합니다.

\[E_2(x) := \mathfrak{R} \int_1^\infty \frac{\exp(-xt)}{t^2} dt\]
double gsl_sf_expint_En(int n, double x)
int gsl_sf_expint_En_e(int n, double x, gsl_sf_result *result)

\(n\) 차 지수 적분 함수 \(E_2(x)\) 를 계산합니다.

\[E_n(x) := \mathfrak{R} \int_1^\infty \frac{\exp(-xt)}{t^n} dt\]

Ei(x)

double gsl_sf_expint_Ei(double x)
int gsl_sf_expint_Ei_e(double x, gsl_sf_result *result)

지수 적분 \(Ei(x)\) 의 값을 계산합니다.

\[Ei(x) = - PV (\int_{-x}^\infty \frac{\exp(-t)}{t} dt)\]

\(PV\) 는 적분 주요값(Principal Value of Integral)입니다 1 .

초기하 적분 (Hyperbolic Integrals)

double gsl_sf_Shi(double x)
int gsl_sf_Shi_e(double x, gsl_sf_result *result)

다음의 적분값을 계산합니다.

\[\text{Shi} (x) := \int_0^x \frac{\sinh(t)}{t} dt\]
double gsl_sf_Chi(double x)
int gsl_sf_Chi_e(double x, gsl_sf_result *result)

다음의 적분값을 계싼합니다.

\[\text{Chi}(x) := \mathfrak{R} [\gamma_E + \log(x) + \int_0^x \frac{\cosh(t) -1}{t} dt]\]

\(\gamma_E\) 는 오일러 상수입니다. 오일러 상수는 매크로 M_EULER 로 라이브러리 내에 있습니다.

Ei_3(x)

double gsl_sf_expint_3(double x)
int gsl_sf_expint_3_e(double x, gsl_sf_result *result)

다음의 \(3\) 차 지수 적분값을 \(x \geq 0\) 에 대해 계산합니다.

\[\text{Ei}_3 (x) := \int_0^x \exp(-t^3) dt\]

삼각 적분 (Trigonometric Integrals)

double gsl_sf_Si(const double x)
int gsl_sf_Si_e(double x, gsl_sf_result *result)

다음의 적분값을 계산합니다.

\[\text{Si} (x) := \int_0^x \frac{\sin(t)}{t} dt\]
double gsl_sf_Ci(const double x)
int gsl_sf_Ci_e(double x, gsl_sf_result *result)

다음의 적분값을 \(x \geq 0\) 에 대해 계산합니다.

\[\text{Ci} (x) := -\int_0^x \frac{\cos(t)}{t} dt\]

역탄젠트 적분 (Arctangent Integral)

double gsl_sf_atanint(double x)
int gsl_sf_atanint_e(double x, gsl_sf_result *result)

다음의 적분값을 계산합니다.

\[\text{AtanInt}(x) := \int_0^x \frac{\text{arctan}}{t} dt\]

각주

1

코시 주요값(Cauchy principal value)이라고도 합니다(*).

페르미 디렉 함수 (Fermi-Dirac Function)

헤더 파일 gsl_sf_fermi_dirac.h 에 정의되어 있습니다.

완비 페르미-디렉 적분 (Complete Fermi-Dirac Integrals)

완비 페르미 디렉 적분 \(F_j(x)\) 는 다음과 같이 정의됩니다.

\[F_j(x) := \frac{1}{\Gamma(j+1)} \int_0^\infty \frac{t^j}{(\exp(t-x)+1)} dt\]

다른 문헌에서 정규화 계수 없이 표현되기도 합니다.

double gsl_sf_fermi_dirac_m1(double x)
int gsl_sf_fermi_dirac_m1_e(double x, gsl_sf_result *result)

\(j=-1\) 인 완비 페르미 디렉 적분 값을 계산합니다. 이 값은 \(F_{-1}(x) = e^x/(1+e^x)\) 로 주어집니다.

double gsl_sf_fermi_dirac_0(double x)
int gsl_sf_fermi_dirac_0_e(double x, gsl_sf_result *result)

\(j=0\) 인 완비 페르미 디렉 적분 값을 계산합니다. 이 값은 \(F_{0}(x) = \ln(1+e^x)\) 로 주어집니다.

double gsl_sf_fermi_dirac_1(double x)
int gsl_sf_fermi_dirac_1_e(double x, gsl_sf_result *result)

\(j=-1\) 인 완비 페르미 디렉 적분 값을 계산합니다. 이 값은 \(F_{1}(x) = \int_0^\infty(t / (\exp(t-x)+1) dt\) 로 주어집니다.

double gsl_sf_fermi_dirac_2(double x)
int gsl_sf_fermi_dirac_2_e(double x, gsl_sf_result *result)

\(j=2\) 인 완비 페르미 디렉 적분 값을 계산합니다. 이 값은 \(F_{2}(x) = (1/2) \int_0^\infty(t^2 / (\exp(t-x)+1) dt\) 로 주어집니다.

double gsl_sf_fermi_dirac_int(int j, double x)
int gsl_sf_fermi_dirac_int_e(int j, double x, gsl_sf_result *result)

일반 완비 페르미 디렉 적분 값을 계산합니다. 이 값은 \(F_{j}(x) = ((1/\Gamma(j+1)) \int_0^\infty (t^j / (\exp(t-x)+1) dt\) 로 주어집니다.

double gsl_sf_fermi_dirac_mhalf(double x)
int gsl_sf_fermi_dirac_mhalf_e(double x, gsl_sf_result *result)

완비 페르미-디렉 적분 \(F_{-1/2}(x)\) 값을 계산합니다.

double gsl_sf_fermi_dirac_half(double x)
int gsl_sf_fermi_dirac_half_e(double x, gsl_sf_result *result)

완비 페르미-디렉 적분 \(F_{1/2}(x)\) 값을 계산합니다.

double gsl_sf_fermi_dirac_3half(double x)
int gsl_sf_fermi_dirac_3half_e(double x, gsl_sf_result *result)

완비 페르미-디렉 적분 \(F_{3/2}(x)\) 값을 계산합니다.

비완비 페르미-디렉 적분 (Incomplete Fermi-Dirac Integrals)

비완비 페르미-디렉 함수 \(F_j(x,b)\) 는 다음과 같이 정의됩니다.

\[F_j(x,b) := \frac{1}{\Gamma (j+1)} \int_b^\infty \frac{t^j}{\exp(t-x)+1} dt\]
double gsl_sf_fermi_dirac_inc_0(double x, double b)
int gsl_sf_fermi_dirac_inc_0_e(double x, double b, gsl_sf_result *result)

\(0\) 차수의 비완비 페르미-디렉 적분, \(F_0 (x,b) = \ln (1+ e^{b-x}) - (b-x)\) 값을 계산합니다.

감마와 베타 함수 (Gamma and Beta Functions)

다음 함수들은 완전/불완전 감마, 베타 함수를 다양한 팩토리얼에 대해 계산합니다. 이 단원의 함수들은 헤더 파일 gsl_sf_gamma.h 에 정의되어 있습니다.

감마 함수 (Gamma Functions)

감마 함수는 다음의 적분으로 정의되어 있습니다.

\[\Gamma (x) = \int_0^\infty t^{x-1} \exp(-t) dt\]

양의 정수 \(n\) 에 대해, 팩토리얼과 \(\Gamma (n) = (n-1)!\) 의 관계를 가집니다. 감마 함수에 관한 더 자세한 설명은 Abramowitz & Stegun, Chapter 6를 참고할 수 있습니다.

double gsl_sf_gamma(double x)
int gsl_sf_gamma_e(double x, gsl_sf_result *result)

\(0\) 이나 음의 정수가 아닌 \(x\) 에 대해, 감마 함수 \(\Gamma(x)\) 의 값을 계산합니다. 실수 Lanczos 방법을 사용합니다. \(x\) 의 최댓값은 \(\Gamma(x)\) 가 오버 플로우 되지 않는 범위로 매크로 GSL_SF_GAMMA_XMAX 로 주어져 있습니다. 이 값은 \(171.0\) 입니다.

double gsl_sf_lngamma(double x)
int gsl_sf_lngamma_e(double x, gsl_sf_result *result)

\(0\) 이나 음의 정수가 아닌 \(x\) 에 대해, 로그 감마 함수 \(\log(\Gamma(x))\) 의 값을 계산합니다. \(x<0\) 에 대해, \(\log(\Gamma(x))\) 값이 반환되며, 이 값은 \(\log(|\Gamma(x)|)\) 의 값과 같습니다. 실수 Lanczos 방법을 사용합니다.

int gsl_sf_lngamma_sgn_e(double x, gsl_sf_result *result_lg, double *sgn)

\(0\) 이나 음의 정수가 아닌 \(x\) 에 대해, 감마 함수의 부호와 그 크기의 로그 값을 계산합니다. 실수 Lanczos 방법을 사용합니다. 감마 함수의 값과 그 오차는 result_lg 두 요소와 다음의 관계를 사용해 계산합니다. \(\Gamma(x) = sgn * \exp(result_{lg})\) .

double gsl_sf_gammastar(double x)
int gsl_sf_gammastar_e(double x, gsl_sf_result *result)

regulated 감마 함수 \(\Gamma^* (x)\)\(x>0\) 에 대해 계산합니다. 다음과 같이 주어집니다.

\[\begin{split}\Gamma^* (x) = \Gamma (x)/(\sqrt{2 \pi}x^{(x-1/2)} \exp(-x))\\ = (1+ \frac{1}{12x}+ \dots) \text{ for } x \rightarrow \infty\end{split}\]
double gsl_sf_gammainv(double x)
int gsl_sf_gammainv_e(double x, gsl_sf_result *result)

실수 Lanczos 방법을 사용해 감마 함수의 역수, \(1/\Gamma(x)\) 값을 계산합니다.

int gsl_sf_lngamma_complex_e(double zr, double zi, gsl_sf_result *lnr, gsl_sf_result *arg)

복소수 Lanczos 방법을 사용해, 주어진 복소수 \(z= z_r + i z_i\) 에 대해, \(\log(\Gamma(z))\) 값을 계산합니다. 이때, \(z\)\(0\) 과 음의 정수가 아닌 복소수입니다. lnr :math:` = log|Gamma(z)|` , arg :math:` = text{arg}(Gamma(z))` 이고, arg\((-\pi, \pi]\) 범위의 값을 가집니다. \(|z|\) 가 매우 클 때, arg 값이 정확히 정해지지 않을 수 있습니다. 이는 \((-\pi, \pi]\) 로 범위를 제약함으로써 생기는 절단점으로 인해 필연적으로 발생합니다. 이러한 상황은 GSL_ELOSS 오류에 속하는 상황입니다. 절대값 lnr 부분은 이러한 정밀도 저하를 격지 않습니다.

팩토리얼

양의 정수 \(n\) 에 대해, 감마 함수를 이용해 \(n! = \Gamma(n+1)\) 팩토리얼을 계산할 수 있습니다. 하지만, 아래의 함수들을 이용하는 것이 작은 \(n\) 값들에 대해 더 효율적입니다. 팩토리얼 값들을 하드 코딩된 테이블에 보관하고 있습니다.

double gsl_sf_fact(unsigned int n)
int gsl_sf_fact_e(unsigned int n, gsl_sf_result *result)

팩토리얼 \(n!\) 를 계산합니다. 팩토리얼은 감마 함수와 \(n! = \Gamma(n+1)\) 의 관계를 가지고 있습니다. \(n\) 의 최댓값은 \(n!\) 이 오버플로우되지 않는 값으로 정해집니다. 이는 매크로 GSL_SF_FACT_NMAX 에 정의되어 있고 \(170\) 입니다.

double gsl_sf_doublefact(unsigned int n)
int gsl_sf_doublefact_e(unsigned int n, gsl_sf_result *result)

더블 팩토리얼 \(n!! = n(n-2)(n-4)\dots\) 을 계산합니다. \(n\) 의 최댓값은 \(n!!\) 이 오버플로 되지 않는 값으로 정해집니다. 이는 매크로 GSL_SF_DOUBLEFACT_NMAX 에 정의되어 있고 \(297\) 입니다.

double gsl_sf_lnfact(unsigned int n)
int gsl_sf_lnfact_e(unsigned int n, gsl_sf_result *result)

\(n\) 팩토리얼의 로그 값, \(\log(n!)\) 값을 계산합니다. 이 알고리즘은 \(n <170\) 에서 \(\ln(\Gamma(n+1))\) 값을 계산하는 gsl_sf_lngamma 보다 빠릅니다. 하지만 큰 \(n\) 에 대해서는 빠르지 않습니다.

double gsl_sf_lndoublefact(unsigned int n)
int gsl_sf_lndoublefact_e(unsigned int n, gsl_sf_result *result)

\(n\) 대해, 더블 팩토리얼의 로그 값 \(\log(n!!)\) 을 계산합니다.

double gsl_sf_choose(unsigned int n, unsigned int m)
int gsl_sf_choose_e(unsigned int n, unsigned int m, gsl_sf_result *result)

조합 계수 n choose m \(= n!/(m!(n-m!))\) 의 값을 계산합니다.

double gsl_sf_lnchoose(unsigned int n, unsigned int m)
int gsl_sf_lnchoose_e(unsigned int n, unsigned int m, gsl_sf_result *result)

조합 계수 \(n choose m\) 의 로그 값을 계산합니다. 이 값은 \(\log(n!) - \log(m!) - \log((n-m)!)\) 과 같습니다.

double gsl_sf_taylorcoeff(int n, double x)
int gsl_sf_taylorcoeff_e(int n, double x, gsl_sf_result *result)

\(x \geq0\) , \(n \geq0\) 에 대해, 테일러 계수 \(x^n/n!\) 값을 계산합니다.

포흐하머 기호

double gsl_sf_poch(double a, double x)
int gsl_sf_poch_e(double a, double x, gsl_sf_result *result)

포흐하머 기호 \((a)_x = \Gamma(a_x)/\Gamma(a)\) 를 계산합니다. 포흐하머 기호는 아펠 기호로도 알려져있으며, \((a,x)\) 로 표기하기도 합니다. \(a\)\(a+x\) 가 음의 정수나 \(0\) 일때, 해당 비의 극한 값이 반환됩니다.

double gsl_sf_lnpoch(double a, double x)
int gsl_sf_lnpoch_e(double a, double x, gsl_sf_result *result)

포흐하머 기호의 로그값 \(\log((a)_x) = \log(\Gamma(a+x)/\Gamma(a))\) 을 계산합니다.

int gsl_sf_lnpoch_sgn_e(double a, double x, gsl_sf_result *result, double *sgn)

포흐하머 기호의 부호와 그 크기의 로그값을 계산합니다. 계산되는 계수들은 result \({}= \log(|(a)_x|)\) 가 오차 값과 함께 계산되고, \((a)_x = \Gamma(a+x)/\Gamma(a)\) 에 대해, sgn \({}= \text{sgn}((a)_x)\) 을 계산합니다.

double gsl_sf_pochrel(double a, double x)
int gsl_sf_pochrel_e(double a, double x, gsl_sf_result *result)

\((a)_x = \Gamma(a+x)/\Gamma(a)\) 에 대해, \(((a)_x -1)/x\) 값을 계산합니다.

불완전 감마 함수

double gsl_sf_gamma_inc(double a, double x)
int gsl_sf_gamma_inc_e(double a, double x, gsl_sf_result *result)

실수 \(a\)\(x \geq 0\) 에 대해, 비정규화된 불완전 감마 함수 \(\Gamma(a,x) = \int_x^\infty t^{(a-1)} \exp(-t) dt\) 값을 계산합니다.

double gsl_sf_gamma_inc_Q(double a, double x)
int gsl_sf_gamma_inc_Q_e(double a, double x, gsl_sf_result *result)

\(a>0\)\(x \leq 0\) 에 대해, 정규화된 불완전 감마 함수 \(Q(a,x) = 1.\Gamma(a) \int_x^\infty t^{(a-1)} \exp(-t) dt\) 의 값을 계산합니다.

double gsl_sf_gamma_inc_P(double a, double x)
int gsl_sf_gamma_inc_P_e(double a, double x, gsl_sf_result *result)

\(a>0\)\(x \geq 0\) 에 대해, \(P(a,x) = 1-Q(a,x) = 1/\Gamma(a) \int_0^x t^{(a-1)} \exp(-t) dt\) 값을 계산합니다.

Abramowtz & Stegun의 6.5단원, 불완전 감마 함수에서 \(P(a,x)\) 표기를 씁니다.

베타 함수

double gsl_sf_beta(double a, double b)
int gsl_sf_beta_e(double a, double b, gsl_sf_result *result)

베타함수 \(B(a,b) = \Gamma(a)\Gamma(b)/\Gamma(a+b)\) 값을 계산합니다. \(a,b\) 는 음의 정수가 아니여야 합니다.

double gsl_sf_lnbeta(double a, double b)
int gsl_sf_lnbeta_e(double a, double b, gsl_sf_result *result)

베타 함수의 로그 값 \(\log(B(a,b))\) 를 계산합니다. \(a,b\) 는 음의 정수가 아니여야 합니다.

불완전 베타 함수

double gsl_sf_beta_inc(double a, double b, double x)
int gsl_sf_beta_inc_e(double a, double b, double x, gsl_sf_result *result)

정규화된 불완전 베타함수 \(I_x (a,b) = B_x(a,b)/ B(a,b)\) 를 계산합니다. \(B_x(a,b)\)\(0 \leq x \leq 1\) 에 대해 다음과 같이 정해집니다.

\[B_x (a,b) = \int_0^x t^{a-1} (1-t)^{b-1} dt\]

이 값은 \(a>0, b>0\) 에 대해, 연속 분수 전개를 이용해 계산됩니다. 다른 경우에는 다음의 관계를 이용해 계산합니다.

\[I_x (a,b,x) = (\frac{1}{a}) x^a \frac{_2F_1 (a, 1-b, a+1, x)}{B(a,b)}\]

구겐바우어 함수 (Gegenbauer Functions)

구겐바우어 다항식은 Abramowitz & Stgun의 22단원에 정의되어 있습니다. 이 다항식은 또 Ultraspherical 다항식으로도 알려져있습니다. 이 함수들은 헤더 파일 gsl_sf_gegenbauer.h 에 정의되어 있습니다.

double gsl_sf_gegenpoly_1(double lambda, double x)
double gsl_sf_gegenpoly_2(double lambda, double x)
double gsl_sf_gegenpoly_3(double lambda, double x)
int gsl_sf_gegenpoly_1_e(double lambda, double x, gsl_sf_result *result)
int gsl_sf_gegenpoly_2_e(double lambda, double x, gsl_sf_result *result)
int gsl_sf_gegenpoly_3_e(double lambda, double x, gsl_sf_result *result)

구겐바우어 다항식 \(C_n^{(\lambda)}(x)\)\(n= 1, 2, 3\) 인 경우에 대해, 정의식을 이용해 계산합니다.

double gsl_sf_gegenpoly_n(int n, double lambda, double x)
int gsl_sf_gegenpoly_n_e(int n, double lambda, double x, gsl_sf_result *result)

구겐바우어 다항식 \(C_n^{(\lambda)}(x)\) 을 주어진 \(n\) , lambda , 그리고 \(x\) 에 대해 계산합니다. 이때, \(\lambda > - \frac{1}{2}. n \geq 0\) 이어야 합니다.

int gsl_sf_gegenpoly_array(int nmax, double lambda, double x, double result_array[])

구겐바우어 다항식 \(C_n^{(\lambda)}(x)\) 배열 값을 계산합니다. \(n = 0, 1, 2, \dots , nmax\) 의 값을 계산하며, \(\lambda > - \frac{1}{2}. nmax \geq 0\) 의 제약을 가집니다.

에르미트 다항식과 함수 (Hermite Polynomials and Functions)

에르미트 다항식과 함수는 Abramowitz & Stegunm Chapter 22 와 Szego, Gabor (1939, 1957, 1967) Orthogonal Polynomials, American Mathematical Society에 기술되어 있습니다. 본 단원의 함수들은 헤더 파일 gsl_sf_hermite.h 에 정의되어 있습니다.

에르미트 다항식(Hermite Polynomials)

에르미트 다항식은 두 가지 형태가 존재합니다. \(H_n(x)\) 는 물리학에서 사용하는 형태이고, \(H_{e_n}(x)\) 는 확률론에서 사용하는 형태입니다.

\[\begin{split}H_n(x) = (-1)^n e^{x^2} (\frac{d}{dx})^n e^{-x^2}\\ H_{e_n}(x) = (-1)^n e^{x^2/2} (\frac{d}{dx})^n e^{-x^2/2}\end{split}\]

이 둘은 다음의 관계를 가지고,

\[\begin{split}H_n(x) = 2^{\frac{n}{2}}H_{e_n}(\sqrt{2}x)\\ H_{e_n}(x) = 2^{-\frac{n}{2}}H_n(\frac{x}{\sqrt{2}})\end{split}\]

다음과 같은 미분 방정식을 만족합니다.

\[\begin{split}H_{n}''(x) -2xH_{n}'(x) + 2nH_{n}(x)=0\\ H_{e_n}''(x) -xH_{e_n}'(x) + nH_{e_n}(x)=0\end{split}\]
double gsl_sf_hermite(const int n, const double x)
int gsl_sf_hermite_e(const int n, const double x, gsl_sf_result *result)

이 함수는 :math”H_n (x) 형태의 에르미트 다항식을 주어진 차수 n 과 변수 x 에 대해 계산합니다. 오버플로우가 감지되면, 오류 관리자를 호출하지 않고 GSL_EOVERFLW 를 반환합니다.

int gsl_sf_hermite_array(const int nmax, const double x, double *result_array)

\(nmax\) 이하의 차수를 가지는 모든 \(H_n(x)\) 형태의 에르미트 다항식을 주어진 변수 \(x\) 에 대해 계산합니다. 결과는 result_array 에 저장됩니다.

double gsl_sf_hermite_series(const int n, const double x, const double *a)
int gsl_sf_hermite_series_e(const int n, const double x, const double *a, gsl_sf_result *result)

\(\sum_{j=0}^n a_j H_j (x)\) 급수 값을 Clenshaw 알고리즘을 이용해 계산합니다.

double gsl_sf_hermite_prob(const int n, const double x)
int gsl_sf_hermite_prob_e(const int n, const double x, gsl_sf_result *result)

\(H_{e_n}(x)\) 형태의 에르미트 다항식을 주어진 차수 \(n\) 변수 \(x\) 대해 계산합니다. 오버플로우가 감지되면, 오류 관리자를 호출하지 않고 GSL_EOVERFLW 를 반환합니다.

int gsl_sf_hermite_prob_array(const int nmax, const double x, double *result_array)

\(nmax\) 이하의 차수를 가지는 모든 \(H_{e_n}(x)\) 형태의 에르미트 다항식을 주어진 변수 \(x\) 에 대해 계산합니다. 결과는 result_array 에 저장됩니다.

double gsl_sf_hermite_prob_series(const int n, const double x, const double *a)
int gsl_sf_hermite_prob_series_e(const int n, const double x, const double *a, gsl_sf_result *result)

\(\sum_{j=0}^n a_j H_{e_j} (x)\) 급수 값을 Clenshaw 알고리즘을 이용해 계산합니다.

에르미트 다항식의 도함수 (Derivatives of Hermite Polynomials)

double gsl_sf_hermite_deriv(const int m, const int n, const double x)
int gsl_sf_hermite_deriv_e(const int m, const int n, const double x, gsl_sf_result *result)

\(n\) 차수의 에르미트 다항식 \(H_n(x)\)\(m\) 차 도함수 값을 주어진 변수 \(x\) 대해 계산합니다.

int gsl_sf_hermite_array_deriv(const int m, const int nmax, const double x, double *result_array)

\(0, \dots, nmax\) 차수의 모든 에르미트 다항식 \(H_n(x)\)\(m\) 차 도함수 값을 주어진 변수 \(x\) 대해 계산합니다. \(d^m / dx^m H_n(x)\) 의 값은 result_array[n] 에 저장됩니다. 계산 결과가 저장되는 result_array 는 최소 \(nmax+1\) 이상의 길이를 가져야 합니다.

int gsl_sf_hermite_deriv_array(const int mmax, const int n, const double x, double *result_array)

\(n\) 차수를 가지는 에르미트 다항식 \(H_n(x)\) 의 모든 \(0, \dots, mmax\) 차 도함수 값을 주어진 변수 \(x\) 대해 계산합니다. \(d^m / dx^m H_n(x)\) 의 값은 result_array[m] 에 저장됩니다. 계산 결과가 저장되는 result_array 는 최소 mmax+1 이상의 길이를 가져야 합니다.

double gsl_sf_hermite_prob_deriv(const int m, const int n, const double x)
int gsl_sf_hermite_prob_deriv_e(const int m, const int n, const double x, gsl_sf_result *result)

\(n\) 차수의 에르미트 다항식 \(H_{e_n}(x)\)\(m\) 차 도함수 값을 주어진 변수 \(x\) 대해 계산합니다.

int gsl_sf_hermite_prob_array_deriv(const int m, const int nmax, const double x, double *result_array)

\(n\) 차수를 가지는 에르미트 다항식 \(H_{e_n}(x)\) 의 모든 \(0, \dots, mmax\) 차 도함수 값을 주어진 변수 \(x\) 대해 계산합니다. \(d^m / dx^m H_{e_n}(x)\) 의 값은 result_array[m] 에 저장됩니다. 계산 결과가 저장되는 result_array 는 최소 mmax+1 이상의 길이를 가져야 합니다.

int gsl_sf_hermite_prob_deriv_array(const int mmax, const int n, const double x, double *result_array)

\(n\) 차수를 가지는 에르미트 다항식 \(H_{e_n}(x)\) 의 모든 \(0, \dots, mmax\) 차 도함수 값을 주어진 변수 \(x\) 대해 계산합니다. \(d^m / dx^m H_{e_n}(x)\) 의 값은 result_array[m] 저장됩니다. 계산 결과가 저장되는 result_array 는 최소 mmax+1 이상의 길이를 가져야 합니다.

에르미트 함수 (Hermite Functions)

에르미트 함수는 다음과 같이 정의됩니다.

\[\psi_n(x) = \frac{1}{(2^n n! \sqrt{\pi})^{\frac{1}{2}}} e^{- \frac{x^2}{2}} H_n(x)\]

그리고 이는 양자 역학에 나오는 슈뢰딩거 방정식의 조화 진동자 형태를 만족합니다.

\[\psi_n''(x) + (2n+1-x^2)\psi_n(x) =0\]

이 들은 서로 직교하는 함수고

\[\int_{-\infty}^{\infty} \psi_m(x)\psi_n(x) \, dx = \delta_{mn}\]

\(L^2 (\mathbb{R})\) 공간의 직교 기저를 형성합니다. 에르미트 함수들은 연속 푸리에 변환의 고유 함수이기도 합니다. GSL은 에르미트 함수를 계산하는 두 가지 방법을 제공합니다. 첫 번째는 수학적으로 정의된 \(3\) 개 항의 재귀 관계를 이용합니다. 이 방법은 \(O(n)\) 의 계산 복잡도를 가지고 가장 정확합니다. 두 번째는 코시 적분 접근 방법을 이용한 방법입니다. 이는 (Bunck, 2009)에 소개 되었으며, \(O(\sqrt{n})\) 의 계산복잡도를 가집니다. 정확도를 조금 희생하지만 \(n\) 값이 클 수록, 기존 방법에 비해 속도에 큰 이점이 있습니다.

double gsl_sf_hermite_func(const int n, const double x)
int gsl_sf_hermite_func_e(const int n, const double x, gsl_sf_result *result)

차수 \(n\) 에르미트 함수 \(\psi_n(x)\) 를 주어진 변수 \(x\) 에 대해 계산합니다. 이 방법은 재귀 관계를 이용하며, \(O(n)\) 의 계산 복잡도를 가집니다.

double gsl_sf_hermite_func_fast(const int n, const double x)
int gsl_sf_hermite_func_fast_e(const int n, const double x, gsl_sf_result *result)

차수 \(n\) 에르미트 함수 \(\psi_n(x)\) 를 주어진 변수 \(x\) 에 대해 계산합니다. 이 방법은 (Bunck, 2009)의 코시 적분을 이용하며, \(O(\sqrt{n})\) 의 계산 복잡도를 가집니다.

int gsl_sf_hermite_func_array(const int nmax, const double x, double *result_array)

\(n=0, \dots ,nmax\) 의 차수를 가지는 에르미트 함수 \(\psi_n(x)\) 를 주어진 변수 \(x\) 대해, 재귀적 방법을 이용해 계산합니다. 계산 결과는 result_array 저장되며 최소 \(nmax+1\) 이상의 길이를 가져야 합니다.

double gsl_sf_hermite_func_series(const int n, const double x, const double *a)
int gsl_sf_hermite_func_series_e(const int n, const double x, const double *a, gsl_sf_result *result)

\(\sum_{j=0}^n a_j \psi_j (x)\) 급수를 계산합니다. \(\psi_j\)\(j\) 의 차수를 가지는 에르미트 함수를 의미하며, Clenshaw 알고리즘을 이용합니다.

에르미트 함수의 도함수 (Derivatives of Hermite Functions)

double gsl_sf_hermite_func_der(const int m, const int n, const double x)
int gsl_sf_hermite_func_der_e(const int m, const int n, const double x, gsl_sf_result *result)

\(n\) 차수의 에르미트 함수 \(\psi_n(x)\)\(m\) 차 도함수를 주어진 \(x\) 대해 계산합니다.

에르미트 함수와 다항식의 근 (Zeros of Hermite Polynomials and Hermite Functions)

이 함수들은 차수 \(n\) 을 가지는 에르미트 함수와 다항식의 \(s\) 번째 근을 계산합니다. 각 근들이 원점을 기준으로 대칭이기 때문에, 양수인 근들만 계산됩니다. 인덱스는 \(1\) 부터 시작해서 오름차순으로 배열됩니다. 홀수 차수의 다항식 만이 \(0\) 에서 \(0\) 번째 근을 가집니다. 해당 값은 항상 \(0\) 입니다.

double gsl_sf_hermite_zero(const int n, const int s)
int gsl_sf_hermite_zero_e(const int n, const int s, gsl_sf_result *result)

\(n\) 차수의 에르미트 다항식 \(H_n(x)\)\(s\) 번째 근을 계산합니다.

double gsl_sf_hermite_prob_zero(const int n, const int s)
int gsl_sf_hermite_prob_zero_e(const int n, const int s, gsl_sf_result *result)

\(n\) 차수의 에르미트 다항식 \(H_{e_n}(x)\)\(s\) 번째 근을 계산합니다.

double gsl_sf_hermite_func_zero(const int n, const int s)
int gsl_sf_hermite_func_zero_e(const int n, const int s, gsl_sf_result *result)

\(n\) 차수의 에르미트 함수 \(\psi_n(x)\)\(s\) 번째 근을 계산합니다.

초기하 함수 (Hypergeometric Functions)

초기하 함수들은 Abramowitz& Stegun의 13, 15 단원의 기술을 기반으로 작성되었습니다. 헤더 파일 gsl_sf_hyperg.h 에 정의 되어있습니다,

double gsl_sf_hyperg_0F1(double c, double x)
int gsl_sf_hyperg_0F1_e(double c, double x, gsl_sf_result *result)

초기하 함수

\[{}_0F_1(c,x)\]

를 계산합니다.

double gsl_sf_hyperg_1F1_int(int m, int n, double x)
int gsl_sf_hyperg_1F1_int_e(int m, int n, double x, gsl_sf_result *result)

합류 초기하 함수(confluent hypergeometric)

\[{}_1F_1(m,n,x) = M(m,n,x)\]

를 정수 인자 \(m\)\(n\) 에 따라 계산합니다.

double gsl_sf_hyperg_1F1(double a, double b, double x)
int gsl_sf_hyperg_1F1_e(double a, double b, double x, gsl_sf_result *result)

합류 초기하 함수

\[_1F_1(a,b,x) = M(a,b,x)\]

를 일반 인자 \(a\) \(b\) 따라 계산합니다.

double gsl_sf_hyperg_U_int(int m, int n, double x)
int gsl_sf_hyperg_U_int_e(int m, int n, double x, gsl_sf_result *result)

합류 초기하 함수 \(U(m,n,x)\) 를 정수 인자 \(m\)\(n\) 에 대해 계산합니다.

int gsl_sf_hyperg_U_int_e10_e(int m, int n, double x, gsl_sf_result_e10 *result)

합류 초기하 함수 \(U(m,n,x)\) 를 정수 인자 \(m\)\(n\) 에 대해 계산하고 확장된 구간에 대해, gsl_sf_result_e10 형의 값을 반환합니다.

double gsl_sf_hyperg_U(double a, double b, double x)
int gsl_sf_hyperg_U_e(double a, double b, double x, gsl_sf_result *result)

합류 초기하 함수 \(U(a,b,x)\) 를 계산합니다.

int gsl_sf_hyperg_U_e10_e(double a, double b, double x, gsl_sf_result_e10 *result)

합류 초기하 함수 \(U(a,b,x)\) 를 계산하고 확장된 구간에 대해 gsl_sf_result_e10 형의 값을 반환합니다.

double gsl_sf_hyperg_2F1(double a, double b, double c, double x)
int gsl_sf_hyperg_2F1_e(double a, double b, double c, double x, gsl_sf_result *result)

구간 \(\|x\|<1\) 에 대해, 가우스 초기하 함수

\[_2F_1(a,b,c,x) = F(a,b,c,x)\]

의 값을 계산합니다. 만약, 인자 \((a,b,c,x)\) 가 특이점(singular point)에 너무 가깝다면, 급수 근사가 너무 느려지게 되고 함수는 오류 값 GSL_EMAXITER 반환합니다. 이러한 지점은 \(x=1\) , \(c-a-b=m, m \in \mathbf{Z}\) 구간에서 발생합니다.

double gsl_sf_hyperg_2F1_conj(double aR, double aI, double c, double x)
int gsl_sf_hyperg_2F1_conj_e(double aR, double aI, double c, double x, gsl_sf_result *result)

구간 \(\|x\|<1\) 에 대해, 가우스 초기하 함수

\[_2F_1 (a_R + i a_I , aR-iaI, c, x)\]

의 복소수 인자 값을 계산합니다.

double gsl_sf_hyperg_2F1_renorm(double a, double b, double c, double x)
int gsl_sf_hyperg_2F1_renorm_e(double a, double b, double c, double x, gsl_sf_result *result)

구간 \(\|x\|<1\) 에 대해, 재규격화 된 가우스 초기하 함수

\[\frac{_2F_1(a,b,c,x)}{\Gamma(c)}\]

의 값을 계산합니다.

double gsl_sf_hyperg_2F1_conj_renorm(double aR, double aI, double c, double x)
int gsl_sf_hyperg_2F1_conj_renorm_e(double aR, double aI, double c, double x, gsl_sf_result *result)

구간 \(\|x\|<1\) 에 대해, 재규격화 된 가우스 초기하 함수

\[\frac{_2F_1(a_R + ia_I,a_R - i a_I,c,x)}{\Gamma(c)}\]

의 값을 계산합니다.

double gsl_sf_hyperg_2F0(double a, double b, double x)
int gsl_sf_hyperg_2F0_e(double a, double b, double x, gsl_sf_result *result)

초기하 함수

\[{}_2F_0(a,b,x)\]

를 계산합니다.

급수 표현은 발산하는 초기하 급수입니다. 하지만, \(x<0\) 이라면 다음을 얻을 수 있습니다.

\[{}_2F_0 (a,b,x) = (-\frac{1}{x})^a U(a, 1+a, -b,- \frac{1}{x})\]

라게르 함수 (Laguerre Functions)

일반화 된 라게르 다항식(다른 이름으로 버금 라게르 다항식이 있습니다)은 합류 초기하 함수(confluent hypergeometric function)로 정의됩니다.

\[L_n^a (x) = \frac{(a+1)_n}{n!} {}_1F_1(-n,a+1,x)\]

\((a)_n\) 는 포흐하머 기호(Pochhammer symbol)입니다. 이들은 일반적인 라게르 다항식 \(L_n(x)\) 과 다음의 관계를 가집니다.

\[\begin{split}L_n^0(x) = L_n(x)\\ L_n^k(x) = (-1)^k (d^k /d x^k)L_{(n+k)}(x)\end{split}\]

더 자세한 정보는 Abramowitz & Strgun, Chapter 22를 참고할 수 있습니다.

이 단원에서 기술된 함수들은 헤더 파일 gsl_sf_laguerre.h 에 정의되어 있습니다.

double gsl_sf_laguerre_1(double a, double x)
double gsl_sf_laguerre_3(double a, double x)
double gsl_sf_laguerre_2(double a, double x)
int gsl_sf_laguerre_1_e(double a, double x, gsl_sf_result *result)
int gsl_sf_laguerre_2_e(double a, double x, gsl_sf_result *result)
int gsl_sf_laguerre_3_e(double a, double x, gsl_sf_result *result)

일반화된 라게르 다항식 \(L_1^a (x), L_2^a (x), L_3^a (x)\) 을 수학 정의식을 이용해 계산합니다.

double gsl_sf_laguerre_n(const int n, const double a, const double x)

일반화된 라게르 다항식 \(L_n^a(x)\)\(a > -1, n \geq 0\) 인 경우를 계산합니다.

람베르트 W 함수 (Lambert W Functions)

람베르트 \(W\) 함수 \(W(x)\)\(W(x)\exp(W(x)) = x\) 방정식의 해로 정의됩니다. 이 함수는 \(x<0\) 에서 다양한 부분 함수들로 나뉘어 집니다. 하지만, \(2\) 개의 실수 함수들이 존재합니다. 일반적으로 \(W_0(x)\) 를 주 함수로 사용합니다. 이 함수는 \(x<0\) 에 대해, \(W>-1\) 값을 가집니다. 그리고, \(w_{-1}(x)\) 는 또다른 실수 함수로 \(x<0\) 에 대해, \(W<-1\) 값을 가집니다.

람베르트 함수들은 헤더 파일 gsl_sf_lambert.h 에 정의되어 있습니다.

double gsl_sf_lambert_W0(double x)
int gsl_sf_lambert_W0_e(double x, gsl_sf_result *result)

주 람베르트 \(W\) 함수 \(W_0(x)\) 값을 계산합니다.

double gsl_sf_lambert_Wm1(double x)
int gsl_sf_lambert_Wm1_e(double x, gsl_sf_result *result)

두 번째 실수 람베르트 \(W\) 함수 \(W_{-1](x)\) 값을 계산합니다.

르장드르 함수와 구면조화 함수 (Legendre Functions and Spherical Harmonics)

르장드르 다항식과 함수들은 Abramowitz & Stegun Chapter 8. 에 기술되어있습니다. gsl_sf_legnedre.h 헤더 파일에 정의되어 있습니다.

르장드르 다항식(Legendre Polynomials)

double gsl_sf_legendre_P1(double x)
double gsl_sf_legendre_P3(double x)
double gsl_sf_legendre_P2(double x)
int gsl_sf_legendre_P1_e(double x, gsl_sf_result *result)
int gsl_sf_legendre_P2_e(double x, gsl_sf_result *result)
int gsl_sf_legendre_P3_e(double x, gsl_sf_result *result)

르장드르 다항식 \(P_l (x)\)\(l= 1, 2, 3\) 에 대해 계산합니다. 이 계산은 해당 함수들의 수학적 정의식을 사용합니다.

double gsl_sf_legendre_Pl(int l, double x)
int gsl_sf_legendre_Pl_e(int l, double x, gsl_sf_result *result)

르장드르 다항식 \(P_l(x)\) 을 주어진 \(l\)\(x\) 에 대해 계산합니다. 이때, \(l\)\(x\)\(l \geq 0\)\(|x| \leq 1\) 을 만족해야 합니다.

int gsl_sf_legendre_Pl_array(int lmax, double x, double result_array[])
int gsl_sf_legendre_Pl_deriv_array(int lmax, double x, double result_array[], double result_deriv_array[])

\(l= 0, \dots lmax\)\(|x| \leq 1\) 에 대해, 르장드르 다항식 \(P_l(x)\) 과 그 도함수 \(dP_l (x)/dx\) 를 주어진 배열 result_array[] 에 계산합니다.

double gsl_sf_legendre_Q0(double x)
int gsl_sf_legendre_Q0_e(double x, gsl_sf_result *result)

\(x>-1\)\(x \neq 1\) 에 대해, 르장드르 함수 \(Q_0(x)\) 를 계산합니다.

double gsl_sf_legendre_Q1(double x)
int gsl_sf_legendre_Q1_e(double x, gsl_sf_result *result)

\(x>-1\)\(x \neq 1\) 에 대해, 르장드르 함수 \(Q_1(x)\) 를 계산합니다.

double gsl_sf_legendre_Ql(int l, double x)
int gsl_sf_legendre_Ql_e(int l, double x, gsl_sf_result *result)

\(x>-1\) , \(x \neq 1\)\(l \geq 0\) 에 대해 ,르장드르 함수 \(Q_l(x)\) 를 계산합니다.

버금 르장드르 함수와 구면 조화 함수 (Associated Legendre Polynomials and Spherical Harmonics)

이 단원의 함수들은 버금 르장드르 함수 \(P_l^m(x)\) 의 값을 계산합니다.

다음 미분 방정식의 해입니다.

\[(1-x^2) \frac{d^2}{d x^2}P_l^m(x) -2x \frac{d}{dx}P_l^m(x) + (l(l+1) - \frac{m^2}{1-x^2})P_l^m(x) =0\]

\(l\)\(m\)\(0 \leq l\)\(0 \leq m \leq l\) 을 만족합니다. 함수 \(P_l^m(x)\) 는 조합적으로 상승하며 \(l\)\(150\) 보다 클 경우 오버플로우가 발생합니다. 이를 대신해서 정규화 된 버금 르장드르 함수를 계산할 수 있습니다. 다양한 종류의 정규화 표현이 존재하고, degree와 order가 \(2700\) 인 함수까지 안정적으로 계산할 수 있습니다. 이 라이브러리에서는 다음의 정규화 표현을 제공합니다.

  • 슈미트 반 정규화 (Schmidt semi-normalization)

    슈미트 반 정규화 버금 르장드르 함수는 자기 교환 계산에 빈번히 사용됩니다. 이는 다음과 같이 계산할 수 있습니다.

    \[\begin{split}S_l^0 (x) = P_l^0 (x)\\ S_l^m (x) = (-1)^m \sqrt{2 \frac{(l-m)!}{(l+m)!}}P_l^m (x), m > 0\end{split}\]

    계수 \((-1)^m\) 은 Condon-Shortley 위상 계수로 불리며, 필요시 함수에서 \(csphase =1\) 인자를 설정해 무시할 수 있습니다.

  • 구면 조화 정규화 (Spherical Harmonic Normalization)

    버금 르장드르 함수는 다음과 같이 구면 조화 함수를 계산하는 데 쓸 수 있습니다.

    \[Y_l^m(x) = (-1)^m \sqrt{\frac{2l+1}{4 \pi}\frac{(l-m)!}{(l+m)!}} P_l^m(x)\]

    계수 \((-1)^m\) 은 필요시 계산에서 제외할 수 있습니다.

  • 완전 정규화 (Full Normalization)

    완전 정규화된 버금 르장드르 다항식은 다음과 같이 정의됩니다.

    \[N_l^m(x) = (-1)^m \sqrt{(1+ \frac{1}{2}\frac{(l-m)!}{(l+m)!}}P_l^m(x)\]

    이 때, 다음과 같은 성질을 가집니다.

    \[\int_{-1}^1 N_l^m(x) \, dx = 1\]

아래에 나올 정규화된 버금 르장드르 함수를 계산하는 함수들은 재귀적 방법을 사용합니다. 이 방법은 degree \(l\) 과 order \(m\) 이 2700 이하라면, 안정적으로 계산할 수 있습니다. 이 값을 넘어서면 계산 함수들은 언더플로우를 일으켜 부정확한 값을 반환합니다. 각 함수들은 \(1\) -계 도함수 \(dP_l^m(x)/dx\)\(2\) -계 도함수 \(d^2 P_l^m(x)/ dx^2\) 을 제공하며, 이와 함께 \(dP_l^m(\cos\theta) / d \theta\)\(d^2 P_l^m(\cos\theta)/d\theta^2\) 도 같이 제공합니다. 이 두 종류의 도함수들은 단순한 비례 관계를 가집니다. \(\theta\) 에 관한 미분은 구면 조화 함수에서 매우 빈번히 사용되기 때문에 이 기능 또한 같이 제공하고 있습니다.

아래의 함수에서 gsl_sf_legendre_t 인자를 이용해 정규화 방법을 선택할 수 있습니다. 가능한 값들은 다음과 같습니다.

type gsl_sf_legendre_t

설명

GSL_SF_LEGENDRE_NONE

비 정규화된 버금 르장드르 다항식 \(P_l^m(x)\)

GSL_SF_LEGENDRE_SCHMIDT

슈미트 반 정규화된 버금 르장드르 다항식 \(S_l^m(x)\)

GSL_SF_LEGENDRE_SPHARM

구면 조화 버금 르장드르 다항식 \(Y_l^m(x)\)

GSL_SF_LEGENDRE_FULL

완전 정규화 버금 르장드르 다항식 \(N_l^m(x)\)

int gsl_sf_legendre_array(const gsl_sf_legendre_t norm, const size_t lmax, const double x, double result_array[])
int gsl_sf_legendre_array_e(const gsl_sf_legendre_t norm, const size_t lmax, const double x, const double csphase, double result_array[])

정규화된 버금 르장드르 다항식들을, \(0 \leq l \leq lmax\) , \(0 \geq m \geq l\) 그리고 \(|x| \leq 1\) 에 대해 계산합니다. norm 인자는 어느 정규화 방법을 사용할지 결정합니다. 정규화된 \(P_l^m(x)\) 값들은 result_array 에 저장됩니다. 이 값은 \(gsl_sf_lengendre_array_n()\) 를 호출해 최소 크기를 결정할 수 있습니다.

배열 \(P_l^m(x)\) 의 지수는 gsl_sf_legendre_array_index(l, m) 를 호출해 얻을 수 있습니다. _e 붇은 함수에서 Condon-Shortly 위상 계수 \((-1)^m\) 의 포함 유무를 조정하려면 csphase\(-1\) 이나 \(1\) 로 설정해 주면 됩니다. 이 계수는 기본적으로 비활성화 되어 있습니다.

int gsl_sf_legendre_deriv_array(const gsl_sf_legendre_t norm, const size_t lmax, const double x, double result_array[], double result_deriv_array[])
int gsl_sf_legendre_deriv_array_e(const gsl_sf_legendre_t norm, const size_t lmax, const double x, const double csphase, double result_array[], double result_deriv_array[])

\(|x| \leq 1\) 값에 대해, 정규화된 버금 르장드르 함수들의 \(1\) 차에서 \(lmax\) 차 까지의 도함수 값을 계산합니다. norm 인자는 어느 정규화 방법을 사용할지 결정합니다. 르장드르 함수 \(P_l^m(x)\)\(d P_l^m(x)/ dx\) 값은 각각 result_arrayresult_deriv_array 에 저장됩니다. _e 붇은 함수에서 Condon-Shortly 위상 계수 \((-1)^m\) 의 포함 유무를 조정하려면 csphase\(-1\) 이나 \(1\) 로 설정해 주면 됩니다. 이 계수는 기본적으로 비활성화 되어 있습니다.

int gsl_sf_legendre_deriv_alt_array(const gsl_sf_legendre_t norm, const size_t lmax, const double x, double result_array[], double result_deriv_array[])
int gsl_sf_legendre_deriv_alt_array_e(const gsl_sf_legendre_t norm, const size_t lmax, const double x, const double csphase, double result_array[], double result_deriv_array[])

\(|x| \leq 1\) 값에 대해, \(1\) 차에서 \(lmax\) 까지의 정규화된 버금 르장드르 함수의 값과 대체된 도함수 값을 계산합니다. 르장드르 함수 \(P_l^m(x)\)\(d P_l^m(\cos(\theta))/ d\theta\) 의 값들은 각각 result_arrayresult_deriv_array 에 저장됩니다. _e 붇은 함수에서 Condon-Shortly 위상 계수 \((-1)^m\) 의 포함 유무를 조정하려면 csphase\(-1\) 이나 \(1\) 로 설정해 주면 됩니다. 이 계수는 기본적으로 비활성화 되어 있습니다.

int gsl_sf_legendre_deriv2_array(const gsl_sf_legendre_t norm, const size_t lmax, const double x, double result_array[], double result_deriv_array[], double result_deriv2_array[])
int gsl_sf_legendre_deriv2_array_e(const gsl_sf_legendre_t norm, const size_t lmax, const double x, const double csphase, double result_array[], double result_deriv_array[], double result_deriv2_array[])

\(|x| \leq 1\) 값에 대해, \(1\) 차에서 \(lmax\) 까지의 정규화된 버금 르장드르 함수들, 그 도함수들과 \(2\) 계 도함수 값들을 계산합니다. norm 인자는 어느 정규화 방법을 사용할지 결정합니다. 르장드르 함수 \(P_l^m(x)\)\(d P_l^m(x)/ dx\) , 그리고 \(2\) 계 도함수 \(d^2 P_l^m(x) / dx^2\) 의 값은 각각 result_array , result_deriv_array , 그리고 result_deriv2_array 에 저장됩니다. _e 붇은 함수에서 Condon-Shortly 위상 계수 \((-1)^m\) 의 포함 유무를 조정하려면 csphase\(-1\) 이나 \(1\) 로 설정해 주면 됩니다. 이 계수는 기본적으로 비활성화 되어 있습니다.

int gsl_sf_legendre_deriv2_alt_array(const gsl_sf_legendre_t norm, const size_t lmax, const double x, double result_array[], double result_deriv_array[], double result_deriv2_array[])
int gsl_sf_legendre_deriv2_alt_array_e(const gsl_sf_legendre_t norm, const size_t lmax, const double x, const double csphase, double result_array[], double result_deriv_array[], double result_deriv2_array[])

\(|x| \leq 1\) 값에 대해, \(1\) 차에서 \(lmax\) 까지의 정규화된 버금 르장드르 함수들, 그 대체 도함수들과 \(2\) 계 도함수 값들을 계산합니다. norm 인자는 어느 정규화 방법을 사용할지 결정합니다. 르장드르 함수 \(P_l^m(x)\)\(d P_l^m(\cos(\theta))/ d\theta\) , 그리고 \(2\) 계 도함수 \(d^2 P_l^m(\cos(\theta)) / d\theta^2\) 의 값은 각각 result_array , result_deriv_array, 그리고 result_deriv2_array 에 저장됩니다. _e 붇은 함수에서 Condon-Shortly 위상 계수 \((-1)^m\) 의 포함 유무를 조정하려면 csphase\(-1\) 이나 \(1\) 로 설정해 주면 됩니다. 이 계수는 기본적으로 비활성화 되어 있습니다.

size_t gsl_sf_legendre_nlm(const size_t lmax)

\(lmax\) 까지의 버금 르장드르 함수 \(P_l^m(x)\) 의 갯수를 반환합니다. 해당 값은 \((lmax_1)* (lmax+2)/2\) 입니다.

size_t gsl_sf_legendre_array_n(const size_t lmax)

최대 차수 \(lmax\) 까지의 배열 버전 버금 르장드르 함수에 필요한 최소 배열의 크기를 반환합니다. 이 값은 \(P_l^m(x)\) 의 최댓값과 재귀식을 계산할 때 필요한, 곱 계수 계산을 위한 공간을 더한 값입니다.

size_t gsl_sf_legendre_array_index(const size_t l, const size_t m)

result_array result_deriv_array 그리고 result_deriv2_array 배열의 색인 값을 반환합니다. 해당 값은 \(P_l^m(x)\) , \(P_l^{'m}(x)\) , 그리고 \(P_l^{''m}(x)\) 에 대응되고 주어진 \(l\) \(m\) 대해, \(l(l_1)//2 +m\) 으로 정해집니다.

\(HAVE_INLINE\) 를 사용하면 인라인 버전의 함수를 사용할 수 있습니다.

double gsl_sf_legendre_Plm(int l, int m, double x)
int gsl_sf_legendre_Plm_e(int l, int m, double x, gsl_sf_result *result)

\(m \geq 0, l \geq m\) 그리고 \(|x| \leq 1\) 에 대해, 버금 르장드르 함수 \(P_l^m(x)\) 의 값을 계산합니다.

double gsl_sf_legendre_sphPlm(int l, int m, double x)
int gsl_sf_legendre_sphPlm_e(int l, int m, double x, gsl_sf_result *result)

구면 조화 함수에서 사용하기 위한, 정규회된 버금 르장드르 다항식 \(\sqrt{(2l+1)/(4\pi)} \sqrt{(l-m)! / (l+m)!} P_l^m(x)\) 값을 계산합니다. 앞의 계수는 \(m \geq 0, l \geq m\) 그리고 \(|x| \leq 1\) 를 만족해야합니다. 표준 정규화과정에서 일어나는 오버플로우를 피할 수 있습니다.

int gsl_sf_legendre_Plm_array(int lmax, int m, double x, double result_array[])
int gsl_sf_legendre_Plm_deriv_array(int lmax, int m, double x, double result_array[], double result_deriv_array[])

현재 비활성화 되어 있고 차후 버전에서 삭제될 예정입니다. gsl_sf_legendre_array()gsl_sf_legendre_deriv_array() 를 참고하길 바랍니다.

int gsl_sf_legendre_sphPlm_array(int lmax, int m, double x, double result_array[])
int gsl_sf_legendre_sphPlm_deriv_array(int lmax, int m, double x, double result_array[], double result_deriv_array[])

현재 비활성화 되어 있고 차후 버전에서 삭제될 예정입니다. gsl_sf_legendre_array()gsl_sf_legendre_deriv_array() 를 참고하길 바랍니다.

int gsl_sf_legendre_array_size(const int lmax, const int m)

현재 비활성화 되어 있고 차후 버전에서 삭제될 예정입니다.

원뿔 함수 (Conial Functions)

원통 함수 \(P_{-(1/2)+i\lambda}^\mu\)\(Q^{\mu}_{-(1/2)+i\lambda}\) 는 Abramowitz & Stegun 8.12 단원에 기술되어 있습니다.

double gsl_sf_conicalP_half(double lambda, double x)
int gsl_sf_conicalP_half_e(double lambda, double x, gsl_sf_result *result)

\(x>-1\) 에 대해, 비정칙 구면 원뿔 함수 \(P_{-1/2+i\lambda}^{1/2} (x)\) 값을 계산합니다.

double gsl_sf_conicalP_mhalf(double lambda, double x)
int gsl_sf_conicalP_mhalf_e(double lambda, double x, gsl_sf_result *result)

\(x>-1\) 에 대해, 정칙 구면 원뿔 함수 \(P_{-1/2+i\lambda}^{1/2} (x)\) 값을 계산합니다.

double gsl_sf_conicalP_0(double lambda, double x)
int gsl_sf_conicalP_0_e(double lambda, double x, gsl_sf_result *result)

\(x>-1\) 에 대해, 원뿔 함수 \(P_{-1/2+i\lambda}^{0} (x)\) 값을 계산합니다.

double gsl_sf_conicalP_1(double lambda, double x)
int gsl_sf_conicalP_1_e(double lambda, double x, gsl_sf_result *result)

\(x>-1\) 에 대해, 원뿔 함수 \(P_{-1/2+i\lambda}^{1} (x)\) 값을 계산합니다.

double gsl_sf_conicalP_sph_reg(int l, double lambda, double x)
int gsl_sf_conicalP_sph_reg_e(int l, double lambda, double x, gsl_sf_result *result)

\(x>-1\) , \(l \geq -1\) 에 대해, 정칙 구면 원뿔 함수 \(P_{1/2+i\lambda}^{-1/2-l} (x)\) 값을 계산합니다.

double gsl_sf_conicalP_cyl_reg(int m, double lambda, double x)
int gsl_sf_conicalP_cyl_reg_e(int m, double lambda, double x, gsl_sf_result *result)

\(x>-1\) , \(m \geq -1\) 에 대해, 정칙 원통 원뿔 함수 \(P_{1/2+i\lambda}^{-m} (x)\) 값을 계산합니다.

쌍곡 공간에서의 방사 함수 (Radial Functions for Hyperbolic Space)

다음의 구면 함수들은 \(3\) 차원 쌍곡 공간 \(H^3\) 속 라플라시안의 고유 함수들인 르장드르 함수들입니다. 특히 중점으로 다루는 부분은 평평한 극한(flat limit)으로 \(\lambda \rightarrow \infty, \eta \rightarrow 0\) 이고, \(\lambda \eta\) 가 상수로 고정된 상황입니다.

double gsl_sf_legendre_H3d_0(double lambda, double eta)
int gsl_sf_legendre_H3d_0_e(double lambda, double eta, gsl_sf_result *result)

\(3\) 차원 쌍곡 공간 라플라시안의 \(0\) 차 고유 함수를 계산합니다. 다음과 같이 정의되어 있습니다. \(\eta \geq 0\) 에 대해,

\[L_0^{H 3d} (\lambda, \eta) := \frac{\sin(\lambda \eta)}{\lambda \sinh(\eta)}\]

평평한 극한값은 \(L_0^{H 3d} (\lambda ,\eta) = j_0(\lambda \eta)\) 의 형태를 가집니다.

double gsl_sf_legendre_H3d_1(double lambda, double eta)
int gsl_sf_legendre_H3d_1_e(double lambda, double eta, gsl_sf_result *result)

\(3\) 차원 쌍곡 공간 라플라시안의 \(1\) 차 고유 함수를 계산합니다. 다음과 같이 정의되어 있습니다. \(\eta \geq 0\) 에 대해,

\[L_1^{H 3d} (\lambda, \eta) := \frac{1}{\sqrt{\lambda^2 +1}}(\frac{\sin(\lambda \eta)}{\lambda \sinh(\eta)}) (\coth(\eta) - \lambda \cot(\lambda \eta))\]

평평한 극한값은 \(L_1^{H 3d} (\lambda ,\eta) = j_1(\lambda \eta)\) 의 형태를 가집니다.

double gsl_sf_legendre_H3d(int l, double lambda, double eta)
int gsl_sf_legendre_H3d_e(int l, double lambda, double eta, gsl_sf_result *result)

\(\eta \geq 0, l \geq 0\) 에 대해, \(l\) 차수의 방사 고유 함수값을 계산합니다. 이 고유 함수들은 \(3\) 차원 쌍곡 공간 속 라플라시안의 고유 함수들입니다. 평평한 극한값은 \(L_l^{H 3d} (\lambda ,\eta) = j_l (\lambda \eta)\) 형태를 가집니다.

int gsl_sf_legendre_H3d_array(int lmax, double lambda, double eta, double result_array[])

\(0 \leq l \leq lmax\) 방사 고유 함수 \(L_l^{H 3d} (\lambda, \eta)\) 의 값을 계산합니다.

마티유 함수 (Mathieu Functions)

이 단원에서는 극 마티유 함수와 방사형 마티유 함수를 계산하는 함수들을 서술합니다. 이 함수들의 특성 값의 계산 기능도 같이 제공합니다. 마트유 함수는 다음 두 미분 방정식의 해로 주어집니다.

\[\begin{split}{{d^2 y}\over{d v^2}}& + (a - 2q\cos 2v)y = 0 \\ {{d^2 f}\over{d u^2}}& - (a - 2q\cosh 2u)f = 0\end{split}\]

극 마티유 함수 \(ce_r(x,q)\)\(se_r(x,q)\) 은 각각 첫번째 미분 방정식의 짝수, 홀수 번째 주기 함수해 입니다. 이 첫번째 미분 방정식은 마티유 방정식으로도 알려저 있습니다. 이 해들은 짝수와 홀수에 대해 각각 특성값의 이산 배열값들 \(a = a_r(q)\)\(a = b_r(q)\) 에 대해서만 존재합니다.

방사형 마티유 함수 \(Mc^{(j)}_{r}(z,q)\)\(Ms^{(j)}_{r}(z,q)\) 는 두번째 미분 방정식의 해로 이 미분 방성직은 수정 마티유 방정식으로 알려져 있습니다. 1, 2, 3, 그리고 4차 방사 마티유 함수들은 계수 \(j\) 로 표현됩니다. \(j\) 는 1, 2, 3, 4 값을 가집니다.

마티유 함수에 대한 더 자세한 정보를 얻고 싶다면 Abramowitz and Stegunn, Chapter 20을 참고하길 바랍니다. 이 함수들은 헤더 파일 gsl_sf_mathieu.h 에 정의되어 있습니다.

마티유 함수 작업 공간 (Mathieu Function Workspace)

마티유 함수들은 단일 차수나 복수 차수의 함수가 있으며, 배열에 기반한 방법을 통해 계산할 수 있습니다. 배열 기반 함수들은 사전 정의된 작업 공간이 필요합니다.

type gsl_sf_mathieu_workspace

배열 기반 방법을 위한 작업 공간입니다.

gsl_sf_mathieu_workspace *gsl_sf_mathieu_alloc(size_t n, double qmax)

마티유 함수들의 배열 기반 방법을 위한 작업 공간을 반환합니다. 인자 nqmax 최대 계수와 공간에서 계산되는 마티유 함수의 \(q\) -값을 결정합니다.

void gsl_sf_mathieu_free(gsl_sf_mathieu_workspace *work)

주어진 작업 공간 work 를 해제합니다.

마티유 함수 특성 값 (Mathieu Function Characteristic Values)

int gsl_sf_mathieu_a(int n, double q)
int gsl_sf_mathieu_a_e(int n, double q, gsl_sf_result *result)
int gsl_sf_mathieu_b(int n, double q)
int gsl_sf_mathieu_b_e(int n, double q, gsl_sf_result *result)

각각, 마티유 함수 \(c e_n (q,x)\)\(s e_n (q,x)\) 의 특성 값 \(a_n(q)\)\(b_n(q)\) 을 계산합니다.

int gsl_sf_mathieu_a_array(int order_min, int order_max, double q, gsl_sf_mathieu_workspace *work, double result_array[])
int gsl_sf_mathieu_b_array(int order_min, int order_max, double q, gsl_sf_mathieu_workspace *work, double result_array[])

마티유 함수의 특성 값 \(a_n(q)\)\(b_n(q)\) 을 주어진 oreder_minorder_max 사이 범위에 있는 n 에 대해 계산합니다. 계산 결과는 result_array 에 저장됩니다.

각 운동량 마티유 함수 (Angular Mathieu Functions)

int gsl_sf_mathieu_ce(int n, double q, double x)
int gsl_sf_mathieu_ce_e(int n, double q, double x, gsl_sf_result *result)
int gsl_sf_mathieu_se(int n, double q, double x)
int gsl_sf_mathieu_se_e(int n, double q, double x, gsl_sf_result *result)

각각, 각 운동량 마티유 함수 \(c e_n (q,x)\)\(s e_n (q,x)\) 를 계산합니다.

int gsl_sf_mathieu_ce_array(int nmin, int nmax, double q, double x, gsl_sf_mathieu_workspace *work, double result_array[])
int gsl_sf_mathieu_se_array(int nmin, int nmax, double q, double x, gsl_sf_mathieu_workspace *work, double result_array[])

각각, 각 운동량 마티유 함수 \(c e_n (q,x)\)\(s e_n (q,x)\) 의 값을 주어진 nminnmax 범위에 있는 n 에 대해 계산합니다. 계산 결과는 result_array 에 저장됩니다.

방사 마티유 함수 (Radial Mathieu Functions)

int gsl_sf_mathieu_Mc(int j, int n, double q, double x)
int gsl_sf_mathieu_Mc_e(int j, int n, double q, double x, gsl_sf_result *result)
int gsl_sf_mathieu_Ms(int j, int n, double q, double x)
int gsl_sf_mathieu_Ms_e(int j, int n, double q, double x, gsl_sf_result *result)

각각, \(j\)\(n\) 차수의 마티유 함수 \(M c_n^{(j)} (q,x)\)\(M s_n^{(j)} (q,x)\) 를 계산합니다.

\(j\) 값은 \(1,2\) 로 한정됩니다. \(j=3, 4\) 는 다음의 관계를 이용해 계산할 수 있습니다. \(M_n^{(j)} =\) \(Mc_n^{(j)}\)\(Ms_n^{(j)}\) 에 대해, \(M_n^{(3)} = M_n^{(1)} + i M_n^{(2)}\)\(M_n^{(4)} = M_n^{(1)} - i M_n^{(2)}\) .

int gsl_sf_mathieu_Mc_array(int j, int nmin, int nmax, double q, double x, gsl_sf_mathieu_workspace *work, double result_array[])
int gsl_sf_mathieu_Ms_array(int j, int nmin, int nmax, double q, double x, gsl_sf_mathieu_workspace *work, double result_array[])

\(j\) 메티유 함수의 값을 주어진 nminnmax 범위에 있는 n 에 대해 계산합니다. 계산 결과는 result_array 에 저장됩니다.

멱함수 (Power Function)

다음의 함수들은 오차 계산과 함께라면, 함수 gsl_pow_int() 같습니다. 이 함수들은 헤더 파일 \(gsl_sf_pow_int.h\) 에 정의되어 있습니다.

double gsl_sf_pow_int(double x, int n)
int gsl_sf_pow_int_e(double x, int n, gsl_sf_result *result)

정수 \(n\) 대해 \(x^n\) 의 값을 계산합니다. 이 때, 거듭 제곱은 곱셈 횟수를 최소화해 계산됩니다. 예를 들어서, \(x^8\) 을 계산하고자 하면, \(3\) 번의 곱셈을 이용해 \(((x^2)^2)^2\) 으로 계산합니다. 효율성을 위해서 이 함수는 오버 플로우나 언더 플로우 조건을 검사하지 않습니다.

다음은 간단한 예제입니다.

#include <gsl/gsl_sf_pow_int.h>
/* compute 3.0**12 */
double y = gsl_sf_pow_int(3.0, 12);

프사이(디감마) 함수 (Psi (Digamma) Function)

차수 \(n\) 의 극 감마 함수는 다음과 같이 정의됩니다.

\[\psi^{(n)}(x) = \left(d \over dx\right)^n \psi(x) = \left(d \over dx\right)^{n+1} \log(\Gamma(x))\]

\(\psi(x) = \Gamma'(x)/\Gamma(x)\) 는 Digamma 함수로 알려져 있습니다. 이 단원의 함수들은 헤더 파일 \(gsl_sf_psi.h\) 에 기술되어 있습니다.

디감마 함수 (Digamma Function)

double gsl_sf_psi_int(int n)
int gsl_sf_psi_int_e(int n, gsl_sf_result *result)

디감마 함수 \(\psi(n)\) 값을 주어진 양의 정수 \(n\) 대해 계산합니다. 디감마 함수는 프사이(Psi) 함수라고도 불립니다.

double gsl_sf_psi(double x)
int gsl_sf_psi_e(double x, gsl_sf_result *result)

디감마 함수 \(\psi(x)\) 를 주어진 값 \(x\) 대해 계산합니다. 이때, \(x \neq 0\) 입니다.

double gsl_sf_psi_1piy(double y)
int gsl_sf_psi_1piy_e(double y, gsl_sf_result *result)

\(1+iy\) 선 위에서의 디감마 함수의 실수부 \(\mathfrak{R}[\psi(1+iy)]\) 를 계산합니다.

Trigamma 함수 (Trigamma Function)

double gsl_sf_psi_1_int(int n)
int gsl_sf_psi_1_int_e(int n, gsl_sf_result *result)

Trigamma 함수 \(\psi'(n)\) 을 주어진 양의 정수 \(n\) 대해 계산합니다.

double gsl_sf_psi_1(double x)
int gsl_sf_psi_1_e(double x, gsl_sf_result *result)

Trigamma 함수 \(\psi'(x)\) 를 주어진 값 \(x\) 대해 계산합니다.

Polygamma 함수 (Polygamma Function)

double gsl_sf_psi_n(int n, double x)
int gsl_sf_psi_n_e(int n, double x, gsl_sf_result *result)

\(n \geq 0, x>0\)\(n,x\) 에 대해, polygamma 함수 \(\psi^{(n)}(x)\) 를 계산합니다.

싱크트론 함수 (Synchrotron Functions)

이 단원의 함수들은 gsl_sf_sychrotron.h 해더 파일에 기술되어 있습니다.

double gsl_sf_synchrotron_1(double x)
int gsl_sf_synchrotron_1_e(double x, gsl_sf_result *result)

\(x \geq 0\) 에 대해, 1차 싱크로트론 함수 \(x \int_x^\infty K_{5/3} \, dt\) 를 계산합니다.

double gsl_sf_synchrotron_2(double x)
int gsl_sf_synchrotron_2_e(double x, gsl_sf_result *result)

\(x \geq 0\) 에 대해, 2차 싱크로트론 함수 \(x \int_x^\infty K_{2/3} \, dt\) 를 계산합니다.

운송 함수 (Transport Functions)

운송 함수 \(J(n,x)\) 는 다음과 같이 적분으로 정의됩니다.

\[ \begin{align}\begin{aligned}J(n,x) = \int_0^x \frac{t^ne^t}{(e^t -1)^2} \, dx\\헤더 파일 :code:`gsl_sf_transport.h` 에 정의되어 있습니다.\end{aligned}\end{align} \]
double gsl_sf_transport_2(double x)
int gsl_sf_transport_2_e(double x, gsl_sf_result *result)

운송 함수 \(J(2,x)\) 의 값을 계산합니다.

double gsl_sf_transport_3(double x)
int gsl_sf_transport_3_e(double x, gsl_sf_result *result)

운송 함수 \(J(3,x)\) 의 값을 계산합니다.

double gsl_sf_transport_4(double x)
int gsl_sf_transport_4_e(double x, gsl_sf_result *result)

운송 함수 \(J(4,x)\) 의 값을 계산합니다.

double gsl_sf_transport_5(double x)
int gsl_sf_transport_5_e(double x, gsl_sf_result *result)

운송 함수 \(J(5,x)\) 의 값을 계산합니다.

삼각 함수 (Trigonometric Functions)

본 라이브러리는 여러 플랫폼에서 일관적인 코드 작성과 신뢰할 수 있는 오류 측정을 위해 자체적인 삼각 함수를 포함하고 있습니다.

헤더파일 gsl_sf_trig.h 에 정의되어 있습니다.

삼각 함수 (Circular Trigonometric Functions)

double gsl_sf_sin(double x)
int gsl_sf_sin_e(double x, gsl_sf_result *result)

sine 함수 \(\sin(x)\) 를 계산합니다.

double gsl_sf_cos(double x)
int gsl_sf_cos_e(double x, gsl_sf_result *result)

cosine 함수 \(\cos(x)\) 를 계산합니다.

double gsl_sf_hypot(double x, double y)
int gsl_sf_hypot_e(double x, double y, gsl_sf_result *result)

오버 플로우와 언더 플로우 없는 \(\sqrt{x^2+y^2}\) 값을 계산합니다.

double gsl_sf_sinc(double x)
int gsl_sf_sinc_e(double x, gsl_sf_result *result)

\(x\) 에 대해, \(\text{sinc}(x) = \sin(\pi x)/(\pi x)\) 값을 계산합니다.

복소 삼각 함수 (Trigonometric Functions for Complex Arguments)

int gsl_sf_complex_sin_e(double zr, double zi, gsl_sf_result *szr, gsl_sf_result *szi)

복소 sine 함수 \(\sin(z_r + i z_i)\) 값을 계산합니다. 계산 값의 실수부는 \(szr\) 에 허수부는 \(szi\) 에 저장됩니다.

int gsl_sf_complex_cos_e(double zr, double zi, gsl_sf_result *czr, gsl_sf_result *czi)

복소 cosine 함수 \(\cos(z_r + i z_i)\) 값을 계산합니다. 계산 값의 실수부는 \(czr\) 에, 허수부는 \(czi\) 에 저장됩니다.

int gsl_sf_complex_logsin_e(double zr, double zi, gsl_sf_result *lszr, gsl_sf_result *lszi)

sine 함수의 로그 값 \(\log(\sin(z_r + i z_i))\) 값을 계산합니다. 계산 값의 실수부는 \(lszr\) 에 허수부는 \(lszi\) 에 저장됩니다.

쌍곡 함수 (Hyperbolic Trigonometric Functions)

double gsl_sf_lnsinh(double x)
int gsl_sf_lnsinh_e(double x, gsl_sf_result *result)

\(\log(\text{sinh}(x))\) 값을, \(x>0\) 에 대해 계산합니다.

double gsl_sf_lncosh(double x)
int gsl_sf_lncosh_e(double x, gsl_sf_result *result)

주어진 값 \(x\) 대해, \(\log(\text{cosh}(x))\) 값을 계산합니다.

좌표 변환 함수 (Conversion Functions)

int gsl_sf_polar_to_rect(double r, double theta, gsl_sf_result *x, gsl_sf_result *y)

극 좌표( \(r\) \(theta\) ) 를 직교 좌표( \(x\) \(y\) )로 변환합니다. \(x = r \cos(\theta), y= r\sin(\theta)\) 관계를 이용합니다.

int gsl_sf_rect_to_polar(double x, double y, gsl_sf_result *r, gsl_sf_result *theta)

직교 좌표( \(x\) \(y\) )를 극 좌표( \(r\) \(theta\) )로 변환합니다. \(x = r \cos(\theta), y= r\sin(\theta)\) 관계를 이용합니다. \(theta\)\([-\pi, \pi]\) 의 범위를 가집니다.

각 제한 함수 (Restriction Functions)

double gsl_sf_angle_restrict_symm(double theta)
int gsl_sf_angle_restrict_symm_e(double *theta)

\(theta\) 값을 \((-\pi, pi]\) 범위 내에 있도록 변환합니다.

참고

실제 \(\pi\) 값은 \(M_PI\) 보다 조금 큽니다.따라서, \(M_PI\) \(-M_PI\) 이 범위에 포함되어 있습니다.

double gsl_sf_angle_restrict_pos(double theta)
int gsl_sf_angle_restrict_pos_e(double *theta)

\(theta\) 값을 \((0, 2pi]\) 범위 내에 있도록 변환합니다.

참고

실제 \(2\pi\) 값은 \(2M_PI\) 보다 조금 큽니다. 따라서, \(2*M_PI\) 이 범위에 포함되어 있습니다.

오차 분석을 포함한 삼각 함수 (Trigonometric Functions With Error Estimates)

int gsl_sf_sin_err_e(double x, double dx, gsl_sf_result *result)

\(x\) 대해, 버금 절대 오차 \(dx\) 가 포함된 sine 값 \(\sin(x \pm dx)\) 을 계산합니다.

참고

오차 전파를 계산하기 위해 제공하는 것이기 때문에 오차 관리 함수 형태만으로 제공됩니다.

int gsl_sf_cos_err_e(double x, double dx, gsl_sf_result *result)

\(x\) 대해, 버금 절대 오차 \(dx\) 포함된 cosine 값 \(\cos(x \pm dx)\) 을 계산합니다.

참고

오차 전파를 계산하기 위해 제공하는 것이기 때문에 오차 관리 함수 형태만으로 제공됩니다.

제타 함수 (Zeta Functions)

리만 제타 함수 (Riemann Zeta Function)

리만 제타 함수는 다음과 같이 무한 급수로 정의됩니다.

\[\zeta (s) = \sum_{k=1}^\infty k^{-s}\]
double gsl_sf_zeta_int(int n)
int gsl_sf_zeta_int_e(int n, gsl_sf_result *result)

정수 \(n\) 대한 리만 제타 함수 \(\zeta(n)\) 의 값을 계산합니다. \(n \neq 1\) 이어야 합니다.

double gsl_sf_zeta(double s)
int gsl_sf_zeta_e(double s, gsl_sf_result *result)

임의의 수 \(s\) 대해 리만 제타 함수 \(\zeta(s)\) 값을 계산합니다. 위 함수와 마찬가지로 \(s \neq 1\) 이어야 합니다.

리만 제타 함수 -1 (Riemann Zeta Function Minus One)

큰 양수 값에 대해, 리만 제타함수는 1로 수렴하게 됩니다. 이 경우에 1이 아닌 부분의 값이 빈번히 사용되므로, 라이브러리에서는 이를 위한 함수를 제공합니다.

double gsl_sf_zetam1_int(int n)
int gsl_sf_zetam1_int_e(int n, gsl_sf_result *result)

정수 \(n\) 대한 \(\zeta(n)-1\) 의 값을 계산합니다. \(n \neq 1\) 이어야 합니다.

double gsl_sf_zetam1(double s)
int gsl_sf_zetam1_e(double s, gsl_sf_result *result)

임의의 수 \(s\) 대해 \(\zeta(s)-1\) 값을 계산합니다. 위 함수와 마찬가지로 \(s \neq 1\) 이어야 합니다.

후르비츠(Hurwitz) 제타 함수 (Hurwitz Zeta Function)

후르비츠 제타함수는 다음과 같이 정의됩니다.

\[\zeta(s,q) = \sum_0^\infty (k+q)^{-s}\]
double gsl_sf_hzeta(double s, double q)
int gsl_sf_hzeta_e(double s, double q, gsl_sf_result *result)

후르비츠 제타 함수 \(\zeta (s,q)\) 의 값을 계산합니다. \(s>1, q>0\) 이어야 합니다.

에타(Eta) 함수 (Eta Function)

에타 함수는 다음과 같이 정의됩니다.

\[\eta(s) = (1-2^{1-s}) \zeta(s)\]
double gsl_sf_eta_int(int n)
int gsl_sf_eta_int_e(int n, gsl_sf_result *result)

정수 \(n\) 대해 에타 함수 \(\eta(n)\) 의 값을 계산합니다.

double gsl_sf_eta(double s)
int gsl_sf_eta_e(double s, gsl_sf_result *result)

임의의 수 \(s\) 대해 에타 함수 \(\eta(s)\) 의 값을 계산합니다.

예제 (Examples)

다음 예제는 오차 관리 형태의 특수 함수 구현체를 사용하는 방법을 보여줍니다. 베셀 함수 \(J_0 (5.0)\) 값을 계산합니다.

#include <stdio.h>
#include <gsl/gsl_errno.h>
#include <gsl/gsl_sf_bessel.h>

int
main (void)
{
  double *x* = 5.0;
  gsl_sf_result result;

  double expected = -0.17759677131433830434739701;

  int status = gsl_sf_bessel_J0_e (x, &result);

  printf ("status  = %s\n", gsl_strerror(status));
  printf ("J0(5.0) = %.18f\n"
          "      +/- % .18f\n",
          result.val, result.err);
  printf ("exact   = %.18f\n", expected);
  return status;
}

다음은 프로그램의 출력 결과입니다.

status  = success
J0(5.0) = -0.    177596771314338264
      +/-  0.    00000000000000019   3
exact   = -0.    177596771314338292

다음 프로그램은 같은 함수(베셀 함수)로 똑같은 값(\(J_0 (5.0)\))을 계산합니다. 이 때, \(result.err\) 반환 상태는 존재하지 않습니다.

#include <stdio.h>
#include <gsl/gsl_sf_bessel.h>

int
main (void)
{
  double *x* = 5.0;
  double expected = -0.  177596771314338304347     39701;

  double y =     gsl_sf_bessel_J0 (x);

  printf ("J0(5.0) = %.  18f\n", y);
  printf ("exact   = %.  18f\n", expected);
  return 0;
}

출력 결과는 다음과 같습니다.

J0(5.0) = -0.177596771314338264
exact   = -0.177596771314338292

참고 문헌과 추가 자료

이 라이브러리는 다음 책의 규약을 따릅니다.

  • Handbook of Mathematical Functions, edited by Abramowitz & Stegun, Dover, ISBN 0486612724.

다음 논문들은 특수 함수들을 계산하기 위한 알고리즘에 관한 내용입니다.

  • Allan J. MacLeod, MISCFUN: A software package to compute uncommon special functions. ACM Trans. Math. Soft., vol.: 22, 1996, 288–301

  • Bunck, B. F., A fast algorithm for evaluation of normalized Hermite functions, BIT Numer. Math, 49: 281-295, 2009.

  • G.N. Watson, A Treatise on the Theory of Bessel Functions, 2nd Edition (Cambridge University Press, 1944).

  • G. Nemeth, Mathematical Approximations of Special Functions, Nova Science Publishers, ISBN 1-56072-052-2

  • B.C. Carlson, Special Functions of Applied Mathematics (1977)

  • N. M. Temme, Special Functions: An Introduction to the Classical Functions of Mathematical Physics (1996), ISBN 978-0471113133.

  • W.J. Thompson, Atlas for Computing Mathematical Functions, John Wiley & Sons, New York (1997).

  • Y.Y. Luke, Algorithms for the Computation of Mathematical Functions, Academic Press, New York (1977).

  • S. A. Holmes and W. E. Featherstone, A unified approach to the Clenshaw summation and the recursive computation of very high degree and order normalised associated Legendre functions, Journal of Geodesy, 76, pg. 279-299, 2002.