FunctionalSpecialFunction(R, F)ΒΆ
combfunc.spad line 443 [edit on github]
- R: Join(Comparable, IntegralDomain) 
- F: FunctionSpace R 
Provides some special functions over an integral domain.
- abs: F -> F
- abs(f)returns the absolute value operator applied to- f.
- airyAi: F -> F
- airyAi(x)returns the Airy- Aifunction applied to- x.
- airyAiPrime: F -> F
- airyAiPrime(x)returns the derivative of Airy- Aifunction applied to- x.
- airyBi: F -> F
- airyBi(x)returns the Airy- Bifunction applied to- x.
- airyBiPrime: F -> F
- airyBiPrime(x)returns the derivative of Airy- Bifunction applied to- x.
- angerJ: (F, F) -> F
- angerJ(v, z)is the Anger- Jfunction.
- belong?: BasicOperator -> Boolean
- belong?(op)returns- trueif- opis a special function operator.
- besselI: (F, F) -> F
- besselI(x, y)returns the Bessel- Ifunction applied to- xand- y.
- besselJ: (F, F) -> F
- besselJ(x, y)returns the Bessel- Jfunction applied to- xand- y.
- besselK: (F, F) -> F
- besselK(x, y)returns the Bessel- Kfunction applied to- xand- y.
- besselY: (F, F) -> F
- besselY(x, y)returns the Bessel- Yfunction applied to- xand- y.
- Beta: (F, F) -> F
- Beta(x, y)returns the beta function applied to- xand- y.
- Beta: (F, F, F) -> F
- Beta(x, a, b)is incomplete Beta function applied to- x, a and- b.
- ceiling: F -> F
- ceiling(x)returns the smallest integer above or equal- x.
- charlierC: (F, F, F) -> F
- charlierC(n, a, z)is the Charlier polynomial.
- conjugate: F -> F
- conjugate(f)returns the conjugate value operator applied to- f.
- digamma: F -> F
- digamma(x)returns the digamma function applied to- x.
- diracDelta: F -> F
- diracDelta(x)is unit mass at zeros of- x.
- ellipticE: (F, F) -> F
- ellipticE(z, m)is the incomplete elliptic integral of the second kind:- ellipticE(z, m) = integrate(sqrt(1-m*t^2)/sqrt(1-t^2), t = 0..z).
- ellipticE: F -> F
- ellipticE(m)is the complete elliptic integral of the second kind:- ellipticE(m) = integrate(sqrt(1-m*t^2)/sqrt(1-t^2), t = 0..1).
- ellipticF: (F, F) -> F
- ellipticF(z, m)is the incomplete elliptic integral of the first kind :- ellipticF(z, m) = integrate(1/sqrt((1-t^2)*(1-m*t^2)), t = 0..z).
- ellipticK: F -> F
- ellipticK(m)is the complete elliptic integral of the first kind:- ellipticK(m) = integrate(1/sqrt((1-t^2)*(1-m*t^2)), t = 0..1).
- ellipticPi: (F, F, F) -> F
- ellipticPi(z, n, m)is the incomplete elliptic integral of the third kind:- ellipticPi(z, n, m) = integrate(1/((1-n*t^2)*sqrt((1-t^2)*(1-m*t^2))), t = 0..z).
- floor: F -> F
- floor(x)returns the largest integer below or equal- x.
- fractionPart: F -> F
- fractionPart(x)returns the fractional part of- x.
- Gamma: (F, F) -> F
- Gamma(a, x)returns the incomplete Gamma function applied to a and- x.
- Gamma: F -> F
- Gamma(f)returns the formal Gamma function applied to- f.
- hahn_p: (F, F, F, F, F) -> F
- hahn_p(n, a, b, bar_a, bar_b, z)is the continuous Hahn polynomial.
- hahnQ: (F, F, F, F, F) -> F
- hahnQ(n, a, b, N, z)is the Hahn polynomial.
- hahnR: (F, F, F, F, F) -> F
- hahnR(n, c, d, N, z)is the dual Hahn polynomial.
- hahnS: (F, F, F, F, F) -> F
- hahnS(n, a, b, c, z)is the continuous dual Hahn polynomial.
- hankelH1: (F, F) -> F
- hankelH1(v, z)is first Hankel function (Bessel function of the third kind).
- hankelH2: (F, F) -> F
- hankelH2(v, z)is the second Hankel function (Bessel function of the third kind).
- hermiteH: (F, F) -> F
- hermiteH(n, z)is the Hermite polynomial.
- hypergeometricF: (List F, List F, F) -> F
- hypergeometricF(la, lb, z)is the generalized hypergeometric function.
- iAiryAi: F -> F
- iAiryAi(x)should be local but conditional.
- iAiryAiPrime: F -> F
- iAiryAiPrime(x)should be local but conditional.
- iAiryBi: F -> F
- iAiryBi(x)should be local but conditional.
- iAiryBiPrime: F -> F
- iAiryBiPrime(x)should be local but conditional.
- iiabs: F -> F
- iiabs(x)should be local but conditional.
- iiAiryAi: F -> F
- iiAiryAi(x)should be local but conditional.
- iiAiryAiPrime: F -> F
- iiAiryAiPrime(x)should be local but conditional.
- iiAiryBi: F -> F
- iiAiryBi(x)should be local but conditional.
- iiAiryBiPrime: F -> F
- iiAiryBiPrime(x)should be local but conditional.
- iiBesselI: List F -> F
- iiBesselI(x)should be local but conditional.
- iiBesselJ: List F -> F
- iiBesselJ(x)should be local but conditional.
- iiBesselK: List F -> F
- iiBesselK(x)should be local but conditional.
- iiBesselY: List F -> F
- iiBesselY(x)should be local but conditional.
- iiBeta: List F -> F
- iiBeta(x)should be local but conditional.
- iiconjugate: F -> F
- iiconjugate(x)should be local but conditional.
- iidigamma: F -> F
- iidigamma(x)should be local but conditional.
- iiGamma: F -> F
- iiGamma(x)should be local but conditional.
- iiHypergeometricF: List F -> F
- iiHypergeometricF(l)should be local but conditional.
- iipolygamma: List F -> F
- iipolygamma(x)should be local but conditional.
- iiPolylog: (F, F) -> F
- iiPolylog(x, s)should be local but conditional.
- iLambertW: F -> F
- iLambertW(x)should be local but conditional.
- jacobiCn: (F, F) -> F
- jacobiCn(z, m)is the Jacobi elliptic- cnfunction, defined by- jacobiCn(z, m)^2 + jacobiSn(z, m)^2 = 1and- jacobiCn(0, m) = 1.
- jacobiDn: (F, F) -> F
- jacobiDn(z, m)is the Jacobi elliptic- dnfunction, defined by- jacobiDn(z, m)^2 + m*jacobiSn(z, m)^2 = 1and- jacobiDn(0, m) = 1.
- jacobiP: (F, F, F, F) -> F
- jacobiP(n, a, b, z)is the Jacobi polynomial.
- jacobiSn: (F, F) -> F
- jacobiSn(z, m)is the Jacobi elliptic- snfunction, defined by the formula- jacobiSn(ellipticF(z, m), m) = z.
- jacobiTheta: (F, F) -> F
- jacobiTheta(z, m)is the Jacobi Theta function in Jacobi notation.
- jacobiZeta: (F, F) -> F
- jacobiZeta(z, m)is the Jacobi elliptic zeta function, defined by- D(jacobiZeta(z, m), z) = jacobiDn(z, m)^2 - ellipticE(m)/ellipticK(m)and- jacobiZeta(0, m) = 0.
- kelvinBei: (F, F) -> F
- kelvinBei(v, z)is the Kelvin bei function defined by equality.- kelvinBei(v, z) = imag(besselJ(v, exp(3*\%pi*\%i/4)*z)). for- zand- vreal.
- kelvinBer: (F, F) -> F
- kelvinBer(v, z)is the Kelvin ber function defined by equality- kelvinBer(v, z) = real(besselJ(v, exp(3*\%pi*\%i/4)*z))for- zand- vreal.
- kelvinKei: (F, F) -> F
- kelvinKei(v, z)is the Kelvin kei function defined by equality- kelvinKei(v, z) = imag(exp(-v*\%pi*\%i/2)*besselK(v, exp(\%pi*\%i/4)*z))for- zand- vreal.
- kelvinKer: (F, F) -> F
- kelvinKer(v, z)is the Kelvin kei function defined by equality- kelvinKer(v, z) = real(exp(-v*\%pi*\%i/2)*besselK(v, exp(\%pi*\%i/4)*z))for- zand- vreal.
- krawtchoukK: (F, F, F, F) -> F
- krawtchoukK(n, p, N, z)is the Krawtchouk polynomial.
- kummerM: (F, F, F) -> F
- kummerM(a, b, z)is the Kummer- Mfunction.
- kummerU: (F, F, F) -> F
- kummerU(a, b, z)is the Kummer- Ufunction.
- laguerreL: (F, F, F) -> F
- laguerreL(n, a, z)is the Laguerre polynomial.
- lambertW: F -> F
- lambertW(x)is the Lambert- Wfunction at- x.
- legendreP: (F, F, F) -> F
- legendreP(nu, mu, z)is the Legendre- Pfunction.
- legendreQ: (F, F, F) -> F
- legendreQ(nu, mu, z)is the Legendre- Qfunction.
- lerchPhi: (F, F, F) -> F
- lerchPhi(z, s, a)is the Lerch Phi function.
- lommelS1: (F, F, F) -> F
- lommelS1(mu, nu, z)is the Lommel- sfunction.
- lommelS2: (F, F, F) -> F
- lommelS2(mu, nu, z)is the Lommel- Sfunction.
- meijerG: (List F, List F, List F, List F, F) -> F
- meijerG(la, lb, lc, ld, z)is the meijerG function.
- meixnerM: (F, F, F, F) -> F
- meixnerM(n, b, c, z)is the Meixner polynomial.
- meixnerP: (F, F, F, F) -> F
- meixnerP(n, phi, lambda, z)is the Meixner-Pollaczek polynomial.
- operator: BasicOperator -> BasicOperator
- operator(op)returns a copy of- opwith the domain-dependent properties appropriate for- F; error if- opis not a special function operator.
- polygamma: (F, F) -> F
- polygamma(x, y)returns the polygamma function applied to- xand- y.
- polylog: (F, F) -> F
- polylog(s, x)is the polylogarithm of order- sat- x.
- racahR: (F, F, F, F, F, F) -> F
- racahR(n, a, b, c, d, z)is the Racah polynomial.
- riemannZeta: F -> F
- riemannZeta(z)is the Riemann Zeta function.
- sign: F -> F
- sign(x)returns the sign of- x.
- struveH: (F, F) -> F
- struveH(v, z)is the Struve- Hfunction.
- struveL: (F, F) -> F
- struveL(v, z)is the Struve- Lfunction defined by the formula- struveL(v, z) = -\%i^exp(-v*\%pi*\%i/2)*struveH(v, \%i*z).
- unitStep: F -> F
- unitStep(x)is 0 for- xless than 0, 1 for- xbigger or equal 0.
- weberE: (F, F) -> F
- weberE(v, z)is the Weber- Efunction.
- weierstrassP: (F, F, F) -> F
- weierstrassP(g2, g3, x)is the Weierstrass- Pfunction.
- weierstrassPInverse: (F, F, F) -> F
- weierstrassPInverse(g2, g3, z)is the inverse of Weierstrass- Pfunction, defined by the formula- weierstrassP(g2, g3, weierstrassPInverse(g2, g3, z)) = z.
- weierstrassPPrime: (F, F, F) -> F
- weierstrassPPrime(g2, g3, x)is the derivative of Weierstrass- Pfunction.
- weierstrassSigma: (F, F, F) -> F
- weierstrassSigma(g2, g3, x)is the Weierstrass Sigma function.
- weierstrassZeta: (F, F, F) -> F
- weierstrassZeta(g2, g3, x)is the Weierstrass Zeta function.
- whittakerM: (F, F, F) -> F
- whittakerM(k, m, z)is the Whittaker- Mfunction.
- whittakerW: (F, F, F) -> F
- whittakerW(k, m, z)is the Whittaker- Wfunction.
- wilsonW: (F, F, F, F, F, F) -> F
- wilsonW(n, a, b, c, d, z)is the Wilson polynomial.