Package 'SMFilter'

Title: Filtering Algorithms for the State Space Models on the Stiefel Manifold
Description: Provides the filtering algorithms for the state space models on the Stiefel manifold as well as the corresponding sampling algorithms for uniform, vector Langevin-Bingham and matrix Langevin-Bingham distributions on the Stiefel manifold.
Authors: Yukai Yang [aut, cre]
Maintainer: Yukai Yang <[email protected]>
License: GPL-3
Version: 1.0.4
Built: 2025-02-11 05:12:28 UTC
Source: https://github.com/yukai-yang/smfilter

Help Index


Compute the squared Frobenius distance between two matrices.

Description

This function Compute the squared Frobenius distance between two matrices.

Usage

FDist2(mX, mY)

Arguments

mX

a p×rp \times r matrix where prp \geq r.

mY

another p×rp \times r matrix where prp \geq r.

Details

The Frobenius distance between two matrices is defined to be

d(X,Y)=tr{AA}d(X, Y) = \sqrt{ \mathrm{tr} \{ A' A \} }

where A=XYA = X - Y.

The Frobenius distance is a possible measure of the distance between two points on the Stiefel manifold.

Value

the Frobenius distance.

Author(s)

Yukai Yang, [email protected]

Examples

FDist2(runif_sm(1,4,2)[1,,], runif_sm(1,4,2)[1,,])

Filtering algorithm for the type one model.

Description

This function implements the filtering algorithm for the type one model. See Details part below.

Usage

FilterModel1(mY, mX, mZ, beta, mB = NULL, Omega, vD, U0,
  method = "max_1")

Arguments

mY

the matrix containing Y_t with dimension T×pT \times p.

mX

the matrix containing X_t with dimension T×q1T \times q_1.

mZ

the matrix containing Z_t with dimension T×q2T \times q_2.

beta

the β\beta matrix.

mB

the coefficient matrix B\boldsymbol{B} before mZ with dimension p×q2p \times q_2.

Omega

covariance matrix of the errors.

vD

vector of the diagonals of DD.

U0

initial value of the alpha sequence.

method

a string representing the optimization method from c('max_1','max_2','max_3','min_1','min_2').

Details

The type one model on Stiefel manifold takes the form:

yt=αtβxt+Bzt+εt\boldsymbol{y}_t \quad = \quad \boldsymbol{\alpha}_t \boldsymbol{\beta} ' \boldsymbol{x}_t + \boldsymbol{B} \boldsymbol{z}_t + \boldsymbol{\varepsilon}_t

αt+1αtML(p,r,αtD)\boldsymbol{\alpha}_{t+1} | \boldsymbol{\alpha}_{t} \quad \sim \quad ML (p, r, \boldsymbol{\alpha}_{t} \boldsymbol{D})

where yt\boldsymbol{y}_t is a pp-vector of the dependent variable, xt\boldsymbol{x}_t and zt\boldsymbol{z}_t are explanatory variables wit dimension q1q_1 and q2q_2, xt\boldsymbol{x}_t and zt\boldsymbol{z}_t have no overlap, matrix B\boldsymbol{B} is the coefficients for zt\boldsymbol{z}_t, εt\boldsymbol{\varepsilon}_t is the error vector.

The matrices αt\boldsymbol{\alpha}_t and β\boldsymbol{\beta} have dimensions p×rp \times r and q1×rq_1 \times r, respectively. Note that rr is strictly smaller than both pp and q1q_1. αt\boldsymbol{\alpha}_t and β\boldsymbol{\beta} are both non-singular matrices. αt\boldsymbol{\alpha}_t is time-varying while β\boldsymbol{\beta} is time-invariant.

Furthermore, αt\boldsymbol{\alpha}_t fulfills the condition αtαt=Ir\boldsymbol{\alpha}_t' \boldsymbol{\alpha}_t = \boldsymbol{I}_r, and therefor it evolves on the Stiefel manifold.

ML(p,r,αtD)ML (p, r, \boldsymbol{\alpha}_{t} \boldsymbol{D}) denotes the Matrix Langevin distribution or matrix von Mises-Fisher distribution on the Stiefel manifold. Its density function takes the form

f(αt+1)=etr{Dαtαt+1}0F1(p2;14D2)f(\boldsymbol{\alpha_{t+1}} ) = \frac{ \mathrm{etr} \left\{ \boldsymbol{D} \boldsymbol{\alpha}_{t}' \boldsymbol{\alpha_{t+1}} \right\} }{ _{0}F_1 (\frac{p}{2}; \frac{1}{4}\boldsymbol{D}^2 ) }

where etr\mathrm{etr} denotes exp(tr())\mathrm{exp}(\mathrm{tr}()), and 0F1(p2;14D2)_{0}F_1 (\frac{p}{2}; \frac{1}{4}\boldsymbol{D}^2 ) is the (0,1)-type hypergeometric function for matrix.

Value

an array aAlpha containing the modal orientations of alpha in the prediction step.

Author(s)

Yukai Yang, [email protected]

Examples

iT = 50
ip = 2
ir = 1
iqx = 4
iqz=0
ik = 0
Omega = diag(ip)*.1

if(iqx==0) mX=NULL else mX = matrix(rnorm(iT*iqx),iT, iqx)
if(iqz==0) mZ=NULL else mZ = matrix(rnorm(iT*iqz),iT, iqz)
if(ik==0) mY=NULL else mY = matrix(0, ik, ip)

alpha_0 = matrix(c(runif_sm(num=1,ip=ip,ir=ir)), ip, ir)
beta = matrix(c(runif_sm(num=1,ip=ip*ik+iqx,ir=ir)), ip*ik+iqx, ir)
mB=NULL
vD = 100

ret = SimModel1(iT=iT, mX=mX, mZ=mZ, mY=mY, alpha_0=alpha_0, beta=beta, mB=mB, vD=vD, Omega=Omega)
mYY=as.matrix(ret$dData[,1:ip])
fil = FilterModel1(mY=mYY, mX=mX, mZ=mZ, beta=beta, mB=mB, Omega=Omega, vD=vD, U0=alpha_0)

Filtering algorithm for the type two model.

Description

This function implements the filtering algorithm for the type two model. See Details part below.

Usage

FilterModel2(mY, mX, mZ, alpha, mB = NULL, Omega, vD, U0,
  method = "max_1")

Arguments

mY

the matrix containing Y_t with dimension T×pT \times p.

mX

the matrix containing X_t with dimension T×q1T \times q_1.

mZ

the matrix containing Z_t with dimension T×q2T \times q_2.

alpha

the α\alpha matrix.

mB

the coefficient matrix B\boldsymbol{B} before mZ with dimension p×q2p \times q_2.

Omega

covariance matrix of the errors.

vD

vector of the diagonals of DD.

U0

initial value of the alpha sequence.

method

a string representing the optimization method from c('max_1','max_2','max_3','min_1','min_2').

Details

The type two model on Stiefel manifold takes the form:

yt=αβtxt+Bzt+εt\boldsymbol{y}_t \quad = \quad \boldsymbol{\alpha} \boldsymbol{\beta}_t ' \boldsymbol{x}_t + \boldsymbol{B}' \boldsymbol{z}_t + \boldsymbol{\varepsilon}_t

βt+1βtML(q1,r,βtD)\boldsymbol{\beta}_{t+1} | \boldsymbol{\beta}_{t} \quad \sim \quad ML (q_1, r, \boldsymbol{\beta}_{t} \boldsymbol{D})

where yt\boldsymbol{y}_t is a pp-vector of the dependent variable, xt\boldsymbol{x}_t and zt\boldsymbol{z}_t are explanatory variables wit dimension q1q_1 and q2q_2, xt\boldsymbol{x}_t and zt\boldsymbol{z}_t have no overlap, matrix B\boldsymbol{B} is the coefficients for zt\boldsymbol{z}_t, εt\boldsymbol{\varepsilon}_t is the error vector.

The matrices α\boldsymbol{\alpha} and βt\boldsymbol{\beta}_t have dimensions p×rp \times r and q1×rq_1 \times r, respectively. Note that rr is strictly smaller than both pp and q1q_1. α\boldsymbol{\alpha} and βt\boldsymbol{\beta}_t are both non-singular matrices. βt\boldsymbol{\beta}_t is time-varying while α\boldsymbol{\alpha} is time-invariant.

Furthermore, βt\boldsymbol{\beta}_t fulfills the condition βtβt=Ir\boldsymbol{\beta}_t' \boldsymbol{\beta}_t = \boldsymbol{I}_r, and therefor it evolves on the Stiefel manifold.

ML(p,r,βtD)ML (p, r, \boldsymbol{\beta}_t \boldsymbol{D}) denotes the Matrix Langevin distribution or matrix von Mises-Fisher distribution on the Stiefel manifold. Its density function takes the form

f(βt+1)=etr{Dβtβt+1}0F1(p2;14D2)f(\boldsymbol{\beta_{t+1}} ) = \frac{ \mathrm{etr} \left\{ \boldsymbol{D} \boldsymbol{\beta}_{t}' \boldsymbol{\beta_{t+1}} \right\} }{ _{0}F_1 (\frac{p}{2}; \frac{1}{4}\boldsymbol{D}^2 ) }

where etr\mathrm{etr} denotes exp(tr())\mathrm{exp}(\mathrm{tr}()), and 0F1(p2;14D2)_{0}F_1 (\frac{p}{2}; \frac{1}{4}\boldsymbol{D}^2 ) is the (0,1)-type hypergeometric function for matrix.

Value

an array aAlpha containing the modal orientations of alpha in the prediction step.

Author(s)

Yukai Yang, [email protected]

Examples

iT = 50
ip = 2
ir = 1
iqx = 4
iqz=0
ik = 0
Omega = diag(ip)*.1

if(iqx==0) mX=NULL else mX = matrix(rnorm(iT*iqx),iT, iqx)
if(iqz==0) mZ=NULL else mZ = matrix(rnorm(iT*iqz),iT, iqz)
if(ik==0) mY=NULL else mY = matrix(0, ik, ip)

alpha = matrix(c(runif_sm(num=1,ip=ip,ir=ir)), ip, ir)
beta_0 = matrix(c(runif_sm(num=1,ip=ip*ik+iqx,ir=ir)), ip*ik+iqx, ir)
mB=NULL
vD = 100

ret = SimModel2(iT=iT, mX=mX, mZ=mZ, mY=mY, alpha=alpha, beta_0=beta_0, mB=mB, vD=vD)
mYY=as.matrix(ret$dData[,1:ip])
fil = FilterModel2(mY=mYY, mX=mX, mZ=mZ, alpha=alpha, mB=mB, Omega=Omega, vD=vD, U0=beta_0)

Sample from the matrix Langevin-Bingham on the Stiefel manifold.

Description

This function draws a sample from the matrix Langevin-Bingham on the Stiefel manifold.

Usage

rmLB_sm(num, mJ, mH, mC, mX, ir)

Arguments

num

number of observations or sample size.

mJ

symmetric ip*ip matrix

mH

symmetric ir*ir matrix

mC

ip*ir matrix

mX

ip*ir matrix, the initial value

ir

ir

Details

The matrix Langevin-Bingham distribution on the Stiefel manifold has the density kernel:

f(X)etr{HXJX+CX}f(X) \propto \mathrm{etr}\{ H X' J X + C' X \}

where XX satisfies XX=IrX'X = I_r, and HH and JJ are symmetric matrices.

Value

an array containing a sample of draws from the matrix Langevin-Bingham on the Stiefel manifold.

Author(s)

Yukai Yang, [email protected]

#' @section References: Hoff, P. D. (2009) "Simulation of the Matrix Bingham—von Mises—Fisher Distribution, With Applications to Multivariate and Relational Data", Journal of Computational and Graphical Statistics, Vol. 18, pp. 438-456.


Sample from the uniform distribution on the Stiefel manifold.

Description

This function draws a sample from the uniform distribution on the Stiefel manifold.

Usage

runif_sm(num, ip, ir)

Arguments

num

number of observations or sample size.

ip

the first dimension pp of the matrix.

ir

the second dimension rr of the matrix.

Details

The Stiefel manifold with dimension pp and rr (prp \geq r) is a space whose points are rr-frames in RpR^p. A set of rr orthonormal vectors in RpR^p is called an rr-frame in RpR^p. The Stiefel manifold is a collection of p×rp \times r full rank matrices XX such that XX=IrX'X = I_r.

Value

an array with dimension num, ip and ir containing a sample of draws from the uniform distribution on the Stiefel manifold.

Author(s)

Yukai Yang, [email protected]

Examples

runif_sm(10,4,2)

Sample from the vector Langevin-Bingham on the Stiefel manifold.

Description

This function draws a sample from the vector Langevin-Bingham on the Stiefel manifold.

Usage

rvlb_sm(num, mA, vc, vx)

Arguments

num

number of observations or sample size.

mA

the matrix A which is symmetric ip*ip matrix.

vc

the vector c with dimension ip.

vx

the vector x, the initial value.

Details

The vector Langevin-Bingham distribution on the Stiefel manifold has the density kernel:

f(X)etr{xAx+cx}f(X) \propto \mathrm{etr}\{ x' A x + c' x \}

where xx satisfies xx=1x'x = 1, and AA is a symmetric matrix.

Value

an array containing a sample of draws from the vector Langevin-Bingham on the Stiefel manifold.

References

Hoff, P. D. (2009) "Simulation of the Matrix Bingham—von Mises—Fisher Distribution, With Applications to Multivariate and Relational Data", Journal of Computational and Graphical Statistics, Vol. 18, pp. 438-456.

Author(s)

Yukai Yang, [email protected]


Simulate from the type one state-space Model on Stiefel manifold.

Description

This function simulates from the type one model on Stiefel manifold. See Details part below.

Usage

SimModel1(iT, mX = NULL, mZ = NULL, mY = NULL, alpha_0, beta,
  mB = NULL, Omega = NULL, vD, burnin = 100)

Arguments

iT

the sample size.

mX

the matrix containing X_t with dimension T×q1T \times q_1.

mZ

the matrix containing Z_t with dimension T×q2T \times q_2.

mY

initial values of the dependent variable for ik-1 up to 0. If mY = NULL, then no lagged dependent variables in regressors.

alpha_0

the initial alpha, p×rp \times r.

beta

the β\beta matrix, iqx+ip*ik, y_1,t-1,y_1,t-2,...,y_2,t-1,y_2,t-2,...

mB

the coefficient matrix B\boldsymbol{B} before mZ with dimension p×q2p \times q_2.

Omega

covariance matrix of the errors.

vD

vector of the diagonals of DD.

burnin

burn-in sample size (matrix Langevin).

Details

The type one model on Stiefel manifold takes the form:

yt=αtβxt+Bzt+εt\boldsymbol{y}_t \quad = \quad \boldsymbol{\alpha}_t \boldsymbol{\beta} ' \boldsymbol{x}_t + \boldsymbol{B} \boldsymbol{z}_t + \boldsymbol{\varepsilon}_t

αt+1αtML(p,r,αtD)\boldsymbol{\alpha}_{t+1} | \boldsymbol{\alpha}_{t} \quad \sim \quad ML (p, r, \boldsymbol{\alpha}_{t} \boldsymbol{D})

where yt\boldsymbol{y}_t is a pp-vector of the dependent variable, xt\boldsymbol{x}_t and zt\boldsymbol{z}_t are explanatory variables wit dimension q1q_1 and q2q_2, xt\boldsymbol{x}_t and zt\boldsymbol{z}_t have no overlap, matrix B\boldsymbol{B} is the coefficients for zt\boldsymbol{z}_t, εt\boldsymbol{\varepsilon}_t is the error vector.

The matrices αt\boldsymbol{\alpha}_t and β\boldsymbol{\beta} have dimensions p×rp \times r and q1×rq_1 \times r, respectively. Note that rr is strictly smaller than both pp and q1q_1. αt\boldsymbol{\alpha}_t and β\boldsymbol{\beta} are both non-singular matrices. αt\boldsymbol{\alpha}_t is time-varying while β\boldsymbol{\beta} is time-invariant.

Furthermore, αt\boldsymbol{\alpha}_t fulfills the condition αtαt=Ir\boldsymbol{\alpha}_t' \boldsymbol{\alpha}_t = \boldsymbol{I}_r, and therefor it evolves on the Stiefel manifold.

ML(p,r,αtD)ML (p, r, \boldsymbol{\alpha}_{t} \boldsymbol{D}) denotes the Matrix Langevin distribution or matrix von Mises-Fisher distribution on the Stiefel manifold. Its density function takes the form

f(αt+1)=etr{Dαtαt+1}0F1(p2;14D2)f(\boldsymbol{\alpha_{t+1}} ) = \frac{ \mathrm{etr} \left\{ \boldsymbol{D} \boldsymbol{\alpha}_{t}' \boldsymbol{\alpha_{t+1}} \right\} }{ _{0}F_1 (\frac{p}{2}; \frac{1}{4}\boldsymbol{D}^2 ) }

where etr\mathrm{etr} denotes exp(tr())\mathrm{exp}(\mathrm{tr}()), and 0F1(p2;14D2)_{0}F_1 (\frac{p}{2}; \frac{1}{4}\boldsymbol{D}^2 ) is the (0,1)-type hypergeometric function for matrix.

Note that the function does not add intercept automatically.

Value

A list containing the sampled data and the dynamics of alpha.

The object is a list containing the following components:

dData

a data.frame of the sampled data

aAlpha

an array of the αt\boldsymbol{\alpha}_{t} with the dimension T×p×rT \times p \times r

Author(s)

Yukai Yang, [email protected]

Examples

iT = 50 # sample size
ip = 2 # dimension of the dependent variable
ir = 1 # rank number
iqx=2 # number of variables in X
iqz=2 # number of variables in Z
ik = 1 # lag length

if(iqx==0) mX=NULL else mX = matrix(rnorm(iT*iqx),iT, iqx)
if(iqz==0) mZ=NULL else mZ = matrix(rnorm(iT*iqz),iT, iqz)
if(ik==0) mY=NULL else mY = matrix(0, ik, ip)

alpha_0 = matrix(c(runif_sm(num=1,ip=ip,ir=ir)), ip, ir)
beta = matrix(c(runif_sm(num=1,ip=ip*ik+iqx,ir=ir)), ip*ik+iqx, ir)
if(ip*ik+iqz==0) mB=NULL else mB = matrix(c(runif_sm(num=1,ip=(ip*ik+iqz)*ip,ir=1)), ip, ip*ik+iqz)
vD = 50

ret = SimModel1(iT=iT, mX=mX, mZ=mZ, mY=mY, alpha_0=alpha_0, beta=beta, mB=mB, vD=vD)

Simulate from the type two state-space Model on Stiefel manifold.

Description

This function simulates from the type two model on Stiefel manifold. See Details part below.

Usage

SimModel2(iT, mX = NULL, mZ = NULL, mY = NULL, beta_0, alpha,
  mB = NULL, Omega = NULL, vD, burnin = 100)

Arguments

iT

the sample size.

mX

the matrix containing X_t with dimension T×q1T \times q_1.

mZ

the matrix containing Z_t with dimension T×q2T \times q_2.

mY

initial values of the dependent variable for ik-1 up to 0. If mY = NULL, then no lagged dependent variables in regressors.

beta_0

the initial beta, iqx+ip*ik, y_1,t-1,y_1,t-2,...,y_2,t-1,y_2,t-2,....

alpha

the α\alpha matrix, p×rp \times r.

mB

the coefficient matrix B\boldsymbol{B} before mZ with dimension p×q2p \times q_2.

Omega

covariance matrix of the errors.

vD

vector of the diagonals of DD.

burnin

burn-in sample size (matrix Langevin).

Details

The type two model on Stiefel manifold takes the form:

yt=αβtxt+Bzt+εt\boldsymbol{y}_t \quad = \quad \boldsymbol{\alpha} \boldsymbol{\beta}_t ' \boldsymbol{x}_t + \boldsymbol{B}' \boldsymbol{z}_t + \boldsymbol{\varepsilon}_t

βt+1βtML(q1,r,βtD)\boldsymbol{\beta}_{t+1} | \boldsymbol{\beta}_{t} \quad \sim \quad ML (q_1, r, \boldsymbol{\beta}_{t} \boldsymbol{D})

where yt\boldsymbol{y}_t is a pp-vector of the dependent variable, xt\boldsymbol{x}_t and zt\boldsymbol{z}_t are explanatory variables wit dimension q1q_1 and q2q_2, xt\boldsymbol{x}_t and zt\boldsymbol{z}_t have no overlap, matrix B\boldsymbol{B} is the coefficients for zt\boldsymbol{z}_t, εt\boldsymbol{\varepsilon}_t is the error vector.

The matrices α\boldsymbol{\alpha} and βt\boldsymbol{\beta}_t have dimensions p×rp \times r and q1×rq_1 \times r, respectively. Note that rr is strictly smaller than both pp and q1q_1. α\boldsymbol{\alpha} and βt\boldsymbol{\beta}_t are both non-singular matrices. βt\boldsymbol{\beta}_t is time-varying while α\boldsymbol{\alpha} is time-invariant.

Furthermore, βt\boldsymbol{\beta}_t fulfills the condition βtβt=Ir\boldsymbol{\beta}_t' \boldsymbol{\beta}_t = \boldsymbol{I}_r, and therefor it evolves on the Stiefel manifold.

ML(p,r,βtD)ML (p, r, \boldsymbol{\beta}_t \boldsymbol{D}) denotes the Matrix Langevin distribution or matrix von Mises-Fisher distribution on the Stiefel manifold. Its density function takes the form

f(βt+1)=etr{Dβtβt+1}0F1(p2;14D2)f(\boldsymbol{\beta_{t+1}} ) = \frac{ \mathrm{etr} \left\{ \boldsymbol{D} \boldsymbol{\beta}_{t}' \boldsymbol{\beta_{t+1}} \right\} }{ _{0}F_1 (\frac{p}{2}; \frac{1}{4}\boldsymbol{D}^2 ) }

where etr\mathrm{etr} denotes exp(tr())\mathrm{exp}(\mathrm{tr}()), and 0F1(p2;14D2)_{0}F_1 (\frac{p}{2}; \frac{1}{4}\boldsymbol{D}^2 ) is the (0,1)-type hypergeometric function for matrix.

Note that the function does not add intercept automatically.

Value

A list containing the sampled data and the dynamics of beta.

The object is a list containing the following components:

dData

a data.frame of the sampled data

aBeta

an array of the βt\boldsymbol{\beta}_t with the dimension T×q1×rT \times q_1 \times r

Author(s)

Yukai Yang, [email protected]

Examples

iT = 50
ip = 2
ir = 1
iqx =3
iqz=2
ik = 1

if(iqx==0) mX=NULL else mX = matrix(rnorm(iT*iqx),iT, iqx)
if(iqz==0) mZ=NULL else mZ = matrix(rnorm(iT*iqz),iT, iqz)
if(ik==0) mY=NULL else mY = matrix(0, ik, ip)

alpha = matrix(c(runif_sm(num=1,ip=ip,ir=ir)), ip, ir)
beta_0 = matrix(c(runif_sm(num=1,ip=ip*ik+iqx,ir=ir)), ip*ik+iqx, ir)
if(ip*ik+iqz==0) mB=NULL else mB = matrix(c(runif_sm(num=1,ip=(ip*ik+iqz)*ip,ir=1)), ip, ip*ik+iqz)
vD = 50

ret = SimModel2(iT=iT, mX=mX, mZ=mZ, mY=mY, alpha=alpha, beta_0=beta_0, mB=mB, vD=vD)

SMFilter: a package implementing the filtering algorithms for the state-space models on the Stiefel manifold.

Description

The package implements the filtering algorithms for the state-space models on the Stiefel manifold. It also implements sampling algorithms for uniform, vector Langevin-Bingham and matrix Langevin-Bingham distributions on the Stiefel manifold.

Details

Two types of the state-space models on the Stiefel manifold are considered.

The type one model on Stiefel manifold takes the form:

yt=αtβxt+Bzt+εt\boldsymbol{y}_t \quad = \quad \boldsymbol{\alpha}_t \boldsymbol{\beta} ' \boldsymbol{x}_t + \boldsymbol{B} \boldsymbol{z}_t + \boldsymbol{\varepsilon}_t

αt+1αtML(p,r,αtD)\boldsymbol{\alpha}_{t+1} | \boldsymbol{\alpha}_{t} \quad \sim \quad ML (p, r, \boldsymbol{\alpha}_{t} \boldsymbol{D})

where yt\boldsymbol{y}_t is a pp-vector of the dependent variable, xt\boldsymbol{x}_t and zt\boldsymbol{z}_t are explanatory variables wit dimension q1q_1 and q2q_2, xt\boldsymbol{x}_t and zt\boldsymbol{z}_t have no overlap, matrix B\boldsymbol{B} is the coefficients for zt\boldsymbol{z}_t, εt\boldsymbol{\varepsilon}_t is the error vector.

The matrices αt\boldsymbol{\alpha}_t and β\boldsymbol{\beta} have dimensions p×rp \times r and q1×rq_1 \times r, respectively. Note that rr is strictly smaller than both pp and q1q_1. αt\boldsymbol{\alpha}_t and β\boldsymbol{\beta} are both non-singular matrices. αt\boldsymbol{\alpha}_t is time-varying while β\boldsymbol{\beta} is time-invariant.

Furthermore, αt\boldsymbol{\alpha}_t fulfills the condition αtαt=Ir\boldsymbol{\alpha}_t' \boldsymbol{\alpha}_t = \boldsymbol{I}_r, and therefor it evolves on the Stiefel manifold.

ML(p,r,αtD)ML (p, r, \boldsymbol{\alpha}_{t} \boldsymbol{D}) denotes the Matrix Langevin distribution or matrix von Mises-Fisher distribution on the Stiefel manifold. Its density function takes the form

f(αt+1)=etr{Dαtαt+1}0F1(p2;14D2)f(\boldsymbol{\alpha_{t+1}} ) = \frac{ \mathrm{etr} \left\{ \boldsymbol{D} \boldsymbol{\alpha}_{t}' \boldsymbol{\alpha_{t+1}} \right\} }{ _{0}F_1 (\frac{p}{2}; \frac{1}{4}\boldsymbol{D}^2 ) }

where etr\mathrm{etr} denotes exp(tr())\mathrm{exp}(\mathrm{tr}()), and 0F1(p2;14D2)_{0}F_1 (\frac{p}{2}; \frac{1}{4}\boldsymbol{D}^2 ) is the (0,1)-type hypergeometric function for matrix.

The type two model on Stiefel manifold takes the form:

yt=αβtxt+Bzt+εt\boldsymbol{y}_t \quad = \quad \boldsymbol{\alpha} \boldsymbol{\beta}_t ' \boldsymbol{x}_t + \boldsymbol{B}' \boldsymbol{z}_t + \boldsymbol{\varepsilon}_t

βt+1βtML(q1,r,βtD)\boldsymbol{\beta}_{t+1} | \boldsymbol{\beta}_{t} \quad \sim \quad ML (q_1, r, \boldsymbol{\beta}_{t} \boldsymbol{D})

where yt\boldsymbol{y}_t is a pp-vector of the dependent variable, xt\boldsymbol{x}_t and zt\boldsymbol{z}_t are explanatory variables wit dimension q1q_1 and q2q_2, xt\boldsymbol{x}_t and zt\boldsymbol{z}_t have no overlap, matrix B\boldsymbol{B} is the coefficients for zt\boldsymbol{z}_t, εt\boldsymbol{\varepsilon}_t is the error vector.

The matrices α\boldsymbol{\alpha} and βt\boldsymbol{\beta}_t have dimensions p×rp \times r and q1×rq_1 \times r, respectively. Note that rr is strictly smaller than both pp and q1q_1. α\boldsymbol{\alpha} and βt\boldsymbol{\beta}_t are both non-singular matrices. βt\boldsymbol{\beta}_t is time-varying while α\boldsymbol{\alpha} is time-invariant.

Furthermore, βt\boldsymbol{\beta}_t fulfills the condition βtβt=Ir\boldsymbol{\beta}_t' \boldsymbol{\beta}_t = \boldsymbol{I}_r, and therefor it evolves on the Stiefel manifold.

ML(p,r,βtD)ML (p, r, \boldsymbol{\beta}_t \boldsymbol{D}) denotes the Matrix Langevin distribution or matrix von Mises-Fisher distribution on the Stiefel manifold. Its density function takes the form

f(βt+1)=etr{Dβtβt+1}0F1(p2;14D2)f(\boldsymbol{\beta_{t+1}} ) = \frac{ \mathrm{etr} \left\{ \boldsymbol{D} \boldsymbol{\beta}_{t}' \boldsymbol{\beta_{t+1}} \right\} }{ _{0}F_1 (\frac{p}{2}; \frac{1}{4}\boldsymbol{D}^2 ) }

where etr\mathrm{etr} denotes exp(tr())\mathrm{exp}(\mathrm{tr}()), and 0F1(p2;14D2)_{0}F_1 (\frac{p}{2}; \frac{1}{4}\boldsymbol{D}^2 ) is the (0,1)-type hypergeometric function for matrix.

Author and Maintainer

Yukai Yang

Department of Statistics, Uppsala University

[email protected]

References

Yang, Yukai and Bauwens, Luc. (2018) "State-Space Models on the Stiefel Manifold with a New Approach to Nonlinear Filtering", Econometrics, 6(4), 48.

Simulation

SimModel1 simulate from the type one state-space model on the Stiefel manifold.

SimModel2 simulate from the type two state-space model on the Stiefel manifold.

Filtering

FilterModel1 filtering algorithm for the type one model.

FilterModel2 filtering algorithm for the type two model.

Sampling

runif_sm sample from the uniform distribution on the Stiefel manifold.

rvlb_sm sample from the vector Langevin-Bingham distribution on the Stiefel manifold.

rmLB_sm sample from the matrix Langevin-Bingham distribution on the Stiefel manifold.

Other Functions

version shows the version number and some information of the package.


Show the version number of some information.

Description

This function shows the version number and some information of the package.

Usage

version()

Author(s)

Yukai Yang, [email protected]