(octave.info)Special Functions


Next: Rational Approximations Prev: Utility Functions Up: Arithmetic
Enter node , (file) or (file)node

17.6 Special Functions
======================

 -- : [A, IERR] = airy (K, Z, OPT)
     Compute Airy functions of the first and second kind, and their
     derivatives.

           K   Function   Scale factor (if "opt" is supplied)
          ---  --------   ---------------------------------------
           0   Ai (Z)     exp ((2/3) * Z * sqrt (Z))
           1   dAi(Z)/dZ  exp ((2/3) * Z * sqrt (Z))
           2   Bi (Z)     exp (-abs (real ((2/3) * Z * sqrt (Z))))
           3   dBi(Z)/dZ  exp (-abs (real ((2/3) * Z * sqrt (Z))))

     The function call ‘airy (Z)’ is equivalent to ‘airy (0, Z)’.

     The result is the same size as Z.

     If requested, IERR contains the following status information and is
     the same size as the result.

       0. Normal return.

       1. Input error, return ‘NaN’.

       2. Overflow, return ‘Inf’.

       3. Loss of significance by argument reduction results in less
          than half of machine accuracy.

       4. Loss of significance by argument reduction, output may be
          inaccurate.

       5. Error—no computation, algorithm termination condition not met,
          return ‘NaN’.

 -- : J = besselj (ALPHA, X)
 -- : J = besselj (ALPHA, X, OPT)
 -- : [J, IERR] = besselj (...)
     Compute Bessel functions of the first kind.

     The order of the Bessel function ALPHA must be real.  The points
     for evaluation X may be complex.

     If the optional argument OPT is 1 or true, the result J is
     multiplied by ‘exp (-abs (imag (X)))’.

     If ALPHA is a scalar, the result is the same size as X.  If X is a
     scalar, the result is the same size as ALPHA.  If ALPHA is a row
     vector and X is a column vector, the result is a matrix with
     ‘length (X)’ rows and ‘length (ALPHA)’ columns.  Otherwise, ALPHA
     and X must conform and the result will be the same size.

     If requested, IERR contains the following status information and is
     the same size as the result.

       0. Normal return.

       1. Input error, return ‘NaN’.

       2. Overflow, return ‘Inf’.

       3. Loss of significance by argument reduction results in less
          than half of machine accuracy.

       4. Loss of significance by argument reduction, output may be
          inaccurate.

       5. Error—no computation, algorithm termination condition not met,
          return ‘NaN’.

     See also: Note: bessely, Note: besseli,
     Note: besselk, Note: besselh.

 -- : Y = bessely (ALPHA, X)
 -- : Y = bessely (ALPHA, X, OPT)
 -- : [Y, IERR] = bessely (...)
     Compute Bessel functions of the second kind.

     The order of the Bessel function ALPHA must be real.  The points
     for evaluation X may be complex.

     If the optional argument OPT is 1 or true, the result Y is
     multiplied by ‘exp (-abs (imag (X)))’.

     If ALPHA is a scalar, the result is the same size as X.  If X is a
     scalar, the result is the same size as ALPHA.  If ALPHA is a row
     vector and X is a column vector, the result is a matrix with
     ‘length (X)’ rows and ‘length (ALPHA)’ columns.  Otherwise, ALPHA
     and X must conform and the result will be the same size.

     If requested, IERR contains the following status information and is
     the same size as the result.

       0. Normal return.

       1. Input error, return ‘NaN’.

       2. Overflow, return ‘Inf’.

       3. Loss of significance by argument reduction results in less
          than half of machine accuracy.

       4. Complete loss of significance by argument reduction, return
          ‘NaN’.

       5. Error—no computation, algorithm termination condition not met,
          return ‘NaN’.

     See also: Note: besselj, Note: besseli,
     Note: besselk, Note: besselh.

 -- : I = besseli (ALPHA, X)
 -- : I = besseli (ALPHA, X, OPT)
 -- : [I, IERR] = besseli (...)
     Compute modified Bessel functions of the first kind.

     The order of the Bessel function ALPHA must be real.  The points
     for evaluation X may be complex.

     If the optional argument OPT is 1 or true, the result I is
     multiplied by ‘exp (-abs (real (X)))’.

     If ALPHA is a scalar, the result is the same size as X.  If X is a
     scalar, the result is the same size as ALPHA.  If ALPHA is a row
     vector and X is a column vector, the result is a matrix with
     ‘length (X)’ rows and ‘length (ALPHA)’ columns.  Otherwise, ALPHA
     and X must conform and the result will be the same size.

     If requested, IERR contains the following status information and is
     the same size as the result.

       0. Normal return.

       1. Input error, return ‘NaN’.

       2. Overflow, return ‘Inf’.

       3. Loss of significance by argument reduction results in less
          than half of machine accuracy.

       4. Complete loss of significance by argument reduction, return
          ‘NaN’.

       5. Error—no computation, algorithm termination condition not met,
          return ‘NaN’.

     See also: Note: besselk, Note: besselj,
     Note: bessely, Note: besselh.

 -- : K = besselk (ALPHA, X)
 -- : K = besselk (ALPHA, X, OPT)
 -- : [K, IERR] = besselk (...)

     Compute modified Bessel functions of the second kind.

     The order of the Bessel function ALPHA must be real.  The points
     for evaluation X may be complex.

     If the optional argument OPT is 1 or true, the result K is
     multiplied by ‘exp (X)’.

     If ALPHA is a scalar, the result is the same size as X.  If X is a
     scalar, the result is the same size as ALPHA.  If ALPHA is a row
     vector and X is a column vector, the result is a matrix with
     ‘length (X)’ rows and ‘length (ALPHA)’ columns.  Otherwise, ALPHA
     and X must conform and the result will be the same size.

     If requested, IERR contains the following status information and is
     the same size as the result.

       0. Normal return.

       1. Input error, return ‘NaN’.

       2. Overflow, return ‘Inf’.

       3. Loss of significance by argument reduction results in less
          than half of machine accuracy.

       4. Complete loss of significance by argument reduction, return
          ‘NaN’.

       5. Error—no computation, algorithm termination condition not met,
          return ‘NaN’.

     See also: Note: besseli, Note: besselj,
     Note: bessely, Note: besselh.

 -- : H = besselh (ALPHA, X)
 -- : H = besselh (ALPHA, K, X)
 -- : H = besselh (ALPHA, K, X, OPT)
 -- : [H, IERR] = besselh (...)
     Compute Bessel functions of the third kind (Hankel functions).

     The order of the Bessel function ALPHA must be real.  The kind of
     Hankel function is specified by K and may be either first (K = 1)
     or second (K = 2).  The default is Hankel functions of the first
     kind.  The points for evaluation X may be complex.

     If the optional argument OPT is 1 or true, the result is multiplied
     by ‘exp (-I*X)’ for K = 1 or ‘exp (I*X)’ for K = 2.

     If ALPHA is a scalar, the result is the same size as X.  If X is a
     scalar, the result is the same size as ALPHA.  If ALPHA is a row
     vector and X is a column vector, the result is a matrix with
     ‘length (X)’ rows and ‘length (ALPHA)’ columns.  Otherwise, ALPHA
     and X must conform and the result will be the same size.

     If requested, IERR contains the following status information and is
     the same size as the result.

       0. Normal return.

       1. Input error, return ‘NaN’.

       2. Overflow, return ‘Inf’.

       3. Loss of significance by argument reduction results in less
          than half of machine accuracy.

       4. Complete loss of significance by argument reduction, return
          ‘NaN’.

       5. Error—no computation, algorithm termination condition not met,
          return ‘NaN’.

     See also: Note: besselj, Note: bessely,
     Note: besseli, Note: besselk.

 -- : beta (A, B)
     Compute the Beta function for real inputs A and B.

     The Beta function definition is

          beta (a, b) = gamma (a) * gamma (b) / gamma (a + b).

     The Beta function can grow quite large and it is often more useful
     to work with the logarithm of the output rather than the function
     directly.  Note: betaln, for computing the logarithm of
     the Beta function in an efficient manner.

     See also: Note: betaln, Note: betainc,
     Note: betaincinv.

 -- : betainc (X, A, B)
 -- : betainc (X, A, B, TAIL)
     Compute the incomplete beta function.

     This is defined as

                         x
                        /
                       |
          I_x (a, b) = | t^(a-1) (1-t)^(b-1) dt
                       |
                       /
                      0

     with real X in the range [0,1].  The inputs A and B must be real
     and strictly positive (> 0).  If one of the inputs is not a scalar
     then the other inputs must be scalar or of compatible dimensions.

     By default, TAIL is "lower" and the incomplete beta function
     integrated from 0 to X is computed.  If TAIL is "upper" then the
     complementary function integrated from X to 1 is calculated.  The
     two choices are related by

     betainc (X, A, B, "upper") = 1 - betainc (X, A, B, "lower").

     ‘betainc’ uses a more sophisticated algorithm than subtraction to
     get numerically accurate results when the "lower" value is small.

     Reference: A. Cuyt, V. Brevik Petersen, B. Verdonk, H. Waadeland,
     W.B. Jones, ‘Handbook of Continued Fractions for Special
     Functions’, ch. 18.

     See also: Note: beta, Note: betaincinv,
     Note: betaln.

 -- : betaincinv (Y, A, B)
 -- : betaincinv (Y, A, B, "lower")
 -- : betaincinv (Y, A, B, "upper")
     Compute the inverse of the normalized incomplete beta function.

     The normalized incomplete beta function is defined as

                         x
                        /
                       |
          I_x (a, b) = | t^(a-1) (1-t)^(b-1) dt
                       |
                       /
                      0

     If two inputs are scalar, then ‘betaincinv (Y, A, B)’ is returned
     for each of the other inputs.

     If two or more inputs are not scalar, the sizes of them must agree,
     and ‘betaincinv’ is applied element-by-element.

     The variable Y must be in the interval [0,1], while A and B must be
     real and strictly positive.

     By default, TAIL is "lower" and the inverse of the incomplete beta
     function integrated from 0 to X is computed.  If TAIL is "upper"
     then the complementary function integrated from X to 1 is inverted.

     The function is computed by standard Newton’s method, by solving

          Y - betainc (X, A, B) = 0

     See also: Note: betainc, Note: beta, Note:
     betaln.

 -- : betaln (A, B)
     Compute the natural logarithm of the Beta function for real inputs
     A and B.

     ‘betaln’ is defined as

          betaln (a, b) = log (beta (a, b))

     and is calculated in a way to reduce the occurrence of underflow.

     The Beta function can grow quite large and it is often more useful
     to work with the logarithm of the output rather than the function
     directly.

     See also: Note: beta, Note: betainc, Note:
     betaincinv, Note: gammaln.

 -- : bincoeff (N, K)
     Return the binomial coefficient of N and K.

     The binomial coefficient is defined as

           /   \
           | n |    n (n-1) (n-2) ... (n-k+1)
           |   |  = -------------------------
           | k |               k!
           \   /

     For example:

          bincoeff (5, 2)
             ⇒ 10

     In most cases, the ‘nchoosek’ function is faster for small scalar
     integer arguments.  It also warns about loss of precision for big
     arguments.

     See also: Note: nchoosek.

 -- : commutation_matrix (M, N)
     Return the commutation matrix K(m,n) which is the unique M*N by M*N
     matrix such that K(m,n) * vec(A) = vec(A') for all m by n matrices
     A.

     If only one argument M is given, K(m,m) is returned.

     See Magnus and Neudecker (1988), ‘Matrix Differential Calculus with
     Applications in Statistics and Econometrics.’

 -- : cosint (X)
     Compute the cosine integral function:

                      +oo
                     /
          Ci (x) = - | (cos (t)) / t dt
                     /
                    x

     An equivalent definition is

                                       x
                                      /
                                      |  cos (t) - 1
          Ci (x) = gamma + log (x) +  | -------------  dt
                                      |        t
                                      /
                                     0

     Reference:

     M. Abramowitz and I.A. Stegun, ‘Handbook of Mathematical Functions’
     1964.

     See also: Note: sinint, Note: expint, Note:
     cos.

 -- : duplication_matrix (N)
     Return the duplication matrix Dn which is the unique n^2 by
     n*(n+1)/2 matrix such that Dn vech (A) = vec (A) for all symmetric
     n by n matrices A.

     See Magnus and Neudecker (1988), ‘Matrix Differential Calculus with
     Applications in Statistics and Econometrics.’

 -- : dawson (Z)
     Compute the Dawson (scaled imaginary error) function.

     The Dawson function is defined as

          (sqrt (pi) / 2) * exp (-z^2) * erfi (z)

     See also: Note: erfc, Note: erf, *note erfcx:
     XREFerfcx, Note: erfi, Note: erfinv, Note:
     erfcinv.

 -- : [SN, CN, DN, ERR] = ellipj (U, M)
 -- : [SN, CN, DN, ERR] = ellipj (U, M, TOL)
     Compute the Jacobi elliptic functions SN, CN, and DN of complex
     argument U and real parameter M.

     If M is a scalar, the results are the same size as U.  If U is a
     scalar, the results are the same size as M.  If U is a column
     vector and M is a row vector, the results are matrices with ‘length
     (U)’ rows and ‘length (M)’ columns.  Otherwise, U and M must
     conform in size and the results will be the same size as the
     inputs.

     The value of U may be complex.  The value of M must be 0 ≤ M ≤ 1.

     The optional input TOL is currently ignored (MATLAB uses this to
     allow faster, less accurate approximation).

     If requested, ERR contains the following status information and is
     the same size as the result.

       0. Normal return.

       1. Error—no computation, algorithm termination condition not met,
          return ‘NaN’.

     Reference: Milton Abramowitz and Irene A Stegun, ‘Handbook of
     Mathematical Functions’, Chapter 16 (Sections 16.4, 16.13, and
     16.15), Dover, 1965.

     See also: Note: ellipke.

 -- : K = ellipke (M)
 -- : K = ellipke (M, TOL)
 -- : [K, E] = ellipke (...)
     Compute complete elliptic integrals of the first K(M) and second
     E(M) kind.

     M must be a scalar or real array with -Inf ≤ M ≤ 1.

     The optional input TOL controls the stopping tolerance of the
     algorithm and defaults to ‘eps (class (M))’.  The tolerance can be
     increased to compute a faster, less accurate approximation.

     When called with one output only elliptic integrals of the first
     kind are returned.

     Mathematical Note:

     Elliptic integrals of the first kind are defined as

                   1
                  /               dt
          K (m) = | ------------------------------
                  / sqrt ((1 - t^2)*(1 - m*t^2))
                 0

     Elliptic integrals of the second kind are defined as

                   1
                  /  sqrt (1 - m*t^2)
          E (m) = |  ------------------ dt
                  /  sqrt (1 - t^2)
                 0

     Reference: Milton Abramowitz and Irene A. Stegun, ‘Handbook of
     Mathematical Functions’, Chapter 17, Dover, 1965.

     See also: Note: ellipj.

 -- : erf (Z)
     Compute the error function.

     The error function is defined as

                                  z
                        2        /
          erf (z) = --------- *  | e^(-t^2) dt
                    sqrt (pi)    /
                              t=0

     See also: Note: erfc, Note: erfcx, *note erfi:
     XREFerfi, Note: dawson, Note: erfinv, Note:
     erfcinv.

 -- : erfc (Z)
     Compute the complementary error function.

     The complementary error function is defined as ‘1 - erf (Z)’.

     See also: Note: erfcinv, Note: erfcx, Note:
     erfi, Note: dawson, Note: erf, Note:
     erfinv.

 -- : erfcx (Z)
     Compute the scaled complementary error function.

     The scaled complementary error function is defined as

          exp (z^2) * erfc (z)

     See also: Note: erfc, Note: erf, *note erfi:
     XREFerfi, Note: dawson, Note: erfinv, Note:
     erfcinv.

 -- : erfi (Z)
     Compute the imaginary error function.

     The imaginary error function is defined as

          -i * erf (i*z)

     See also: Note: erfc, Note: erf, *note erfcx:
     XREFerfcx, Note: dawson, Note: erfinv,
     Note: erfcinv.

 -- : erfinv (X)
     Compute the inverse error function.

     The inverse error function is defined such that

          erf (Y) == X

     See also: Note: erf, Note: erfc, *note erfcx:
     XREFerfcx, Note: erfi, Note: dawson, Note:
     erfcinv.

 -- : erfcinv (X)
     Compute the inverse complementary error function.

     The inverse complementary error function is defined such that

          erfc (Y) == X

     See also: Note: erfc, Note: erf, *note erfcx:
     XREFerfcx, Note: erfi, Note: dawson, Note:
     erfinv.

 -- : expint (X)
     Compute the exponential integral.

     The exponential integral is defined as:

                     +oo
                    /
                    | exp (-t)
          E_1 (x) = | -------- dt
                    |    t
                    /
                   x

     Note: For compatibility, this function uses the MATLAB definition
     of the exponential integral.  Most other sources refer to this
     particular value as E_1 (x), and the exponential integral as

                      +oo
                     /
                     | exp (-t)
          Ei (x) = - | -------- dt
                     |    t
                     /
                   -x

     The two definitions are related, for positive real values of X, by
     ‘E_1 (-x) = -Ei (x) - i*pi’.

     References:

     M. Abramowitz and I.A. Stegun, ‘Handbook of Mathematical
     Functions’, 1964.

     N. Bleistein and R.A. Handelsman, ‘Asymptotic expansions of
     integrals’, 1986.

     See also: Note: cosint, Note: sinint, Note:
     exp.

 -- : gamma (Z)
     Compute the Gamma function.

     The Gamma function is defined as

                       infinity
                      /
          gamma (z) = | t^(z-1) exp (-t) dt.
                      /
                   t=0

     Programming Note: The gamma function can grow quite large even for
     small input values.  In many cases it may be preferable to use the
     natural logarithm of the gamma function (‘gammaln’) in calculations
     to minimize loss of precision.  The final result is then ‘exp
     (RESULT_USING_GAMMALN).’

     See also: Note: gammainc, Note: gammaln,
     Note: factorial.

 -- : gammainc (X, A)
 -- : gammainc (X, A, TAIL)
     Compute the normalized incomplete gamma function.

     This is defined as

                                          x
                                 1       /
          gammainc (x, a) = ---------    | exp (-t) t^(a-1) dt
                            gamma (a)    /
                                      t=0

     with the limiting value of 1 as X approaches infinity.  The
     standard notation is P(a,x), e.g., Abramowitz and Stegun (6.5.1).

     If A is scalar, then ‘gammainc (X, A)’ is returned for each element
     of X and vice versa.

     If neither X nor A is scalar then the sizes of X and A must agree,
     and ‘gammainc’ is applied element-by-element.  The elements of A
     must be non-negative.

     By default, TAIL is "lower" and the incomplete gamma function
     integrated from 0 to X is computed.  If TAIL is "upper" then the
     complementary function integrated from X to infinity is calculated.

     If TAIL is "scaledlower", then the lower incomplete gamma function
     is multiplied by gamma(a+1)*exp(x)/(x^a).  If TAIL is
     "scaledupper", then the upper incomplete gamma function is
     multiplied by the same quantity.

     References:

     M. Abramowitz and I.A. Stegun, ‘Handbook of mathematical
     functions’, Dover publications, Inc., 1972.

     W. Gautschi, ‘A computational procedure for incomplete gamma
     functions’, ACM Trans.  Math Software, pp.  466–481, Vol 5, No.  4,
     2012.

     W. H. Press, S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery,
     ‘Numerical Recipes in Fortran 77’, ch. 6.2, Vol 1, 1992.

     See also: Note: gamma, *note gammaincinv:
     XREFgammaincinv, Note: gammaln.

 -- : gammaincinv (Y, A)
 -- : gammaincinv (Y, A, TAIL)
     Compute the inverse of the normalized incomplete gamma function.

     The normalized incomplete gamma function is defined as

                                          x
                                 1       /
          gammainc (x, a) = ---------    | exp (-t) t^(a-1) dt
                            gamma (a)    /
                                      t=0

     and ‘gammaincinv (gammainc (X, A), A) = X’ for each non-negative
     value of X.  If A is scalar then ‘gammaincinv (Y, A)’ is returned
     for each element of Y and vice versa.

     If neither Y nor A is scalar then the sizes of Y and A must agree,
     and ‘gammaincinv’ is applied element-by-element.  The variable Y
     must be in the interval [0,1] while A must be real and positive.

     By default, TAIL is "lower" and the inverse of the incomplete gamma
     function integrated from 0 to X is computed.  If TAIL is "upper",
     then the complementary function integrated from X to infinity is
     inverted.

     The function is computed with Newton’s method by solving

          Y - gammainc (X, A) = 0

     Reference: A. Gil, J. Segura, and N. M. Temme, ‘Efficient and
     accurate algorithms for the computation and inversion of the
     incomplete gamma function ratios’, SIAM J. Sci.  Computing, pp.
     A2965–A2981, Vol 34, 2012.

     See also: Note: gammainc, Note: gamma,
     Note: gammaln.

 -- : L = legendre (N, X)
 -- : L = legendre (N, X, NORMALIZATION)
     Compute the associated Legendre function of degree N and order M =
     0 ... N.

     The value N must be a real non-negative integer.

     X is a vector with real-valued elements in the range [-1, 1].

     The optional argument NORMALIZATION may be one of "unnorm", "sch",
     or "norm".  The default if no normalization is given is "unnorm".

     When the optional argument NORMALIZATION is "unnorm", compute the
     associated Legendre function of degree N and order M and return all
     values for M = 0 ... N.  The return value has one dimension more
     than X.

     The associated Legendre function of degree N and order M:

           m         m      2  m/2   d^m
          P(x) = (-1) * (1-x  )    * ----  P(x)
           n                         dx^m   n

     with Legendre polynomial of degree N:

                    1    d^n   2    n
          P(x) = ------ [----(x - 1) ]
           n     2^n n!  dx^n

     ‘legendre (3, [-1.0, -0.9, -0.8])’ returns the matrix:

           x  |   -1.0   |   -0.9   |   -0.8
          ------------------------------------
          m=0 | -1.00000 | -0.47250 | -0.08000
          m=1 |  0.00000 | -1.99420 | -1.98000
          m=2 |  0.00000 | -2.56500 | -4.32000
          m=3 |  0.00000 | -1.24229 | -3.24000

     When the optional argument NORMALIZATION is "sch", compute the
     Schmidt semi-normalized associated Legendre function.  The Schmidt
     semi-normalized associated Legendre function is related to the
     unnormalized Legendre functions by the following:

     For Legendre functions of degree N and order 0:

            0      0
          SP(x) = P(x)
            n      n

     For Legendre functions of degree n and order m:

            m      m         m    2(n-m)! 0.5
          SP(x) = P(x) * (-1)  * [-------]
            n      n              (n+m)!

     When the optional argument NORMALIZATION is "norm", compute the
     fully normalized associated Legendre function.  The fully
     normalized associated Legendre function is related to the
     unnormalized associated Legendre functions by the following:

     For Legendre functions of degree N and order M

            m      m         m    (n+0.5)(n-m)! 0.5
          NP(x) = P(x) * (-1)  * [-------------]
            n      n                  (n+m)!

 -- : gammaln (X)
 -- : lgamma (X)
     Return the natural logarithm of the gamma function of X.

     See also: Note: gamma, Note: gammainc.

 -- : psi (Z)
 -- : psi (K, Z)
     Compute the psi (polygamma) function.

     The polygamma functions are the Kth derivative of the logarithm of
     the gamma function.  If unspecified, K defaults to zero.  A value
     of zero computes the digamma function, a value of 1, the trigamma
     function, and so on.

     The digamma function is defined:

          psi (z) = d (log (gamma (z))) / dx

     When computing the digamma function (when K equals zero), Z can
     have any value real or complex value.  However, for polygamma
     functions (K higher than 0), Z must be real and non-negative.

     See also: Note: gamma, Note: gammainc,
     Note: gammaln.

 -- : sinint (X)
     Compute the sine integral function:

                     x
                    /
          Si (x) =  | sin (t) / t dt
                    /
                   0

     Reference: M. Abramowitz and I.A. Stegun, ‘Handbook of Mathematical
     Functions’, 1964.

     See also: Note: cosint, Note: expint, Note:
     sin.


automatically generated by info2www version 1.2.2.9