[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7 Funções Especiais

Esse capítulo descreve a biblioteca de funções especias da GSL. A biblioteca de funções especias inclui rotinas para o cálculo dos valores das funções de Airy, das funções de Bessel, das funções de Clausen, das funções de onda de Coulomb, dos coeficientes de Coupling, da função de Dawson, das funções de Debye, dos Dilogaritmos, de integrais Elípticas, das funções elípticas de Jacobi, das funções de Erro, das integrais de exponenciais, das funções de Fermi-Dirac, das funções Gama, das funções de Gegenbauer, das funções Hipergeométricas, das funções de Laguerre, das funções de Legendre e Harmônicas esféricas, da Função Psi (Digama), Das funções do sincrotron, das funções de Transporte, das funções trigonométricas e das Funções Zeta. Cada rotina também calcula uma estimativa do erro numérico no valor calculado da função.

As funções nesse capítulo são declaradas em arquivos de cabeçalho individuais, tais como ‘gsl_sf_airy.h’, ‘gsl_sf_bessel.h’, etc. O conjunto completo de cabeçalhos pode ser incluído usando o arquivo ‘gsl_sf.h’.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.1 Utilização

As funções especias estão disponíveis em duas convenções de chamada, uma forma natural que retorna o valor numérico da função e uma forma controladora de erro que retorna um código de erro. Os dois tipos de função fornecem caminhos alternativos de acesso o mesmo código básico.

A forma natural retorna somente o valor da função e pode ser usado diretamente em expressões matemáticas. Por exemplo, A seguinte chamada de função irá calcular o vaor da função de Bessel J_0(x),

double y = gsl_sf_bessel_J0 (x);

Não existe forma de acessar um código de erro ou de estimar o erro usando esse método. Para permitir o acesso a essa informação a forma alternativa controladora de erro armazena o valor da função e o erro em um argumento modificável,

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

As funções controladoras de erro possuem o sufixo _e. O valor da situação atual indica condições de erro tais como overflow, underflow ou perda de precisão. Caso não hajam erros as funções controladoras de erro retornam GSL_SUCCESS.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.2 A estrutura gsl_sf_result

A forma controladora de erro das funções especiais sempre calculam um erro estimado juntamente com o valor do resultado. Portanto, estruturas são fornecidas para amalgamar um valor e um erro estimado. Essas estruturas são declaradas no arquivo de cabeçalho ‘gsl_sf_result.h’.

A estrutura gsl_sf_result contém campos valor e erro.

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

O campo val contém o valor e o campo err contém uma estimativa de erro absoluto para o valor.

Em alguns casos, um overflow ou underflow pode ser detectado e manipulado pela função controladora de erro. Nesse caso, pode ser possível retornar um expoente dimensionado bem como um par erro/valor com o objetivo de evitar que o resultado exceda o intervalo dinâmico dos tipos da dado construídos dentro da função. A estrutura gsl_sf_result_e10 contém campos valor e erro bem como um campo expoente tal que o resultado atual é obtido como result * 10^(e10).

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

[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.3 Modos

O objetivo da biblioteca é encontrar grande exatidão na precisão dupla sempre que possível. Todavia o custo de avaliação de algumas funções especiais em precisão dupla pode ser significante, particularmente onde são requeridos termos de ordem muito alta. Nesses casos um argumento mode permite que a exatidão do termo seja reduzida com o objetivo de melhorar o desempenho. Os seguintes níveis de exatidão estão disponíveis para o argumento mode,

GSL_PREC_DOUBLE

Precisão dupla, uma exatidão relativa de aproximadamente 2 * 10^-16.

GSL_PREC_SINGLE

Precisão simples, uma exatidão relativa de aproximadamente 10^-7.

GSL_PREC_APPROX

Valores aproximados, uma exatidão relativa de aproximadamente 5 * 10^-4.

O modo de valores aproximados fornecem a avaliação mais rápida com a exatidão mais baixa.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.4 Funções de Airy e Derivadas

As funções de Airy Ai(x) e Bi(x) são definidas através das representações de integral,

Ai(x)
= 1

π



0 
cos(t3/3 + xt )  dt,
Bi(x)
= 1

π



0 
(e−t3/3 + xt + sin(t3/3 + xt))  dt.
Para informações adicionais veja Abramowitz & Stegun, Seção 10.4. As funções de Airy são definidas no arquivo de cabeçalho ‘gsl_sf_airy.h’.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.4.1 Funções Airy

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

Essas rotinas calculam a função de Airy Ai(x) com uma exatidão especificada através de mode.

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

Essas rotinas calculam a função de Airy Bi(x) com uma exatidão especificada através de mode.

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

Essas rotinas calculam uma versão ajustada proporcionalmente da função de Airy S_A(x) Ai(x). Para x>0 o fator de proporcionalidade S_A(x) é \exp(+(2/3) x^(3/2)), e é 1 para x<0.

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

Essas rotinas calculam uma versão alterada proporcionalmente da função de Airy S_B(x) Bi(x). Para x>0 o fator de proporcionalidade S_B(x) é exp(-(2/3) x^(3/2)), e é 1 para x<0.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.4.2 Derivadas das Funções de Airy

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

Essas rotinas calculam a derivada da função de Airy Ai'(x) com uma exatidão especificada através de mode.

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

Essas rotinas calculam a derivada da função de Airy Bi'(x) com uma exatidão especificada através de mode.

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

Essas rotinas calculam a derivada da função de Airy alterada proporcionalmente S_A(x) Ai'(x). Para x>0 o fator de proporcionalidade S_A(x) é \exp(+(2/3) x^(3/2)), e é 1 para x<0.

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

Essas rotinas calculam a derivada da função de Airy alterada proporcionalmente S_B(x) Bi'(x). Para x>0 o fator de proporcionalidade S_B(x) é exp(-(2/3) x^(3/2)), e é 1 para x<0.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.4.3 Zeros das Funções de Airy

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

Essas rotinas calculam a localização da s-ésima raíz da função de Airy Ai(x).

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

Essas rotinas calculam a localização da s-ésima raíz da função de Airy Bi(x).


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.4.4 Zeros de Derivadas de Funções de Airy

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

Essas rotinas calculam a localização da s-ésima raíz da derivada da função de Airy Ai'(x).

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

Essas rotinas calculam a localização da s-ésima raíz da derivada da função de Airy Bi'(x).


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.5 Funções de Bessel

As rotinas descritas nessa seção calculam as funções de Bessel Cilíndricas J_n(x), Y_n(x), funções de Bessel Cilíndricas Modificadas I_n(x), K_n(x), funções de Bessel Esféricas j_l(x), y_l(x), e funções de Bessel Esféricas Modificadas i_l(x), k_l(x). Para mais informação veja Abramowitz & Stegun, Capítulos 9 e 10. As funções de Bessel são definidas no arquivo de cabeçalho ‘gsl_sf_bessel.h’.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.5.1 Funções de Bessel Cilíndricas Regulares

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

Essas rotinas calculam a função de Bessel cilíndrica regular de zero-ésima ordem, J_0(x).

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

Essas rotinas calculam a função de Bessel cilíndrica regular primeira ordem, J_1(x).

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

Essas rotinas calculam a função de Bessel cilíndrica regular de ordem n, J_n(x).

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

Essa rotina calcula os valores das funções de Bessel cilíndricas regulares J_n(x) para n variando de nmin a nmax inclusive, armazenando os resultados no vetor estático result_array. Os valores são calculados usando relações de recorrência por eficiência, e portanto podem diferir levemente dos valores exatos.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.5.2 Funções de Bessel Cilíndricas Irregulares

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

Essas rotinas calculam a função de Bessel cilíndrica irregular da zero-ésima ordem, Y_0(x), for x>0.

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

Essas rotinas calculam a função de Bessel cilíndrica irregular da primeira ordem, Y_1(x), for x>0.

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

Essas rotinas calculam a função de Bessel cilíndrica irregular de ordem n, Y_n(x), para x>0.

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

Essa rotina calcula os valores das funções de Bessel cilíndricas irregulares Y_n(x) para n variando de nmin a nmax inclusive, armazenando os resultados no vetor estático result_array. O domínio da função é x>0. Os valores são calculados usando relações de recorrência por eficiência, e portanto podem diferir levemente dos valores exatos.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.5.3 Funções de Bessel Cilíndricas Modificadas Regulares

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

Essas rotinas calculam a função de Bessel cilíndrica modificada regular de zero-ésima ordem, I_0(x).

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

Essas rotinas calculam a função de Bessel cilíndrica modificada regular de primeira ordem, I_1(x).

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

Essas rotinas calculam a função de Bessel cilíndrica modificada regular de ordem n, I_n(x).

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

Essa rotina calcula os valores das funções de Bessel cilíndricas modificadas regulares I_n(x) para n variando de nmin a nmax inclusive, armazenando os resultados no vetor estático result_array. O início do intervalo nmin deve ser positivo ou zero. Os valores são calculados usando relações de recorência por eficiência, e portanto podem diferir levemente dos valores exatos.

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

Essas rotinas calculam a função de Bessel cilíndrica modificada regular de forma proporcional da zero-ésima ordem \exp(-|x|) I_0(x).

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

Essas rotinas calculam a função de Bessel cilíndrica modificada regular de forma proporcional de primeira ordem \exp(-|x|) I_1(x).

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

Essas rotinas calculam a função de Bessel cilíndrica modificada regular de forma proporcional de ordem n, \exp(-|x|) I_n(x)

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

essa rotina calcula os valores das funções de Bessel cilíndricas regulares alteradas proporcionalmente \exp(-|x|) I_n(x) para n variando de nmin a nmax inclusive, armazenando os resultados no vetor estático result_array. O início do intervalo nmin deve ser positivo ou zero. Os valores são calculados usando relações de recorrência por eficiência, e portanto podem diferir levemente dos valores exatos.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.5.4 Funções de Bessel Cilíndricas Modificadas Irregulares

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

Essas rotinas calculam a função de Bessel cilíndrica modificada irregular de zero-ésima ordem, K_0(x), for x > 0.

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

Essas rotinas calculam a função de Bessel cilídrica modificada irregular de primeira ordem, K_1(x), para x > 0.

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

Essas rotinas calculam a função de Bessel cilíndrica modificada irregular de ordem n, K_n(x), para x > 0.

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

Essa rotina calcula os valores das funções de Bessel cilíndricas modificadas irregulares K_n(x) para n variando de nmin a nmax inclusive, armazenando os resultados no vetor estático result_array. O início do intervalo nmin deve ser positivo ou zero. O domínio da função é x>0. Os valores são calculados usando relações de recorrência por eficiência, e portanto podem diferir levemente dos valores exatos.

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

Essas rotinas calculam a função de Bessel cilíndrica modificada irregular ajustada proporcionalmente de zero-ésima ordem \exp(x) K_0(x) para x>0.

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

Essas rotinas calculam a função de Bessel cilíndrica modificada irregular de primeira ordem \exp(x) K_1(x) para x>0.

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

Essas rotinas calculam a função de Bessel cilíndrica modificada irregular ajustada proporcionalmente de ordem n, \exp(x) K_n(x), para x>0.

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

Essa rotina calcula os valores das funções de Bessel cilíndricas irregulares ajustadas proporcionalmente \exp(x) K_n(x) para n variando de nmin a nmax inclusive, armazenando os resultados no vetor estático result_array. O início do intervalo nmin deve ser positivo ou zero. O domínio da função é x>0. Os valores são calculados usando relações de recorrência por eficiência, e portanto podem diferir levemente dos valores exatos.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.5.5 Funções de Bessel Esférica Regular

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

Essas rotinas calculam a função de Bessel esférica regular de zero-ésima ordem, j_0(x) = \sin(x)/x.

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

Essas rotinas calculam a função de Bessel esférica de primeira ordem, j_1(x) = (\sin(x)/x - \cos(x))/x.

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

Essas rotinas calculam a função de Bessel esférica regular de segunda ordem, j_2(x) = ((3/x^2 - 1)\sin(x) - 3\cos(x)/x)/x.

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

Essas rotinas calculam a função de Bessel esférica regular de ordem l, j_l(x), para l >= 0 e x >= 0.

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

Essa rotina calcula os valores das funções de Bessel esféricas regulares j_l(x) para l variando de 0 a lmax inclusive para lmax >= 0 e x >= 0, armazenando os resultados no vetor estático result_array. Os valores são calculados usando relações de recoreência por eficiência e portanto podem diferir levemente dos valores exatos.

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

Essa rotina usa o método de Steed para calcular os valores das funções de Bessel esféricas regulares j_l(x) para l variando de 0 a lmax inclusive para lmax >= 0 e x >= 0, armazenando os resultados no vetor estático result_array. O algoritmo de Steed/Barnett está descrito em Comp. Phys. Comm. 21, 297 (1981). O método de Steed é mais estável que a recorrência usada em outras funções mas é também mais lento.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.5.6 Funções de Bessel Esféricas Regulares

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

Essa rotina calcula a função de Bessel esférica irregula de zero-ésima ordem, y_0(x) = -\cos(x)/x.

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

Essa rotina calcula a função de Bessel esférica irregular de primeira ordem, y_1(x) = -(\cos(x)/x + \sin(x))/x.

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

Essa rotina calcula a função de Bessel esférica irregular de segunda ordem, y_2(x) = (-3/x^3 + 1/x)\cos(x) - (3/x^2)\sin(x).

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

Essa rotina calcula a função de Bessel esférica irregular de ordem l, y_l(x), para l >= 0.

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

Essa rotina calcula os valores das funções de Bessel esféricas irregulares y_l(x) para l variado de 0 a lmax inclusive para lmax >= 0, armazenando os resultados no vetor estático result_array. Os valores são calculados usando relações de recorrência por eficiência, e portanto podem diferir levemente dos valores exatos.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.5.7 Funções de Bessel Esféricas Modificadas Regulares

As funções de Bessel esféricas modificadas regulares i_l(x) estão relacionadas às funções de Bessel modificadas de ordem fracionária, i_l(x) = \sqrt{\pi/(2x)} I_{l+1/2}(x)

Function: double gsl_sf_bessel_i0_scaled (double x)
Function: int gsl_sf_bessel_i0_scaled_e (double x, gsl_sf_result * result)

Essa rotina calcula a função de Bessel esférica modificada ajustada proporcionalmente de zero-ésima ordem, \exp(-|x|) i_0(x).

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

Essa rotina calcula a função de Bessel esférica modificada ajustada proporcionalmente de primeira ordem, \exp(-|x|) i_1(x).

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

Essa rotina calcula a função de Bessel esférica modificada ajustada proporcionalmente de segunda ordem, \exp(-|x|) i_2(x)

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

Essa rotina calcula função de Bessel esférica modificada ajustada proporcionalmente de ordem l, \exp(-|x|) i_l(x)

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

Essa rotina calcula os valores das funções de Bessel cilíndricas modificadas regulares ajustadas proporcionalmente \exp(-|x|) i_l(x) para l variando de 0 a lmax inclusive para lmax >= 0, armazenando os resultados no vetor estático result_array. Os valores são calculados usando relações de recorrência por eficiência, e portanto podem diferir levemente do valor exato.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.5.8 Funções de Bessel esféricas Modificadas Irregulares

As funções de Bessel esféricas modificadas irregulares k_l(x) estão relacionadas às funções de Bessel modificadas irregulares de ordem fracionária, k_l(x) = \sqrt{\pi/(2x)} K_{l+1/2}(x).

Function: double gsl_sf_bessel_k0_scaled (double x)
Function: int gsl_sf_bessel_k0_scaled_e (double x, gsl_sf_result * result)

Essa rotina calcula a função de Bessel esférica modificada irregular ajustada proporcionalmente de zero-ésima ordem, \exp(x) k_0(x), for x>0.

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

Essa rotina calcula a função de Bessel esférica modificada irregular ajustada proporcionalmente de primeira ordem, \exp(x) k_1(x), para x>0.

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

Essa rotina calcula a função de Bessel esférica modificada irregular ajustada proporcionalmente de segunda ordem, \exp(x) k_2(x), for x>0.

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

Essa rotina calcula a função de Bessel esférica modificada irregular ajustada proporcionalmente de ordem l, \exp(x) k_l(x), para x>0.

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

Essa rotina calcula os valores das funções de Bessel esféricas modificadas irregulares ajustadas proporcionalmente \exp(x) k_l(x) para l variando de 0 a lmax inclusive para lmax >= 0 e x>0, armazenando os resultados no vetor estático result_array. Os valores são calculados usando recorrência por eficiência, e portanto podem diferir levemente dos valores exatos.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.5.9 Regular Bessel Function—Fractional Order

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

Essas rotinas calculam a função de Bessel cilíndrica regular de ordem fracionária \nu, J_\nu(x).

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

Essa função calcula a função de Bessel cilíndrica regular de ordem fracionária \nu, J_\nu(x), avaliada em uma série de valores de x. O vetor estático v de comprimento size contém os valores de x. Eles são assumidos para serem estritamente ordenados e positivos. O vetor estático é sobregravado com o valores de J_\nu(x_i).


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.5.10 Irregular Funções de Bessel—Fractional Order

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

Essa rotina calcula a função de Bessel cilíndrica irregular de ordem fracionária \nu, Y_\nu(x).


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.5.11 Regular Modified Funções de Bessel—Ordem Fracionária

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

Essa rotina calcula a função de Bessel modificada regular de ordem fracionária \nu, I_\nu(x) for x>0, \nu>0.

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

Essa rotina calcula a função de Bessel modificada regular ajustada proporcionalmente de ordem fracionária \nu, \exp(-|x|)I_\nu(x) for x>0, \nu>0.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.5.12 Funções de Bessel Modificadas Irregulares—Ordem Fracionária

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

Essa rotina calcula a função de Bessel modificada irregular de ordem fracionária \nu, K_\nu(x) for x>0, \nu>0.

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

Essa rotina calcula o logaritmo da função de Bessel modificada irregular de ordem fracionária \nu, \ln(K_\nu(x)) for x>0, \nu>0.

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

Essa rotina calcula a função de Bessel modificada irregular ajustada proporcionalmente de ordem fracionária \nu, \exp(+|x|) K_\nu(x) for x>0, \nu>0.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.5.13 Zeros de Funções de Bessel Regulares

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

Essa rotina calcula a localização da s-ésima raíz positiva da função de Bessel J_0(x).

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

Essa rotina calcula a localização da s-ésima raíz positiva da função de Bessel J_1(x).

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

Essa rotina calcula a localização da s-ésima raíz positiva da função de Bessel J_\nu(x). A implementação atual não suporta valores negativos de nu.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.6 Funções de Clausen

A função de Clausen é definida pela seguinte integral,

Cl2(x) = −
x

0 
dt log(2 sin(t/2))
A função de Clausem está relacionada a dilogaritmo por Cl_2(\theta) = \Im Li_2(\exp(i\theta)). As funções de Clausen são declaradas no arquivo d cabeçalho ‘gsl_sf_clausen.h’.

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

Essas rotinas calculam a integral de Clausen Cl_2(x).


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.7 Funções de Coulomb

Os protótipos das funções de Coulomb são declarados no arquivo de cabeçalho ‘gsl_sf_coulomb.h’. Ambos estado de limite e soluções de dispersão estão disponíveis.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.7.1 Estado de Limite de Hidrogênio Normalizado

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

Essas rotinas calculam a função de onda radial do estado de limite do hidrogênio normalizado de menor ordem R_1 := 2Z \sqrt{Z} \exp(-Z r).

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

Essas rotinas calculam a n-ésima função de onda radial do estado de limite do hidrogênio normalizado.

Rn : = 2 Z3/2

n2

2Z r

n

l

 
  ⎛


(n−l−1)!

(n+l)!
 
exp(−Z r/n) L2l+1n−l−1(2Z r / n).
onde L^a_b(x) é o polinômio generalizado de Laguerre (veja seção Funções de Laguerre). A normalização é escolhida de forma que a função de onda \psi seja fornecida por \psi(n,l,r) = R_n Y_{lm}.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.7.2 Funções de Onda de Coulomb

As funções de onda de Coulomb F_L(\eta,x), G_L(\eta,x) são descritas em Abramowitz & Stegun, Capítulo 14. Devido ao fato de poder existir um grande intervalo dinâmico de valores para essas funções, overflows são tratados normalmente. Se um overflow ocorre, GSL_EOVRFLW é sinalizado e expoente(s) é(são) retornado(s) apesar dos parâmetros modificáveis exp_F, exp_G. A solução completa pode ser reconstruída a partir das seguintes relações,

FL(η,x)
= fc[kL] * exp(expF)
GL(η,x)
= gc[kL] * exp(expG)
FL′(η,x)
= fcp[kL] * exp(expF)
GL′(η,x)
= gcp[kL] * exp(expG)

Function: 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)

Essa função calcula as funções de onda de Coulomb F_L(\eta,x), G_{L-k}(\eta,x) e suas derivadas F'_L(\eta,x), G'_{L-k}(\eta,x) em relação a x. Os parâmetros estão restritos a L, L-k > -1/2, x > 0 e ao inteiro k. Note que L em si mesmo não está restrito aos inteiros. Os resultados são armazenados nos parâmetros F, G para os valores da função e Fp, Gp para os valores das derivadas. Se um overflow ocorrer, GSL_EOVRFLW é retornado e expoentes ajustados proporcionalmente são armazenados nos parâmetros modificáveis exp_F, exp_G.

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

Essa função calcula a função de onda de Coulomb F_L(\eta,x) para L = Lmin \dots Lmin + kmax, armazenando os resultados em fc_array. No caso de overflow o expoente é armazenado em F_exponent.

Function: 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)

