cf2DistGP | R Documentation |
cf2DistGP(cf, x, prob, options)
calcuates the CDF/PDF/QF from the Characteristic Function CF
by using the Gil-Pelaez inversion formulae:
cdf(x) = 1/2 + (1/\pi) * Integral_0^inf imag(exp(-1i*t*x)*cf(t)/t)*dt,
pdf(x) = (1/\pi) * Integral_0^inf real(exp(-1i*t*x)*cf(t))*dt.
The FOURIER INTEGRALs are calculated by using the simple TRAPEZOIDAL QUADRATURE method, see below.
cf2DistGP(cf, x, prob, options)
cf |
function handle for the characteristic function CF. |
x |
vector of values where the CDF/PDF is computed. |
prob |
vector of values from |
options |
structure (list) with the following default parameters:
|
If options.DIST
is provided, then cf2DistGP
evaluates CDF/PDF based on
this information only (it is useful, e.g., for subsequent evaluation of
the quantiles). options.DIST
is created automatically after first call.
This is supposed to be much faster, bacause there is no need for further
evaluations of the characteristic function. In fact, in such case the
function handle of the CF is not required, i.e. in such case do not specify cf
.
The required integrals are evaluated approximately by using the simple
trapezoidal rule on the interval (0,T)
, where T = N * dt
is a sufficienly
large integration upper limit in the frequency domain.
If the optimum values of N
and T
are unknown, we suggest, as a simple
rule of thumb, to start with the application of the six-sigma-rule for
determining the value of dt = (2*\pi)/(xMax-xMin)
, where xMax = xMean +
6*xStd
, and xMin = xMean - 6*xStd
, see [1]
.
Please note that THIS (TRAPEZOIDAL) METHOD IS AN APPROXIMATE METHOD: Frequently, with relatively low numerical precision of the results of the calculated PDF/CDF/QF, but frequently more efficient and more precise than comparable Monte Carlo methods.
However, the numerical error (truncation error and/or the integration error) could be and should be properly controled!
CONTROLING THE PRECISION:
Simple criterion for controling numerical precision is as follows: Set N
and T = N*dt
such that the value of the integrand function
Imag(e^(-1i*t*x) * cf(t)/t)
is sufficiently small for all t > T
, i.e.
PrecisionCrit = abs(cf(t)/t) <= tol
,
for pre-selected small tolerance value, say tol = 10^-8
. If this
criterion is not valid, the numerical precission of the result is
violated, and the method should be improved (e.g. by selecting larger N
or considering other more sofisticated algorithm - not considered here).
result |
structure (list) with CDF/PDF/QF amd further details. |
Ver.: 16-Sep-2018 17:59:07 (consistent with Matlab CharFunTool v1.3.0, 22-Sep-2017 11:11:11).
[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, 2001, 94, 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, 2015, 96, 223-231.
[4] WITKOVSKY, V.: Numerical inversion of a characteristic function: An alternative tool to form the probability distribution of output quantity in linear measurement models. Acta IMEKO, 2016, 5(3), 32-44.
[5] WITKOVSKY, V., WIMMER, G., DUBY, T. Computing the aggregate loss distribution based on numerical inversion of the compound empirical characteristic function of frequency and severity. ArXiv preprint, 2017, arXiv:1701.08299.
For more details see: https://arxiv.org/pdf/1701.08299.pdf.
Other CF Inversion Algorithm:
cf2DistFFT()
,
cf2PMF_FFT()
## EXAMPLE1
# Calculate CDF/PDF of N(0,1) by inverting its CF
cf <- function(t)
exp(-t ^ 2 / 2)
result <- cf2DistGP(cf)
## EXAMPLE2
# PDF/CDF of the compound Binomial-Exponential distribution
n <- 25
p <- 0.3
lambda <- 5
cfX <- function(t)
cfX_Exponential(t, lambda)
cf <- function(t)
cfN_Binomial(t, n, p, cfX)
x <- seq(0, 5, length.out = 101)
prob <- c(0.9, 0.95, 0.99)
options <- list()
options$isCompound <- TRUE
result <- cf2DistGP(cf, x, prob, options)
## 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(0, 8, length.out = 101)
prob <- c(0.9, 0.95, 0.99)
options <- list()
options$isCompound <- TRUE
result <- cf2DistGP(cf, x, prob, options)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.