GQD.density: Generate the Transition Density of a Scalar Generalized...

Description Usage Arguments Details Value Warning Author(s) References See Also Examples

Description

GQD.density approximates the transition density of a scalar generalized quadratic diffusion model (GQD). Given an initial value for the diffusion, Xs, the approximation is evaluated for all Xt at equispaced time-nodes given by splitting [s, t] into subintervals of length delt. GQD.density() approximates transitional densities of diffusions of the form:

ScalarEqn1.png

where

ScalarEqn2.png

and

ScalarEqn3.png

Usage

1
2
3
GQD.density(Xs, Xt, s, t, delt=1/100, Dtype='Saddlepoint', Trunc=c(4,4),
            P=100, alpha=0, lower=0, upper=50, print.output=TRUE,
            eval.density=TRUE)

Arguments

Xs

Initial value of the process at time s.

Xt

Vector of values at which the transition density is to be evaluated over the trajectory of the transition density from time s to t.

s

The starting time of the process.

t

The time horizon up to and including which the transitional density is evaluated.

delt

Size of the time increments at which successive evaluations are made.

Dtype

Character string indicating the type of density approximation (see details) to use. Types: 'Saddlepoint', 'Normal', 'Gamma', 'InvGamma' and 'Beta' are supported (default = 'Saddlepoint').

Trunc

Vector of length 2 containing the cumulant truncation order and the density truncation order respectively. May take on values 4, 6 and 8 with the constraint that Trunc[1] >= Trunc[2]. Default is c(4,4).

P

Normalization parameter indicating the number of points to use when normalizing members of the Pearson system (see details)

alpha

Normalization parameter controlig the mesh concentration when normalizing members of the Pearson system (see details). Increasing alpha decreases concentration around the mean and vice versa (default alpha = 0).

lower,upper

Lower and upper bounds for the normalization range.

print.output

If TRUE information about the model and algorithm is printed to the console.

eval.density

If TRUE, the density is evaluated in addition to calculating the moment eqns.

Details

GQD

GQD.density constructs an approximate transition density for a class of quadratic diffusion models. This is done by first evaluating the trajectory of the cumulants/moments of the diffusion numerically as the solution of a system of ordinary differential equations over a time horizon [s,t] split into equi-distant points delt units apart. Subsequently, the resulting cumulants/moments are carried into a density approximant (by default, a saddlepoint approximation) in order to evaluate the transtion surface.

Value

density

A matrix giving the density over the spatio-temporal mesh whose vertices are defined by paired permutations of the elements of X_t and time

Xt

A vector of points defining the state space at which the density was evaluated(recycled from input).

time

A vector of time points at which the density was evaluated.

cumulants

A matrix giving the cumulants of the diffusion. Row i gives the i-th cumulant.

moments

A matrix giving the moments of the diffusion. Row i gives the i-th cumulant.

mesh

A matrix giving the mesh used for normalization of the density.

Warning

Warning [1]: The system of ODEs that dictate the evolution of the cumulants do so approximately. Thus, although it is unlikely such cases will be encountered in inferential contexts, it is worth checking (by simulation) whether cumulants accurately replicate those of the target GQD. Furthermore, it may in some cases occur that the cumulants are indeed accurate whilst the density approximation fails. This can again be verified by simulation after which alternate density approximants may be specified through the variable Dtype.

Warning [2]: The parameter delt is also used as the stepsize for solving a system of ordinary differential equations (ODEs) that govern the evolution of the cumulants of the diffusion. As such delt is required to be small for highly non-linear models in order to ensure sufficient accuracy.

Author(s)

Etienne A.D. Pienaar: etiannead@gmail.com

References

Updates available on GitHub at https://github.com/eta21.

Daniels, H.E. 1954 Saddlepoint approximations in statistics. Ann. Math. Stat., 25:631–650.

Feagin, T. 2007 A tenth-order Runge-Kutta method with error estimate. In Proceedings of the IAENG Conf. on Scientifc Computing.

Varughese, M.M. 2013 Parameter estimation for multivariate diffusion systems. Comput. Stat. Data An., 57:417–428.

See Also

See GQD.mcmc and GQD.mle for likelihood based inference procedures for GQDs.

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
#===============================================================================
# Generate the transition density of a CIR process with time dependant
# drift and volatility.
#-------------------------------------------------------------------------------
  
  # Remove any existing coefficients
  GQD.remove()         
  
  # Define drift Coefficients. Note that the limiting mean is sinusoidal.
  G0 <- function(t){2*(10+sin(2*pi*(t-0.5)))}    
  G1 <- function(t){-2}
  
  # Define sinusoidal diffusion coefficient with `faster' oscillation.
  Q1 <- function(t){0.25*(1+0.75*(sin(4*pi*t)))}
  
  states     <-  seq(5,15,1/10) # State values
  initial    <-  8              # Starting value of the process
  Tmax       <-  5              # Time horizon
  Tstart     <-  1              # Time starts at 1
  increment  <-  1/100          # Incremental time steps
  
  # Generate the transitional density
  M <- GQD.density(Xs=initial,Xt=states,s=Tstart,t=Tmax,delt=increment)
  
  if(require(rgl, quietly = TRUE))
  {
    open3d(windowRect=c(50,50,640+50,50+640),zoom=0.95)
    persp3d(x=M$Xt,y=M$time,z=M$density,col=3,box=FALSE,xlab='State (X_t)',
            ylab='Time (t)',zlab='Density f(X_t|X_s)')
    play3d(spin3d(axis=c(0,0,1), rpm=3), duration=10)
  }else
  {
    persp(x=M$Xt,y=M$time,z=M$density,col=3,xlab='State (X_t)',ylab='Time (t)',
          zlab='Density f(X_t|X_s)',border=NA,shade=0.5,theta=145)
  }

#===============================================================================
# Generate the transition density for a diffusion process with restricted domain.
# The diffusion has reflective boundaries at 0 and 1.
#-------------------------------------------------------------------------------
  
  GQD.remove()         
  
  G0 <- function(t){0.4*(0.5+sin(2*pi*t))}    
  G1 <- function(t){-0.4}
  Q1 <- function(t){0.25}
  Q2 <- function(t){-0.25}
  
  states     <-  seq(0.005,0.995,1/200) 
  initial    <-  0.5              
  Tmax       <-  5              
  increment  <-  1/50          
  
  # Generate the transitional density
  M <- GQD.density(Xs=initial,Xt=states,s=0,t=Tmax,delt=increment,
                   Dtype='Beta',Trunc=c(8,8))
  
  if(require(rgl, quietly = TRUE))
  { 
    open3d(windowRect=c(50,50,640+50,50+640),zoom=0.95)
    persp3d(x=M$Xt,y=M$time,z=M$density,col='steelblue',box=FALSE,
            xlab='State (X_t)',ylab='Time (t)',zlab='Density f(X_t|X_s)')
    play3d(spin3d(axis=c(0,0,1), rpm=3), duration=10)
  }else
  {
    persp(x=M$Xt,y=M$time,z=M$density,col=3,xlab='State (X_t)',ylab='Time (t)',
          zlab='Density f(X_t|X_s)',border=NA,shade=0.5,theta=145)
  }
#===============================================================================

eta21/DiffusionRgqd documentation built on May 16, 2019, 8:53 a.m.