Essa função calcula as funções F_L(\eta,x), G_L(\eta,x) para L = Lmin \dots Lmin + kmax armazenando os resultados em fc_array e gc_array. No caso de overflow os expoentes são armazenados em F_exponent e G_exponent.

Function: 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)

Essa função calcula as funções F_L(\eta,x), G_L(\eta,x) e suas derivadas F'_L(\eta,x), G'_L(\eta,x) para L = Lmin \dots Lmin + kmax armazenando os resultados em fc_array, gc_array, fcp_array e gcp_array. No caso de overflow os expoentes são armazenados em F_exponent e G_exponent.

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

Essa função calcula a função de onda de Coulomb dividida pelo argumento F_L(\eta, x)/x para L = Lmin \dots Lmin + kmax, armazenando os resultados em fc_array. No caso de overflow o expoente é armazenado em F_exponent. Essa função se reduz a funções de Bessel no limite \eta \to 0.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.7.3 Constante de Normalização da Função de Onda de Coulomb

A constante de normalização da função de onda de Coulomb é definida em Abramowitz 14.1.7.

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

Essa função calcula a constante de normalização da função de onda de Coulomb C_L(\eta) for L > -1.

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

Essa função calcula a constante de normalização da função de onda de Coulomb C_L(\eta) for L = Lmin \dots Lmin + kmax, Lmin > -1.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.8 Acoplando Coeficientes

