cf2DistFFT: Evaluating CDF/PDF/QF (quantiles) from the characteristic...

View source: R/cf2DistFFT.R

cf2DistFFTR Documentation

Evaluating CDF/PDF/QF (quantiles) from the characteristic function CF of a (continuous) distribution F by using the Fast Fourier Transform (FFT) algorithm

Description

TEST VERSION ! cf2DistFFT(cf, x, prob, options) evaluates the approximate values CDF(x), PDF(x), and/or the quantiles QF(prob) for given x and prob, by interpolation from the PDF-estimate computed by the numerical inversion of the given characteristic function CF by using the FFT algorithm.

Usage

cf2DistFFT(cf, x, prob, options)

Arguments

cf

function handle for the characteristic function CF.

x

vector of values from domain of the distribution F, if x is left empty, cf2DistFFT automatically selects vector x from the domain.

prob

vector of values from [0,1] for which the quantiles will be estimated, if prob is left empty, cf2DistFFT automatically selects vector prob = c(0.9, 0.95, 0.975, 0.99, 0.995, 0.999).

options

structure (list) with the following default parameters:

  • options$isCompound = FALSE treat the compound distributions, of the RV Y = X_1 + ... + X_N, where N is discrete RV and X\ge 0 are iid RVs from nonnegative continuous distribution,

  • options$N = 2^10 N points used by FFT,

  • options$xMin = -Inf set the lower limit of x,

  • options$xMax = Inf set the upper limit of x,

  • options$xMean = NULL set the MEAN value of x,

  • options$xStd = NULL set the STD value of x,

  • options$dt = NULL set grid step dt = 2*\pi/xRange,

  • options$T = NULL set upper limit of (0,T), T = N*dt,

  • options$SixSigmaRule = 6 set the rule for computing domain,

  • options$tolDiff = 1e-4 tol for numerical differentiation,

  • options$isPlot = TRUE plot the graphs of PDF/CDF,

  • options$DIST list with information for future evaluations, options$DIST is created automatically after first call:

    • options$DIST$xMin = xMin the lower limit of x,

    • options$DIST$xMax = XMax the upper limit of x,

    • options$DIST$xMean = xMean the MEAN value of x,

    • options$DIST$cft = cft CF evaluated at t_j : cf(t_j).

Details

The outputs of the algorithm cf2DistFFT are approximate values! The precission of the presented results depends on several different factors:
- application of the FFT algorithm for numerical inversion of the CF
- selected number of points used by the FFT algorithm (by default options$N = 2^10),
- estimated/calculated domain [A,B] of the distribution F,
used with the FFT algorithm. Optimally, [A,B] covers large part of the distribution domain, say more than 99\%. However, the default automatic procedure for selection of the domain [A,B] may fail. It is based on the 'SixSigmaRule': A = MEAN - SixSigmaRule * STD, and B = MEAN + SixSigmaRule * STD. Alternatively, change the options$SixSigmaRule to different value, say 12, or use the options$xMin and options$xMax to set manually the values of A and B.

Value

result

structure (list) with with the following results values:

result$x = x;
result$cdf = cdf;
result$pdf = pdf;
result$prob = prob;
result$qf = qf;
result$xFTT = xFFT;
result$pdfFFT = pdfFFT;
result$cdfFFT = cdfFFT;
result$SixSigmaRule = options.SixSigmaRule;
result$N = N;
result$dt = dt;
result$T = t[length(t)];
result$PrecisionCrit = PrecisionCrit;
result$myPrecisionCrit = options.crit;
result$isPrecisionOK = isPrecisionOK;
result$xMean = xMean;
result$xStd = xStd;
result$xMin = xMin;
result$xMax = xMax;
result$cf = cf;
result$options = options;
result$tictoc = toc.

References

[1] WITKOVSKY, V.: On the exact computation of the density and of the quantiles of linear combinations of t and F random variables. Journal of Statistical Planning and Inference 94 (2001), 1-13.

[2] WITKOVSKY, V.: Matlab algorithm TDIST: The distribution of a linear combination of Student's t random variables. In COMPSTAT 2004 Symposium (2004), J. Antoch, Ed., Physica-Verlag/Springer 2004, Heidelberg, Germany, pp. 1995-2002.

[3] WITKOVSKY, V., WIMMER, G., DUBY, T.: Logarithmic Lambert W x F random variables for the family of chi-squared distributions and their applications. Statistics & Probability Letters 96 (2015), 223-231.

