dstable | R Documentation |
These functions provide information about the stable distribution
with the location, the dispersion, the skewness and the tail thickness
respectively modelled by the parameters loc
, disp
,
skew
and tail
. These differ from those of 'stabledist' (i.e., Nolan's 0-parameterization and
1-parameterization as detailed in Nolan (2020)). See the README for how to make equivalent calls
to those of 'stabledist' (i.e., Nolan's 0-parameterization and
1-parameterization as detailed in Nolan (2020)).
dstable( x, loc = 0, disp = 1/sqrt(2), skew = 0, tail = 2, npt = 501, up = 10, eps = 1e-06, integration = "Romberg" ) pstable(q, loc = 0, disp = 1/sqrt(2), skew = 0, tail = 2, eps = 1e-06) qstable(p, loc = 0, disp = 1/sqrt(2), skew = 0, tail = 2, eps = 1e-06) rstable(n = 1, loc = 0, disp = 1/sqrt(2), skew = 0, tail = 2, eps = 1e-06) hstable(x, loc = 0, disp = 1/sqrt(2), skew = 0, tail = 2, eps = 1e-06)
x, q |
vector of quantiles. |
loc |
vector of (real) location parameters. |
disp |
vector of (positive) dispersion parameters. |
skew |
vector of skewness parameters (in [-1,1]). |
tail |
vector of parameters (in [0,2]) related to the tail thickness. |
npt, up, integration |
As detailed herein – only available when using |
eps |
scalar giving the required precision in computation. |
p |
vector of probabilites. |
n |
number of observations. |
dstable
, pstable
, qstable
and hstable
compute the density, the distribution, the quantile and the hazard functions
of a stable variate. rstable
generates random deviates with
the prescribed stable distribution.
loc
is a location parameter in the same way as the mean
in the normal distribution: it can take any real value.
disp
is a dispersion parameter in the same way as the standard
deviation in the normal distribution: it can take any positive value.
skew
is a skewness parameter: it can take any value in (-1,1).
The distribution is right-skewed, symmetric and left-skewed
when skew
is negative, null or positive respectively.
tail
is a tail parameter (often named the characteristic exponent):
it can take any value in (0,2) (with tail=1
and tail=2
yielding the Cauchy and the normal distributions respectively
when symmetry holds).
If loc
, disp
, skew
, or tail
are not
specified they assume the default values of 0, 1/sqrt(2),
0 and 2 respectively. This corresponds to a normal
variate with mean=0 and variance=1/2 disp^2.
The stable characteristic function is given by
phi(t) = i loc t - disp |t|^tail [1+i skew sign(t) omega(t,tail)]
where
omega(t,tail) = (2/pi) log|t|
when tail=1
, and
omega(t,tail) = tan(pi alpha / 2)
otherwise.
The characteristic function is inverted using Fourier's transform
to obtain the corresponding stable density. This inversion requires the
numerical evaluation of an integral from 0 to infinity.
Two algorithms
are proposed for this. The default is Romberg's method
(integration
="Romberg") which is used to evaluate the integral
with an error bounded by eps
.
The alternative method is Simpson's integration
(integration
="Simpson"): it approximates the
integral from 0 to infinity by an integral
from 0 to up
with npt
points subdividing (O, up).
These three extra arguments – integration
, up
and
npt
– are only available when using dstable
.
The other functions are all based on Romberg's algorithm.
[Swihart 2022 update:]
See the README for how to make equivalent calls to those of 'stabledist' (i.e., Nolan's 0-parameterization and 1-parameterization as detailed in Nolan (2020)). See github for Lambert and Lindsey 1999 JRSS-C journal article, which details the parameterization of the Buck (1995) stable distribution which allowed a Fourier inversion to arrive at a form similar to but not exactly the $g_d$ function as detailed in Nolan (2020), Abdul-Hamid and Nolan (1998) and Nolan (1997). The Nolan (2020) reference is a textbook that provides an accessible and comprehensive summary of stable distributions in the 25 years or so since the core of this R package was made and put on CRAN. The Buck (1995) parameterization most closely resembles the Zolotarev B parameterization outlined in Definition 3.6 on page 93 of Nolan (2020) – except that Buck (1995) did not allow the scale parameter to multiply with the location parameter. This explains why the 'Zolotarev B' entry in Table 3.1 on page 97 of Nolan (2020) has the location parameter being multiplied by the scale parameter whereas in converting the Lindsey and Lambert (1999) to Nolan 1-parameterization the location parameter stays the same.
To be clear, stable::dstable
and stable::pstable
are evaluated
by numerically integrating the inverse Fourier transform. The code works
reasonably for small and moderate values of x, but will have numerical issues
in some cases (such as values from stable::pstable
being greater than
1 or or not being monotonic). The arguments npt
, up
,
integration
, and eps
can be adjusted
to improve accuracy at the cost of speed, but will still have limitations.
Functions that better handle these problems are available in other packages
(such as stabledist
and stable
) that use an alternative method
(as detailed in Nolan 1997)
distinct from directly numerically integrating the Fourier inverse transform.
See last example in the README.
dstable
: density
pstable
: cdf
qstable
: quantiles
rstable
: random deviates
hstable
: hazard
Philippe Lambert (Catholic University of Louvain, Belgium, phlambert@stat.ucl.ac.be)
Jim Lindsey
Lambert, P. and Lindsey, J.K. (1999) Analysing financial returns using regression models based on non-symmetric stable distributions. Applied Statistics, 48, 409-424.
Nolan, John P. Univariate stable distributions. Berlin/Heidelberg, Germany: Springer, 2020.
Nolan, John P. "Numerical calculation of stable densities and distribution functions." Communications in statistics. Stochastic models 13.4 (1997): 759-774.
Abdul-Hamid, Husein, and John P. Nolan. "Multivariate stable densities as functions of one dimensional projections." Journal of multivariate analysis 67.1 (1998): 80-89.
stablereg
to fit generalized nonlinear regression models
for the stable distribution parameters.
R packages stabledist
and libstableR
provide [dpqr] functions.
par(mfrow=c(2,2)) x <- seq(-5,5,by=0.1) # Influence of loc (location) plot(x,dstable(x,loc=-2,disp=1/sqrt(2),skew=-0.8,tail=1.5), type="l",ylab="",main="Varying LOCation") lines(x,dstable(x,loc=0,disp=1/sqrt(2),skew=-0.8,tail=1.5)) lines(x,dstable(x,loc=2,disp=1/sqrt(2),skew=-0.8,tail=1.5)) # Influence of disp (dispersion) plot(x,dstable(x,loc=0,disp=0.5,skew=0,tail=1.5), type="l",ylab="",main="Varying DISPersion") lines(x,dstable(x,loc=0,disp=1/sqrt(2),skew=0,tail=1.5)) lines(x,dstable(x,loc=0,disp=0.9,skew=0,tail=1.5)) # Influence of skew (skewness) plot(x,dstable(x,loc=0,disp=1/sqrt(2),skew=-0.8,tail=1.5), type="l",ylab="",main="Varying SKEWness") lines(x,dstable(x,loc=0,disp=1/sqrt(2),skew=0,tail=1.5)) lines(x,dstable(x,loc=0,disp=1/sqrt(2),skew=0.8,tail=1.5)) # Influence of tail (tail) plot(x,dstable(x,loc=0,disp=1/sqrt(2),skew=0,tail=0.8), type="l",ylab="",main="Varying TAIL thickness") lines(x,dstable(x,loc=0,disp=1/sqrt(2),skew=0,tail=1.5)) lines(x,dstable(x,loc=0,disp=1/sqrt(2),skew=0,tail=2))
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.