O símbolos de Wigner 3-j, 6-j e 9-j fornecem os coeficiente de acoplamento para vetores de momentum angular combinado. Uma vez que os argumentos das funções de coeficiente de acoplamento padronizadas são inteiras ou meio inteiras, os argumentos das seguintes funções são, por convenção, inteiros iguais a duas vezes os valor atual do giro. Para informação sobre coeficientes 3-j veja Abramowitz & Stegun, Seção 27.9. A função descrita nessa seção está declarada no arquivo de cabeçalho ‘gsl_sf_coupling.h’.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.8.1 Símbolos 3-j

Function: double gsl_sf_coupling_3j (int two_ja, int two_jb, int two_jc, int two_ma, int two_mb, int two_mc)
Function: 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)

Essa rotina calcula o coeficiente 3-j de Wigner,




ja
jb
jc
ma
mb
mc



onde os argumentos são fornecidos em unidades de meio inteiro, ja = two_ja/2, ma = two_ma/2, etc.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.8.2 Símbolos 6-j

Function: double gsl_sf_coupling_6j (int two_ja, int two_jb, int two_jc, int two_jd, int two_je, int two_jf)
Function: 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)

Essas rotina cacula o coeficiente 6-j de Wigner,





ja
jb
jc
jd
je
jf




onde os argumentos são dados em unidades de meio inteiro, ja = two_ja/2, ma = two_ma/2, etc.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.8.3 Símbolos 9-j

Function: 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)
Function: 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)

Essa rotina calcula o coeficiente 9-j de Wigner,





ja
jb
jc
jd
je
jf
jg
jh
ji




onde os argumentos são dados em unidades de meio inteiro, ja = two_ja/2, ma = two_ma/2, etc.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.9 Funções de Dawson

A integral de Dawson é definida por \exp(-x^2) \int_0^x dt \exp(t^2). Uma tabela da integral de Dawson pode ser encontrada em Abramowitz & Stegun, Tabela 7.5. As funções de Dawson são declaradas no arquivo de cabeçalho ‘gsl_sf_dawson.h’.

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

Essas rotinas calculam o valor da integral de Dawson para x.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.10 Funções de Debye

As fuções de Debye D_n(x) são definidas pela seguinte integral,

Dn(x) = n

xn

x

0 
dt tn

et − 1
Para informaçcões adicionais veja Abramowitz & Stegun, Seção 27.1. As funções de Debye são declaradas no arquivo de cabeçalho ‘gsl_sf_debye.h’.

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

Essas otinas calculam a função de Debye de primeira ordem D_1(x) = (1/x) \int_0^x dt (t/(e^t - 1)).

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

Essas rotinas calculam a função de Debye de segunda ordem D_2(x) = (2/x^2) \int_0^x dt (t^2/(e^t - 1)).

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

Essas rotinas calculam a função de Debye de terceira ordem D_3(x) = (3/x^3) \int_0^x dt (t^3/(e^t - 1)).

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

Essas rotinas calculam a função de Debye de quarta ordem D_4(x) = (4/x^4) \int_0^x dt (t^4/(e^t - 1)).

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

Essas rotinas calculam a função de Debye de quinta ordem D_5(x) = (5/x^5) \int_0^x dt (t^5/(e^t - 1)).

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

Essas rotinas calculam a função de Debye de sexta ordem D_6(x) = (6/x^6) \int_0^x dt (t^6/(e^t - 1)).


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.11 Dilogaritmo

As funções descritas nessa seção são declaradas no arquivo de cabeçalho ‘gsl_sf_dilog.h’.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.11.1 Argumento Real

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

Essas rotinas calculam o dilogaritmo para um argumento real. Na notação de Lewin isso é Li_2(x), a parte real do dilogaritmo de um real x. O dilogaritmo é definido através da representação da integral Li_2(x) = - \Re \int_0^x ds \log(1-s) / s. Note que \Im(Li_2(x)) = 0 for x <= 1, e -\pi\log(x) por x > 1.

Note que Abramowitz & Stegun referem-se à integral de Spence S(x)=Li_2(1-x) como o dilogaritmo em lugar de Li_2(x).


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.11.2 Argumento Complexo

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

Essa função calcula dilogaritmo para valores complexos com ambas as partes real e imaginária para o argumento complexo z = r \exp(i \theta). As partes reais e imaginárias do resultado são retornadas em result_re, result_im.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.12 Operações Elementares

As seguintes funções permitem estudar a propagação de erros quando combinando quantidades por multiplicação. As funções são declaradas no arquivo de cabeçalho ‘gsl_sf_elementary.h’.

Function: int gsl_sf_multiply_e (double x, double y, gsl_sf_result * result)

Essa função multiplica x e y armazenando o produto e seu erro associado em result.

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