[4] WITKOVSKY, V. (2016): Numerical inversion of a chracteristic function: An alternative tool to form the probability distribution of output quantity in linear measurement models. Acta IMEKO, 5(3), 32-34.

[5] WITKOVSKY, V., WIMMER G., DUBY T. (2016): Computing the aggregate loss distribution based on numerical inversion of the compound empirical characteristic function of frequency and severity. Preprint submitted to Insurance: Mathematics and Economics.

[6] DUBY, T., WIMMER, G., WITKOVSKY, V. (2016): MATLAB toolbox CRM for computing distributions of collective risk models. Preprint submitted to Journal of Statistical Software.

See Also

For more details see: https://arxiv.org/pdf/1701.08299.pdf.

Other CF Inversion Algorithm: cf2DistGP(), cf2PMF_FFT()

Examples

## EXAMPLE 1
## DISTRIBUTION OF A LINEAR COMBINATION OF THE INDEPENDENT RVs
## (Normal, Student's t, Rectangular, Triangular & Arcsine distribution)
## Y = X_{N} + X_{t} + 5*X_{R} + X_{T} + 10*X_{U}
## CFs: Normal, Student's t, Rectangular, Triangular, and Arcsine
cf_N  <- function(t) exp(-t^2/2)
cf_t <- function(t, nu) {pmin(1, Bessel::BesselK(abs(t) * sqrt(nu), nu / 2, TRUE) *
                                      exp(-abs(t) * sqrt(nu)) *
                                      (sqrt(nu) * abs(t))^(nu / 2) / 2^(nu / 2 - 1)/ gamma(nu / 2))
                        }
#cf_t2 <- function(t) {min(1, Bessel::BesselK(abs(t) * sqrt(1), 1 / 2, TRUE) *
                                                                 #exp(-abs(t) * sqrt(1)) *
                                                                   #   (sqrt(1) * abs(t))^(1 / 2) / 2^(1 / 2 - 1)/ gamma(1 / 2))
#}
cf_R <- function(t) pmin(1, sin(t) / t)
cf_T <- function(t) pmin(1, (2 - 2 * cos(t)) / t^2)
cf_U <- function(t) Bessel::BesselJ(t, 0)
## Characteristic function of the linear combination Y
c <- c(1, 1, 5, 1, 10)
nu <- 1
cf_Y <- function(t) {cf_N(c[1] * t) * cf_t(c[2] * t, nu) * cf_R(c[3] * t) *
                cf_T(c[4] * t) * cf_U(c[5] * t)}
options <- list()
options$N <- 2^10
options$xMin <- -50
options$xMax <- 50
result <- cf2DistFFT(cf = cf_Y, options = options)
# title('CDF of Y = X_{N} + X_{t} + 5*X_{R} + X_{T} + 10*X_{U}')
# problems with BesselK() function...

## EXAMPLE 2
## DISTRIBUTION OF A LINEAR COMBINATION OF THE INDEPENDENT CHI2 RVs
## (Chi-squared RVs with 1 and 10 degrees of freedom)
## Y = 10*X_{Chi2_1} + X_{Chi2_10}
## Characteristic functions of X_{Chi2_1} and X_{Chi2_10}
df1 <- 1
df2 <- 10
cfChi2_1 <- function(t) (1 - 2i * t)^(-df1 / 2)
cfChi2_10 <- function(t) (1 - 2i * t)^(-df2 / 2)
cf_Y <- function(t) cfChi2_1(10 * t) * cfChi2_10(t)
options <- list()
options$xMin <- 0
result <- cf2DistFFT(cf = cf_Y, options = options)
# title('CDF of Y = 10*X_{\chi^2_{1}} + X_{\chi^2_{10}}')

## EXAMPLE3 (PDF/CDF of the compound Poisson-Exponential distribution)
lambda1 <- 10
lambda2 <- 5
cfX <- function(t) cfX_Exponential(t, lambda2)
cf <- function(t) cfN_Poisson(t, lambda1, cfX)
x <- seq(from = 0, to = 8, length.out = 101)
prob <- c(0.9, 0.95, 0.99)
options <- list()
options$isCompound <- 1
result <- cf2DistFFT(cf, x, prob, options)

gajdosandrej/CharFunToolR documentation built on June 3, 2024, 7:46 p.m.