Essa função multiplica x e y com erros absolutos associados dx e dy. O produto xy +/- xy \sqrt((dx/x)^2 +(dy/y)^2) é armazenado result.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.13 Integrais Elípticas

As funções descritas nessa seção são declaradas no arquivo de cabeçalho ‘gsl_sf_ellint.h’. Informações adicionais sobre as integrais elípticas podem ser encontradas em Abramowitz & Stegun, Capítulo 17.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.13.1 Formas de Legendre - Definição

As formas de Legendre das integrais elípticas F(\phi,k), E(\phi,k) e \Pi(\phi,k,n) são definidas por,

F(ϕ,k)
=
ϕ

0 
dt 1




(1 − k2 sin2(t))
E(ϕ,k)
=
ϕ

0 
dt

 

(1 − k2 sin2(t))
 
Π(ϕ,k,n)
=
ϕ

0 
dt 1

(1 + n sin2(t))


1 − k2 sin2(t)
As formas de Legendre completas são denotadas por K(k) = F(\pi/2, k) e E(k) = E(\pi/2, k).

A notação usada aqui é baseada em Carlson, Numerische Mathematik 33 (1979) 1 e difere levemente daquela usada por Abramowitz & Stegun, onde as funções são fornecidas em termos de parâmetro m = k^2 e n é substituído por -n.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.13.2 Formas de Carlson - Definição

As formas simétricas de Carlson das integrais elípticas RC(x,y), RD(x,y,z), RF(x,y,z) e RJ(x,y,z,p) são definidas por,

RC(x,y)
= 1/2


0 
dt (t+x)−1/2 (t+y)−1
RD(x,y,z)
= 3/2


0 
dt (t+x)−1/2 (t+y)−1/2 (t+z)−3/2
RF(x,y,z)
= 1/2


0 
dt (t+x)−1/2 (t+y)−1/2 (t+z)−1/2
RJ(x,y,z,p)
= 3/2


0 
dt (t+x)−1/2 (t+y)−1/2 (t+z)−1/2 (t+p)−1


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.13.3 Forma de Legendre de Integrais Elípticas Completas

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

Essas rotinas calculam a integral elíptica completa K(k) para a exatidão especificada pela varável de modo mode. Note que Abramowitz & Stegun definem essa função em termos do parâmetro m = k^2.

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

Essas rotinas calculam a integral elíptica completa E(k) para a exatidão especificada pela variável de modo mode. Note que Abramowitz & Stegun definem essa função em termos do parâmetro m = k^2.

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

Essas rotinas calculam a integral elíptica completa \Pi(k,n) para a exatidão especificada pela variável de modo mode. Note que Abramowitz & Stegun definem essa função em termos dos parâmetros m = k^2 e \sin^2(\alpha) = k^2, com a mudança de sinal n \to -n.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.13.4 Forma de Legendre de Integrais Elípticas Incompletas

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

Essas rotinas calculam a integral elíptica incompleta F(\phi,k) para a exatidão especificada pela variável de modo mode. Note que Abramowitz & Stegun definem essa função em termos do parâmetro m = k^2.

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

Essas rotinas calculam a integral elíptica incompleta E(\phi,k) para a exatidão especificada pela variável de modo mode. Note que Abramowitz & Stegun definem essa função em termos do parâmetro m = k^2.

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

Essas rotinas calculam a integral elíptica incompleta \Pi(\phi,k,n) para a exatidão especificada pela variável de modo mode. Note que Abramowitz & Stegun definem essa função em termos dos parâmetros m = k^2 e \sin^2(\alpha) = k^2, com a mudança de sinal n \to -n.

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

Essas funções calculam a integral elíptica incompleta D(\phi,k) que é definida através da forma de Carlson RD(x,y,z) pela seguinte relação,

D(ϕ,k,n) = 1

3
(sinϕ)3 RD (1−sin2(ϕ), 1−k2 sin2(ϕ), 1).
O argumento n não é usado e irá ser removido em uma versão futura.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.13.5 Formas de Carlson

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

Essas rotinas calculam a integral elíptica imcompleta RC(x,y) para a exatidão especificada pela variável de modo mode.

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

Essas rotinas calculam a integral elíptica incompleta RD(x,y,z) para a exatidão especificada pela variável de modo mode.

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

Essas rotinas calculam a integral elíptica incompleta RF(x,y,z) para a exatidão especificada pela variável de modo mode.

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

Essas rotinas calculam a integral elíptica incompleta RJ(x,y,z,p) para a exatidão especificada pela variável de modo mode.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.14 Funções Elípticas de Jacobi

As funções elípticas de Jacobi são definidas em Abramowitz & Stegun, Capítulo 16. As funções são declaradas no arquivo de cabeçalho ‘gsl_sf_elljac.h’.

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

Essa função calcula as funções elípticas Jacobianas sn(u|m), cn(u|m), dn(u|m) através das transformações descendentes de Landen.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.15 Funções de Erro

A função error está descrita em Abramowitz & Stegun, Caítulo 7. As funções nessa seção são declaradas no arquivo de cabeçalho ‘gsl_sf_erf.h’.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.15.1 Função de Erro

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

Essas rotinas calculam a função de erro erf(x), onde erf(x) = (2/\sqrt(\pi)) \int_0^x dt \exp(-t^2).


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.15.2 Função de Erro Complementar

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

Essas rotinas calculam a função de erro complementar erfc(x) = 1 - erf(x) = (2/\sqrt(\pi)) \int_x^\infty \exp(-t^2).


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.15.3 Registros Complementares a Funções de Erro

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

Essas rotinas calculam o logaritmo da função de erro complementar \log(\erfc(x)).


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.15.4 Funções de probabilidade

As funções de probalidade para a distribuição Normal ou distribuição de Gauss são descritas em Abramowitz & Stegun, Seção 26.2.

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

Essas rotinas calculam a função densidade de probabilidade de Gauss Z(x) = (1/\sqrt{2\pi}) \exp(-x^2/2).

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

Essas rotinas calculam o limite superior da função de probabilidade de Gauss Q(x) = (1/\sqrt{2\pi}) \int_x^\infty dt \exp(-t^2/2).

A função de risco para a distribuição normal, também conhecida como a razão de Mills inversa, é definida como,

h(x) = Z(x)

Q(x)
=   ⎛


2

π
 
exp(−x2 / 2)

\erfc(x/√2)
A função de risco decresce rapidamente à medida que x tende a -\infty e tende assintóticamente para h(x) \sim x à medida que x aproxíma-se de +\infty.

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

Essas rotinas calculam a função de risco para a distribuição normal.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.16 Funções Exponenciais

As funções descritas nessa seção são declaradas no arquivo de cabeçalho ‘gsl_sf_exp.h’.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.16.1 Função Exponencial

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

Essas rotinas fornecem uma função exponencial \exp(x) usando a semântica da GSL e verificação de erro.

Function: int gsl_sf_exp_e10_e (double x, gsl_sf_result_e10 * result)

Essa função calcula a exponencial \exp(x) usando o tipo gsl_sf_result_e10 para retornar um resultado com intervalo extendido. Essa função pode ser útil se o valor de \exp(x) puder vir a extrapolar o intervalo numérico de trabalho do tipo de dado double.

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

Essa rotina calcula a base neperiana “e” elevada ao expoente x e multiplica o resultado pelo fator y para retornar o produto y \exp(x).

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

Essa função calcula o produto y \exp(x) usando o tipo gsl_sf_result_e10 para retornar um resultado com intervalo numérico extendido.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.16.2 Funções Exponenciais Relativas

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

Essas rotinas calculam a quantidade \exp(x)-1 usando um algoritmo que é bastante aproximado para pequenos valores de x.

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

Essas rotinas calculam a quantidade (\exp(x)-1)/x usando um algoritmo que é bastante aproximado para equenos valores de x. Para pequenos x o algoritmo é baseado na expansão (\exp(x)-1)/x = 1 + x/2 + x^2/(2*3) + x^3/(2*3*4) + \dots.

Function: double gsl_sf_exprel_2 (double x)
Function: int gsl_sf_exprel_2_e (double x, gsl_sf_result * result)

Essas rotinas calculam a quantidade 2(\exp(x)-1-x)/x^2 usando um algoritmo que é bastante aproximado para pequenos valores de x. Para pequenos valores de x o algoritmo é baseado na expansão 2(\exp(x)-1-x)/x^2 = 1 + x/3 + x^2/(3*4) + x^3/(3*4*5) + \dots.

Function: double gsl_sf_exprel_n (int n, double x)
Function: int gsl_sf_exprel_n_e (int n, double x, gsl_sf_result * result)

Essas rotinas calculam a exponencial N-relativa, que é a generalização n-ésima das funções gsl_sf_exprel e gsl_sf_exprel_2. A exponencial N-relativa é fornecida por,

exprelN(x)
= N!/xN
exp(x) − N−1

k=0 
xk/k!
= 1 + x/(N+1) + x2/((N+1)(N+2)) + ...
= 1F1(1,1+N,x)


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.16.3 Exponenciação com Erro Estimado

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

Essa função calcula o resultado da base “e” elevada ao expoente x com um erro absoluto associado dx.

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

Essa função calcula o resultado da base “e” elevada ao expoente x com um erro absoluto associado dx usando o tipo gsl_sf_result_e10 para retornar um resultado com intervalo extendido.

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

Essa rotina calcula o produto y \exp(x) para as quantidades x, y com erros absolutos associados dx, dy.

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

Essa rotina calcula o produto y \exp(x) para as quantidades x, y com erros absolutos associados dx, dy usando o tipo gsl_sf_result_e10 para retornar um resultado com intervalo extendido.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.17 Integrais Exponenciais

Informações sobre integrais de exponenciais podem ser encontradas em Abramowitz & Stegun, Capítulo 5. Essas funções são declaradas no arquivo de cabeçalho ‘gsl_sf_expint.h’.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.17.1 Integral de Exponencial

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

Essas rotinas calculam a integral de exponencial E_1(x),

E1(x) : = ℜ


1 
dt exp(−xt)/t.

Function: double gsl_sf_expint_E2 (double x)
Function: int gsl_sf_expint_E2_e (double x, gsl_sf_result * result)

Essas rotinas calculam a integral de exponencial de segunda ordem E_2(x),

E2(x) : = ℜ


1 
dt exp(−xt)/t2.

Function: double gsl_sf_expint_En (int n, double x)
Function: int gsl_sf_expint_En_e (int n, double x, gsl_sf_result * result)

Essas rotinas calculam a integral de exponencial E_n(x) de ordem n,

En(x) : = ℜ


1 
dt exp(−xt)/tn.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.17.2 Ei(x)

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

Essas rotinas calculam a integral de exponencial Ei(x),

Ei(x) : = − PV



−x 
dt exp(−t)/t
onde PV denota o valor principal da integral.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.17.3 Integrais Hiperbólicas

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

Essas rotinas calculam a integral Shi(x) = \int_0^x dt \sinh(t)/t.

Function: double gsl_sf_Chi (double x)
Function: int gsl_sf_Chi_e (double x, gsl_sf_result * result)

Essas rotinas calculam a integral Chi(x) := \Re[ \gamma_E + \log(x) + \int_0^x dt (\cosh(t)-1)/t] , onde \gamma_E é a constante de Euler (disponível como a macro M_EULER).


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.17.4 Ei_3(x)

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

Essas rotinas calculam a integral de exponencial de terceira ordem Ei_3(x) = \int_0^xdt \exp(-t^3) for x >= 0.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.17.5 Integrais Trigonométricas

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

Essas rotinas calculam a integral do Seno Si(x) = \int_0^x dt \sin(t)/t.

Function: double gsl_sf_Ci (const double x)
Function: int gsl_sf_Ci_e (double x, gsl_sf_result * result)

Essas rotinas calculam a integral do Cosseno Ci(x) = -\int_x^\infty dt \cos(t)/t para x > 0.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.17.6 Integral de arco tangente

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

Essas rotinas calculam a integral do Arco-tangente, que é definida como AtanInt(x) = \int_0^x dt \arctan(t)/t.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.18 Funções de Fermi-Dirac

As funções descritas nessa seção são declaradas no arquivo de cabeçallho ‘gsl_sf_fermi_dirac.h’.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.18.1 Integrais de Fermi-Dirac Completas

A integral completa de Fermi-Dirac F_j(x) é fornecida por,

Fj(x) : = 1

Γ(j+1)



0 
dt tj

(exp(t−x) + 1)
Note que a integral de Fermi-Dirac é algumas vezes definida sem a normalização em outros textos.

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

Essas rotinas calculam a integral completa de Fermi-Dirac com um índice de -1. Essa integral é fornecida por F_{-1}(x) = e^x / (1 + e^x).

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

Essas rotinas calculam a integral completa de Fermi-Dirac com um índice de 0. Essa integral é fornecida por F_0(x) = \ln(1 + e^x).

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

Essas rotinas calculam a integral completa de Fermi-Dirac com um índice de 1, F_1(x) = \int_0^\infty dt (t /(\exp(t-x)+1)).

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

Essas rotinas calculam a integral completa de Fermi-Dirac com um índice de 2, F_2(x) = (1/2) \int_0^\infty dt (t^2 /(\exp(t-x)+1)).

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

Essas rotinas calculam a integral completa de Fermi-Dirac com um índice inteiro de j, F_j(x) = (1/\Gamma(j+1)) \int_0^\infty dt (t^j /(\exp(t-x)+1)).

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

Essas rotinas calculam a integral completa de Fermi-Dirac F_{-1/2}(x).

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

Essas rotinas calculam a integral completa de Fermi-Dirac F_{1/2}(x).

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

Essas rotinas calculam a integral completa de Fermi-Dirac F_{3/2}(x).


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.18.2 Integral de Fermi-Dirac Incompleta

A integral incompleta de Fermi-Dirac F_j(x,b) é fornecida por,

Fj(x,b) : = 1

Γ(j+1)



b 
dt tj

(exp(t−x) + 1)

Function: double gsl_sf_fermi_dirac_inc_0 (double x, double b)
Function: int gsl_sf_fermi_dirac_inc_0_e (double x, double b, gsl_sf_result * result)

Essas rotinas calculam a integral incompleta de Fermi-Dirac com um índice de zero, F_0(x,b) = \ln(1 + e^{b-x}) - (b-x).


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.19 Funções Gama e Beta

Essas rotinas adiante calculam as funções gama e beta em suas duas formas, completa e incompleta, bem como vários tipos de fatoriais. As funções descritas nessa seção são declaradas no arquivo de cabeçalho ‘gsl_sf_gamma.h’.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.19.1 Funções Gama

A função gama é definida através da seguinte integral,

Γ(x) =


0 
dt  tx−1 exp(−t)
A função gama é relacionada com a função fatorial através de \Gamma(n)=(n-1)! para inteiros positivos n. Informação adicional sobre a função Gama pode ser encontrado em Abramowitz & Stegun, Capítulo 6.

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

Essas rotinas calculam a função Gama \Gamma(x), na condição de que x não seja um inteiro negativo ou zero. A função é calculada usando o método real de Lanczos. O valor máximo de x tal que \Gamma(x) não é considerado um overflow é fornecido pela macro GSL_SF_GAMMA_XMAX e é 171.0.

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

Essas rotinas calculam o logaritmo da função Gama, \log(\Gamma(x)), na condição de que x não seja um inteiro negativo ou zero. Para x<0 a parte real de \log(\Gamma(x)) é retornada, a qual é equivalente a \log(|\Gamma(x)|). A função é calculada usando o método Lanczos real.

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

Essa rotina calcula o sinal da função gama e o logaritmo de seu módulo, na condição de que x não seja um inteiro negativo ou zero. A função é calculada usando o método de Lanczos real. O valor da função gama e seu erro pode ser reconstruído usando a relação \Gamma(x) = sgn * \exp(result\_lg), levando em consideração as duas componentes de result_lg.

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

Essas rotinas calculam a Função Gama regulada \Gamma^*(x) para x > 0. A função gama regulada é fornecida por,

Γ*(x)
= Γ(x)/(

 


 
x(x−1/2) exp(−x))
=
1 + 1

12x
+ ...
   for  x→ ∞
e é uma sugestão útil de Temme (20).

Function: double gsl_sf_gammainv (double x)
Function: int gsl_sf_gammainv_e (double x, gsl_sf_result * result)

Essas rotinas calculam o recíproco da função gama, 1/\Gamma(x) usando o método de Lanczos real.

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

Essa rotina calcula \log(\Gamma(z)) para o complexo z=z_r+i z_i e z sendo inteiro estritamente positivo, usando o método de Lanczos complexo. Os parâmetros retornados são lnr = \log|\Gamma(z)| e arg = \arg(\Gamma(z)) dentro do intervalo (-\pi,\pi]. Note que a parte de fase (arg) não é bem determinado quando |z| é muito grande, devido a arredondamentos inevitáveis na restrição (-\pi,\pi]. Irá resultar em um erro GSL_ELOSS quando isso ocorrer. A parte de valor absoluto (lnr), todavia, nunca sofre perda de precisão.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.19.2 Fatoriais

Embora fatoriais possam ser calculados a partir da função Gama, usando a relação n! = \Gamma(n+1) para o inteiro não negativo n, esse cálculo de fatorial é mais eficiente se for feito por meio das funções nessa seção, particularmente para pequenos valores de n, os quais valores fatoriais são mantidos em tabelas codificadas.

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

Essas rotinas calculam o fatorial n!. O fatorial é relacionado com a função Gama através de n! = \Gamma(n+1). O valor máximo de n é tal que n! não seja considerado um overflow é fornecido pela macro GSL_SF_FACT_NMAX e é 170.

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

Essas rotinas calculam o fatorial duplo n!! = n(n-2)(n-4) \dots. o valor máximo de n tal que n!! não seja considerado um overflow é fornecido pela macro GSL_SF_DOUBLEFACT_NMAX e é 297.

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

Essas rotinas calculam o logaritmo do fatorial de n, \log(n!). O algoritmo é mais rápido que calcular \ln(\Gamma(n+1)) via gsl_sf_lngamma para n < 170, mas esse comportamento muda para grandes valores de n.

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

Essas rotinas calculam o logaritmo do fatorial duplo de n, \log(n!!).

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

Essas rotinas calculam a combinação simples entre n escolha m = n!/(m!(n-m)!)

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

Essas rotinas calculam o logaritmo da combinação simples. Essas rotinas equivalem ao somatório \log(n!) - \log(m!) - \log((n-m)!).

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

Essas rotinas calculam o coeficiente de Taylor x^n / n! para x >= 0, n >= 0.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.19.3 Símbolo de Pochhammer

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

Essas rotinas calculam o símbolo de Pochhammer (a)_x = \Gamma(a + x)/\Gamma(a). O símbolo de Pochhammer é também conhecido como símbolo de Apell e algumas vezes escrito como (a,x). Quando a e a+x são inteiros negativos ou zero, o valor limitante da razão é retornado.

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

Essas rotinas calculam o logaritmo do símbolo de Pochhammer, \log((a)_x) = \log(\Gamma(a + x)/\Gamma(a)).

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

Essas rotinas calculam o sinal do símbolo de Pochhammer e o logaritmo de seu módulo. Os parâmetros calculados são result = \log(|(a)_x|) com um correspondente termo de erro, e sgn = \sgn((a)_x) onde (a)_x = \Gamma(a + x)/\Gamma(a).

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

Essas rotinas calculam símbolo de Pochhammer relativo ((a)_x - 1)/x onde (a)_x = \Gamma(a + x)/\Gamma(a).


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.19.4 Funções Gamma Incompletas

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

Essas funções calculam a Função Gama incompleta não normalizada \Gamma(a,x) = \int_x^\infty dt t^{a-1} \exp(-t) para a real e x >= 0.

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

Essas rotinas calculam Função Gama incompleta normalizada Q(a,x) = 1/\Gamma(a) \int_x^\infty dt t^{a-1} \exp(-t) para a > 0, x >= 0.

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

Essas rotinas calculam a Função Gama incompleta normalizada complementar P(a,x) = 1 - Q(a,x) = 1/\Gamma(a) \int_0^x dt t^{a-1} \exp(-t) for a > 0, x >= 0.

Note que Abramowitz & Stegun chamam P(a,x) à função gama incompleta (seção 6.5).


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.19.5 Funções Beta

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

Essas rotinas calculam a Função Beta, B(a,b) = \Gamma(a)\Gamma(b)/\Gamma(a+b) na condição de que a e b não sejam inteiros negativos.

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

Essas rotinas calculam o logaritmo da Função Beta, \log(B(a,b)) na condição que a e b não sejam inteiros negativos.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.19.6 Função Beta Incompleta

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

Essas rotinas calculam função Beta incompleta normalizada I_x(a,b)=B_x(a,b)/B(a,b) onde B_x(a,b) = \int_0^x t^{a-1} (1-t)^{b-1} dt para 0 <= x <= 1. Para a > 0, b > 0 o valor é calculado usando uma expansão contínua de fração. Para todos os outros valores a função Beta incompleta normalizada é calculada usando a relação I_x(a,b,x) = (1/a) x^a 2F1(a,1-b,a+1,x)/B(a,b).


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.20 Funções de Gegenbauer

Os polinômios de Gegenbauer são definidos em Abramowitz & Stegun, Capítulo 22, onde eles são conhecidos polinômios Ultraesféricos. As funções descritas nessa seção são declaradas no arquivo de cabeçalho ‘gsl_sf_gegenbauer.h’.

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

Essas funções avaliam os polinômios de Gegenbauer C^{(\lambda)}_n(x) usando representações explícitas para n =1, 2, 3.

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

Essas funções avaliam o polinômio de Gegenbauer C^{(\lambda)}_n(x) para um específico valor de n, lambda, x na condição que \lambda > -1/2, n >= 0.

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

Essa função calcula um vetor estático dos polinômios de Gegenbauer C^{(\lambda)}_n(x) for n = 0, 1, 2, \dots, nmax, na condição que \lambda > -1/2, nmax >= 0.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.21 Funções Hipergeométricas

Funções Hipergeométricas são descritas em Abramowitz & Stegun, Capítulos 13 e 15. Essas funções são declaradas no arquivo de cabeçalho ‘gsl_sf_hyperg.h’.

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

Essas rotinas calculam a função hipergeométrica 0F1(c,x).

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

Essas rotinas calculam a função hipergeométrica confluente 1F1(m,n,x) = M(m,n,x) for integer parameters m, n.

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

Essas rotinas calculam a função hipergeométrica confluente 1F1(a,b,x) = M(a,b,x) para parâmetros gerais a, b.

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

Essas rotinas calculam a função hipergeométrica confluente U(m,n,x) para parâmetros inteiros m, n.

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

Essas rotinas calculam a função hipergeométrica confluente U(m,n,x) para parâmetros inteiros m, n usando o tipo gsl_sf_result_e10 para retornar um resultado com intervalo extendido.

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

Essas rotinas calculam a função hipergeométrica confluente U(a,b,x).

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

Essa rotina calcula a função hipergeométrica confluente U(a,b,x) usando o tipo gsl_sf_result_e10 para retornar um resultado com intervalo extendido.

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

Essas rotinas calculam a função hipergeométrica de Gauss 2F1(a,b,c,x) = F(a,b,c,x) para |x| < 1.

Se os argumentos (a,b,c,x) forem muito próximos a uma singularidade então a função pode retornar o código de erro GSL_EMAXITER quando a série de aproximação converge muito lentamente. Isso ocorre na região de x=1, c - a - b = m para o inteiro m.

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

Essas rotinas calculam a função hipergeométrica de Gauss 2F1(a_R + i a_I, a_R - i a_I, c, x) com parâmetros complexos para |x| < 1.

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

Essas rotinas calculam a função hipergeométrica de Gauss renormalizada 2F1(a,b,c,x) / \Gamma(c) para |x| < 1.

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

Essas rotinas calculam a função hipergeométrica de Gauss renormalizada 2F1(a_R + i a_I, a_R - i a_I, c, x) / \Gamma(c) for |x| < 1.

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

Essas rotinas calculam a função hipergeométrica 2F0(a,b,x). A série de representação é uma série hipergeométrica convergente. Todavia, para x < 0 temos 2F0(a,b,x) = (-1/x)^a U(a,1+a-b,-1/x)


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.22 Funções de Laguerre

Os polinômios de Laguerre generalizados são definidos em termos de funções hipergeométricas confluentes como L^a_n(x) = ((a+1)_n / n!) 1F1(-n,a+1,x), e são algumas vezes referidos como os polinômios de Laguerre associados. Eles são relacionado aos polinômios de Laguerre planos L_n(x) através de L^0_n(x) = L_n(x) e de L^k_n(x) = (-1)^k (d^k/dx^k) L_(n+k)(x). Para mais informações veja Abramowitz & Stegun, Capítulo 22.

As funções descritas nessa seção são declaradas no arquivo de cabeçalho ‘gsl_sf_laguerre.h’.

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

Essas rotinas avaliam os polinômios de Laguerre generalizados L^a_1(x), L^a_2(x), L^a_3(x) usando representações explícitas.

Function: double gsl_sf_laguerre_n (const int n, const double a, const double x)
Function: int gsl_sf_laguerre_n_e (int n, double a, double x, gsl_sf_result * result)

Essas rotinas avaliam os polinômios de Laguerre generalizados L^a_n(x) para a > -1, n >= 0.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.23 Funçãoes W de Lambert

As funções W de lambert, W(x), são definidas para serem soluções da equação W(x) \exp(W(x)) = x. Essa função tem múltiplos ramos para x < 0; todavia, tem somente dois ramos com valores reais. Definimos W_0(x) para ser o ramo principal, onde W > -1 para x < 0, e W_{-1}(x) para ser o outros ramo real, onde W < -1 para x < 0. As funções de Lambert são declaradas no arquivo de cabeçalho ‘gsl_sf_lambert.h’.

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

Essas funções calculam o ramo principal da função W de Lambert, W_0(x).

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

Essas funções calculam o ramo secundário de valores reais da função W de Lambert, W_{-1}(x).


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.24 Funções de Legendre e Harmônicos Esféricos

As Funções de Legendre e os Polinômios de Legendre são descritos em Abramowitz & Stegun, Capítulo 8. Essas funções são declaradas no arquivo de cabeçalho ‘gsl_sf_legendre.h’.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.24.1 Polinômios de Legendre

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

Essas funções avaliam os polinômios de Legendre P_l(x) usando representações explícitas para l=1, 2, 3.

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

Essas funções avaliam o polinômio de Legendre P_l(x) para o valor específico de l, x na condição de que l >= 0, |x| <= 1

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

Essas funções calculam vetores estáticos de polinômios de Legendre P_l(x) e derivadas dP_l(x)/dx, for l = 0, \dots, lmax, |x| <= 1

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

Essas rotina calcula a função de Legendre Q_0(x) para x > -1, x != 1.

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

Essas rotinas calculam a função de Legendre Q_1(x) para x > -1, x != 1.

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

Essas rotinas calculam a função de Legendre Q_l(x) para x > -1, x != 1 e l >= 0.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.24.2 Polinômios Associados de Legendre e Harmônicos Esféricos

As seguintes funções calculam os polinômios de Legendre associados P_l^m(x). Note que essa função cresce combinatorialmente com l e pode causar overflow em l maior que 150. Não existe nenhum problema para pequenos valores de m, mas overflow ocorre quando m e l são ambos grandes. Em lugar de permitir overflows, essas funções recusam-se a calcular P_l^m(x) e retornam GSL_EOVRFLW quando elas percebem que l e m são muito grandes.

Se você deseja calcular uma harmônica esférica, então não use essas funções. Ao invés de usá-las use gsl_sf_legendre_sphPlm mostrada adiante, que utiliza uma recursão similar, mas com as funções normalizadas.

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

Essas rotinas calculam o polinômio de Legendre associado P_l^m(x) for m >= 0, l >= m, |x| <= 1.

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

Essas funções calculam vetores estáticos de polinômios de Legendre P_l^m(x) e derivadas dP_l^m(x)/dx, for m >= 0, l = |m|, ..., lmax, |x| <= 1.

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

Essas rotinas calculam o polinômios de Legendre associado normalizado \sqrt{(2l+1)/(4\pi)} \sqrt{(l-m)!/(l+m)!} P_l^m(x) adequado para ser usado em harmônicas esféricas. Os parâmetros devem satisfazer m >= 0, l >= m, |x| <= 1. Essas rotinas evitam os overflows que ocorrem para a normalização padrão de P_l^m(x).

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

Essas funções calculam vetores estáticos de funções de Legendre associadas normalizadas \sqrt{(2l+1)/(4\pi)} \sqrt{(l-m)!/(l+m)!} P_l^m(x), e derivadas, para m >= 0, l = |m|, ..., lmax, |x| <= 1.0

Function: int gsl_sf_legendre_array_size (const int lmax, const int m)

Essa função retorna o tamanho de result_array[] necessário para versões do vetor estático de P_l^m(x), lmax - m + 1. Uma versão modificada dessa função é usada quando HAVE_INLINE for definida.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.24.3 Funções Canônicas

As Funções Canônicas P^\mu_{-(1/2)+i\lambda}(x) e Q^\mu_{-(1/2)+i\lambda} são descritas em Abramowitz & Stegun, Seção 8.12.

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

Essas rotinas calculam a Função Canônica Esférica irregular P^{1/2}_{-1/2 + i \lambda}(x) for x > -1.

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

Essas rotinas calculam a Função Canônica Esférica regular P^{-1/2}_{-1/2 + i \lambda}(x) for x > -1.

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

Essas rotinas calculam a função canônica P^0_{-1/2 + i \lambda}(x) para x > -1.

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

Essas rotinas calculam a função canônica P^1_{-1/2 + i \lambda}(x) for x > -1.

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

Essas rotinas calculam a Função Canônica Esférica Regular P^{-1/2-l}_{-1/2 + i \lambda}(x) para x > -1, l >= -1.

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

Essas rotinas calculam a Função Canônica Cilíndrica Regular P^{-m}_{-1/2 + i \lambda}(x) for x > -1, m >= -1.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.24.4 Funções Radiais no Espaço Hiperbólico

As seguintes funções esféricas são especializações das funções de Legendre as quais fornecem autofunções regulares do Laplaciano sobre um espaço hiperbólico tridimencional H3d. De particular interesse é o limite uniforme, \lambda \to \infty, \eta \to 0, \lambda\eta fixado.

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

Essas rotinas calculam a autofunção de zero-ésimo radial do Laplaciano sobre o espaço hiperbólico tridimensional, L^{H3d}_0(\lambda,\eta) := \sin(\lambda\eta)/(\lambda\sinh(\eta)) para \eta >= 0. No limite plano isso toma a forma L^{H3d}_0(\lambda,\eta) = j_0(\lambda\eta).

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

Essas rotinas calculam a primeira autofunção radial do Laplaciano sobre o espaço hiperbólico tridimensional, L^{H3d}_1(\lambda,\eta) := 1/\sqrt{\lambda^2 + 1} \sin(\lambda \eta)/(\lambda \sinh(\eta)) (\coth(\eta) - \lambda \cot(\lambda\eta)) para \eta >= 0. No limite plano isso toma a forma L^{H3d}_1(\lambda,\eta) = j_1(\lambda\eta).

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

Essas rotinas calculam a l-ésima autofunção radial do Laplaciano sobre o espaço hiperbólico tridimencional \eta >= 0, l >= 0. No limite plano isso toma a forma L^{H3d}_l(\lambda,\eta) = j_l(\lambda\eta).

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

Essa função calcula um vetor estático de autofunções radiais L^{H3d}_l(\lambda, \eta) for 0 <= l <= lmax.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.25 Logaritmos e Funções Relacionadas

Informações sobre as propriedades da função Logaritmica podem ser encontradas em Abramowitz & Stegun, Capítulo 4. As funções descritas nessa seção são declaradas no arquivo de cabeçalho ‘gsl_sf_log.h’.

Function: double gsl_sf_log (double x)
Function: int gsl_sf_log_e (double x, gsl_sf_result * result)

Essas rotinas calculam o logaritmo de x, \log(x), para x > 0.

Function: double gsl_sf_log_abs (double x)
Function: int gsl_sf_log_abs_e (double x, gsl_sf_result * result)

Essas rotinas calculam o logaritmo do módulo de x, \log(|x|), for x \ne 0.

Function: int gsl_sf_complex_log_e (double zr, double zi, gsl_sf_result * lnr, gsl_sf_result * theta)

Essa rotina calcula o logaritmo complexo de z = z_r + i z_i. Os resultados são retornados como lnr, theta tal que \exp(lnr + i \theta) = z_r + i z_i, onde \theta localiza-se no intervalo [-\pi,\pi].

Function: double gsl_sf_log_1plusx (double x)
Function: int gsl_sf_log_1plusx_e (double x, gsl_sf_result * result)

Essas rotinas calculam \log(1 + x) para x > -1 usando um algoritmo que tem grande exatidão para pequenos valores de x.

Function: double gsl_sf_log_1plusx_mx (double x)
Function: int gsl_sf_log_1plusx_mx_e (double x, gsl_sf_result * result)

Essas rotinas calculam \log(1 + x) - x para x > -1 usando um algoritmo que tem grande exatidão para pequenos valores de x.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.26 Funções de Mathieu

As rotinas descritas nessa seção calculam a Função de Mathieu angular e a Função de Mathieu radial, e seus valores característicos. As funções de Mathieu são as soluções das duas seguintes equações diferenciais:

d2 y

d v2
+ (a − 2qcos2v)y = 0,
d2 f

d u2
− (a − 2qcosh2u)f = 0.
As funções angulares de Mathieu ce_r(x,q), se_r(x,q) são as soluções periódicas pares e ímpares da primeira equação, que é conhecida como equação de Mathieu. Essas soluções existem somente para a sequência discreta de valores característicos a=a_r(q) (periódicos pares) e a=b_r(q) (periódicos ímpares).

O radial das funções de Mathieu Mc^{(j)}_{r}(z,q), Ms^{(j)}_{r}(z,q) são as soluções da segunda equação, que é referenciada como equação de Mathieu modificada. O radial das funções de Mathieu de primeiro, segundo, terceiro e quarto tipo são denotadas pelo parâmetro j, que recebe os valores 1, 2, 3 ou 4.

Para maiores informações sobre as funções de Mathieu, veja Abramowitz e Stegun, Capítulo 20. Essas funções são definidas no arquivo de cabeçalho ‘gsl_sf_mathieu.h’.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.26.1 Função de Mathieu - espaço de trabalho

As funções de Mathieu podem ser calculadas para uma única ordem ou para múltiplas ordens, usando rotinas baseadas em vetores estáticos. As rotinas baseadas em vetores estáticos requerem um espaço de trabalho pré-alocado.

Function: gsl_sf_mathieu_workspace * gsl_sf_mathieu_alloc (size_t n, double qmax)

Essa função retorna um espaço de trabalho para as versões baseadas em vetores estáticos das rotinas de Mathieu. Os argumentos n e qmax especificam a maior ordem e o valor q das funções de Mathieu que pode ser calculado com esse espaço de trabalho.

Function: void gsl_sf_mathieu_free (gsl_sf_mathieu_workspace * work)

Essa função libera o espaço de trabalho work.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.26.2 Função de Mathieu - Valores Característicos

Function: int gsl_sf_mathieu_a (int n, double q, gsl_sf_result * result)
Function: int gsl_sf_mathieu_b (int n, double q, gsl_sf_result * result)

Essas rotinas calculam os valores característicos a_n(q), b_n(q) das funções de Mathieu ce_n(q,x) and se_n(q,x), respectivamente.

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

Essas rotinas calculam uma série de valores característicos de Mathieu a_n(q), b_n(q) para n variando de order_min a order_max inclusive, armazenando os resultados no vetor estático result_array.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.26.3 Funções angulares de Mathieu

Function: int gsl_sf_mathieu_ce (int n, double q, double x, gsl_sf_result * result)
Function: int gsl_sf_mathieu_se (int n, double q, double x, gsl_sf_result * result)

Essas rotinas calculam as funções de Mathieu angulares ce_n(q,x) e se_n(q,x), respectivamente.

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

Essas rotinas calculam uma série de funções de Mathieu angulares ce_n(q,x) e se_n(q,x) de ordem n com n variando de nmin a nmax inclusive, armazenando os resultados no veotr estático result_array.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.26.4 Funções Radiais de Mathieu

Function: int gsl_sf_mathieu_Mc (int j, int n, double q, double x, gsl_sf_result * result)
Function: int gsl_sf_mathieu_Ms (int j, int n, double q, double x, gsl_sf_result * result)

Essas rotinas calculam as funções de Mathieu radial de j-ésimo tipo Mc_n^{(j)}(q,x) e Ms_n^{(j)}(q,x) de ordem n.

Os valores permitidos de j são 1 e 2. As funções de j = 3,4 podem ser calculadas como M_n^{(3)} = M_n^{(1)} + iM_n^{(2)} e M_n^{(4)} = M_n^{(1)} - iM_n^{(2)}, onde M_n^{(j)} = Mc_n^{(j)} or Ms_n^{(j)}.

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

Essas rotinas calculam uma série de funções de Mathieu radial de tipo j, com ordem variando de nmin a nmax inclusive, armazenado os resultados no vetor estático result_array.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.27 Função Potência

As seguintes funções são equivalentes à função gsl_pow_int (veja seção Potência de pequenos inteiros) com uma estimativa de erro. Essas funções são declaradas no arquivo de cabeçalho ‘gsl_sf_pow_int.h’.

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

Essas rotinas calculam a potência x^n para o inteiro n. A potência é calculada usando o menor número de multiplicações. Por exemplo, x^8 é calculado como ((x^2)^2)^2, requerendo somente 3 multiplicações. Por razões de eficiência, essas funções não verificam por condições de overflow ou underflow.

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

[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.28 Função Psi (Digama)

As funções poligama de ordem n são definidas por

ψ(n)(x) =
d

dx

n

 
ψ(x) =
d

dx

n+1

 
log(Γ(x))
onde \psi(x) = \Gamma'(x)/\Gamma(x) é conhecida como função digama. Essas funções são declaradas no arquivo de cabeçalho ‘gsl_sf_psi.h’.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.28.1 Função Digama

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

Essas rotinas calculam a função digama \psi(n) para inteiros positivos n. A função digama é também chamada de função Psi.

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

Essas rotinas calculam a função digama \psi(x) para valores de x quaisquer, x \ne 0.

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

Essas rotinas calculam a parte real da função digama sobre a linha 1+i y, \Re[\psi(1 + i y)].


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.28.2 Função Trigama

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

Essas rotinas calculam a função Trigama \psi'(n) para o inteiro positivo n.

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

Essas rotinas calculam função Trigamma \psi'(x) para valores de x quaisquer.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.28.3 Função Poligama

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

Essas rotinas calculam a função poligama \psi^{(n)}(x) para n >= 0, x > 0.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.29 Funções Sincrotron

As funções descritas nessa seção são declaradas no arquivo de cabeçalho ‘gsl_sf_synchrotron.h’.

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

Essas rotinas calculam a primeira função síncrotron x \int_x^\infty dt K_{5/3}(t) para x >= 0.

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

Essas rotinas calculam a segunda função síncrotron x K_{2/3}(x) for x >= 0.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.30 Funções de Transporte

As funções de transporte J(n,x) são definidas através das representações de integral J(n,x) := \int_0^x dt t^n e^t /(e^t - 1)^2. As funções de transporte são declaradas no arquivo de cabeçalo ‘gsl_sf_transport.h’.

Function: double gsl_sf_transport_2 (double x)
Function: int gsl_sf_transport_2_e (double x, gsl_sf_result * result)

Essas rotinas calculam a função de transporte J(2,x).

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

Essas rotinas calculam a função de transporte J(3,x).

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

Essas rotinas calculam a função de transporte J(4,x).

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

Essas rotinas calculam a função de transporte J(5,x).


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.31 Funções Trigonométricas

A biblioteca inclui suas próprias funções trigonométricas com o objetivo de fornecer consistência independente da plataforma de instalação e estmativas de eroo confiáveis. Essas funções são declaradas no arquivo de cabeçalho ‘gsl_sf_trig.h’.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.31.1 Funções Trigonométricas Circulares

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

Essas rotinas calculam a função seno \sin(x).

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

Essas rotinas calculam a função cosseno \cos(x).

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

Essas rotinas calculam função hipotenusa \sqrt{x^2 + y^2} evitando overflow e underflow.

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

Essas rotinas calculam \sinc(x) = \sin(\pi x) / (\pi x) para qualquer valor de x.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.31.2 Funções Trigonométricas para Argumentos Complexos

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

Essa função calcula o seno complexo, \sin(z_r + i z_i) armazenando as partes real e imaginária em szr, szi.

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

Essa função calcula o cosseno complexo, \cos(z_r + i z_i) armazenando as partes real e imaginária em czr, czi.

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

Essa função calcula o seno de logaritmo complexo, \log(\sin(z_r + i z_i)) armazenando as partes real e imaginária em lszr, lszi.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.31.3 Funções Trigonométricas Hiperbólicas

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

Essas rotinas calculam \log(\sinh(x)) para x > 0.

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

Essas rotinas calculam \log(\cosh(x)) para qualquer x.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.31.4 Funções de Conversão

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

Essa função converte as coordenadas polares (r,theta) para coordenadas retangulares (x,y), x = r\cos(\theta), y = r\sin(\theta).

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

Essa função converte as coordenadas retangulares (x,y) para coordenadas polares (r,theta), tais que x = r\cos(\theta), y = r\sin(\theta). O argumento theta encontra-se no intervalo [-\pi, \pi].


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.31.5 Funções de Restrição

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

Essas rotinas forçam o ângulo theta a localizar-se no intervalo (-\pi,\pi].

Note que o valor matemático de \pi é ligeiramente maior que M_PI, de forma que os números programados no computador M_PI e -M_PI estejam incluídos no intervalo.

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

Essas rotinas forçam o ângulo theta a localizar-se no intervalo [0, 2\pi).

Note que o valore matemático de 2\pi é ligeiramente maior que 2*M_PI, de forma que o número programado no computador 2*M_PI esteja incluído no intervalo.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.31.6 Funções Trigonométricas Com Estimativas de Erro

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

Essa rotina calcula o seno do ângulo x com um erro abosoluto associado dx, \sin(x \pm dx). Note que essa função é fornecida na forma de controlador de erro somente uma vez que seu propósito é calcular o erro propagado.

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

Essa rotina calcula o cosseno de um ângulo x com um erro absoluto associado dx, \cos(x \pm dx). Note que essa função é fornecida na forma de controlador de erro somente uma vez que seu propósito é calcular o erro propagado.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.32 Funções Zeta

A função zeta de Riemann é definida em Abramowitz & Stegun, Seção 23.2. As funções descritas nessa seção são declaradas no arquivo de cabeçalho ‘gsl_sf_zeta.h’.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.32.1 Função Zeta de Riemann

A função zeta de Riemann é definida através do somatório infinito \zeta(s) = \sum_{k=1}^\infty k^{-s}.

Function: double gsl_sf_zeta_int (int n)
Function: int gsl_sf_zeta_int_e (int n, gsl_sf_result * result)

Essas rotinas calculam a função zeta de Riemann \zeta(n) para o inteiro n, n \ne 1.

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

Essas rotinas calculam a função zeta de Riemann \zeta(s) para um s arbitrário, s \ne 1.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.32.2 Função Zeta de Riemann Menos Um

Para argumentos positivos elevados, a função zeta de Riemann aproxima-se da unidade. Nessa região a parte fracionária é interessante, e portanto precisamos de uma função para avaliar essa parte fracionária explicitamente.

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

Essas rotinas calculam \zeta(n) - 1 para o inteiro n, n \ne 1.

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

Essas rotinas calculam \zeta(s) - 1 para valores arbitrários de s, s \ne 1.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.32.3 Função Zeta de Hurwitz

A função zeta de Hurwitz é definida através de \zeta(s,q) = \sum_0^\infty (k+q)^{-s}.

Function: double gsl_sf_hzeta (double s, double q)
Function: int gsl_sf_hzeta_e (double s, double q, gsl_sf_result * result)

Essas rotinas calculam função zeta de Hurwitz \zeta(s,q) para s > 1, q > 0.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.32.4 Função Eta

A função eta é definida através de \eta(s) = (1-2^{1-s}) \zeta(s).

Function: double gsl_sf_eta_int (int n)
Function: int gsl_sf_eta_int_e (int n, gsl_sf_result * result)

Essas rotinas calculam a função eta \eta(n) para o inteiro n.

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

Essas rotinas calculam a função eta \eta(s) para um valor arbitrário de s.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.33 Exemplos

Os sequintes exemplos demonstram o uso da forma controladora de erro das funções especiais, no caso para calcular função de Bessel 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;
}

Aqui está o resultado da execução do programa,

$ ./a.out 
status  = success
J0(5.0) = -0.177596771314338292 
      +/-  0.000000000000000193
exact   = -0.177596771314338292

O programa seguinte calcula a mesma quantidade usando a forma natural da função. Nesse caso o termo de erro result.err e a situação de retorno não estão acessíveis.

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

int
main (void)
{
  double x = 5.0;
  double expected = -0.17759677131433830434739701;
  
  double y = gsl_sf_bessel_J0 (x);

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

Os resultados da função são os mesmos,

$ ./a.out 
J0(5.0) = -0.177596771314338292
exact   = -0.177596771314338292

[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

7.34 Referências e Leituras Adicionais

A biblioteca segue as convenções de Abramowitz & Stegun onde possível,

Os seguintes artigos possuem informação sobre o algoritmo usado para calcular as funções especiais,


[ << ] [ >> ]           [Topo] [Conteúdo] [Índice] [ ? ]

Esse documento foi gerado em 23 de Julho de 2013 usando texi2html 5.0.