Package 'multiwave'

Title: Estimation of Multivariate Long-Memory Models Parameters
Description: Computation of an estimation of the long-memory parameters and the long-run covariance matrix using a multivariate model (Lobato (1999) <doi:10.1016/S0304-4076(98)00038-4>; Shimotsu (2007) <doi:10.1016/j.jeconom.2006.01.003>). Two semi-parametric methods are implemented: a Fourier based approach (Shimotsu (2007) <doi:10.1016/j.jeconom.2006.01.003>) and a wavelet based approach (Achard and Gannaz (2016) <doi:10.1111/jtsa.12170>).
Authors: Sophie Achard [aut, cre], Irene Gannaz [aut]
Maintainer: Sophie Achard <[email protected]>
License: GPL (>= 2)
Version: 1.4
Built: 2025-03-03 05:33:14 UTC
Source: https://github.com/cran/multiwave

Help Index


Estimation of multivariate long-memory models parameters: memory parameters and long-run covariance matrix (also called fractal connectivity).

Description

This package computes an estimation of the long-memory parameters and the long-run covariance matrix using a multivariate model (Lobato, 1999; Shimotsu 2007). Two semi-parametric methods are implemented: a Fourier based approach (Shimotsu 2007) and a wavelet based approach (Achard and Gannaz 2014).

Details

Package: multiwave
Type: Package
Version: 1.0
Date: 2015-09-17
License: GPL (>= 2)

Author(s)

Sophie Achard and Irene Gannaz

Maintainer: Sophie Achard <[email protected]>, Irene Gannaz <[email protected]>

References

S. Achard, I. Gannaz (2016) Multivariate wavelet Whittle estimation in long-range dependence. Journal of Time Series Analysis, Vol 37, N. 4, pages 476-512. http://arxiv.org/abs/1412.0391.

S. Achard, I Gannaz (2019) Wavelet-Based and Fourier-Based Multivariate Whittle Estimation: multiwave. Journal of Statistical Software, Vol 89, N. 6, pages 1-31.

Examples

rho<-0.4
cov <- matrix(c(1,rho,rho,1),2,2)
d<-c(0.4,0.2)
J <- 9
N <- 2^J

resp <- fivarma(N, d, cov_matrix=cov)

x <- resp$x
long_run_cov <- resp$long_run_cov

#### Compute wavelets this is also included in the functions without _wav
res_filter <- scaling_filter('Daubechies',8);
filter <- res_filter$h
M <- res_filter$M
alpha <- res_filter$alpha

LU <- c(1,11)

if(is.matrix(x)){
    N <- dim(x)[1]
    k <- dim(x)[2]
}else{
    N <- length(x)
    k <- 1
}
mat_x <- as.matrix(x,dim=c(N,k))

## Wavelet decomposition
xwav <- matrix(0,N,k)
    for(j in 1:k){
        xx <- mat_x[,j]   
        resw <- DWTexact(xx,filter)
        xwav_temp <- resw$dwt
        index <- resw$indmaxband
        Jmax <- resw$Jmax
        xwav[1:index[Jmax],j] <- xwav_temp;
    }
## we free some memory
new_xwav <- matrix(0,min(index[Jmax],N),k)
    if(index[Jmax]<N){
        new_xwav[(1:(index[Jmax])),] <- xwav[(1:(index[Jmax])),]
    }
    xwav <- new_xwav
    index <- c(0,index)


##### Compute the wavelet functions 
res_psi <- psi_hat_exact(filter,J)
psih<-res_psi$psih
grid<-res_psi$grid


##### Estimate using Fourier #############

m <- floor(N^{0.65}) ## default value of Shimotsu
res_mfw <- mfw(x,m)
res_d_mfw<-res_mfw$d
res_rho_mfw<-res_mfw$cov[1,2]

### Eval MFW

res_mfw_eval <- mfw_eval(d,x,m)
res_mfw_cov_eval <- mfw_cov_eval(d,x,m)

###### Estimate using Wavelets #############

## Using xwav 

if(dim(xwav)[2]==1) xwav<-as.vector(xwav)
res_mww_wav <- mww_wav(xwav,index,psih,grid,LU)

### Eval MWW_wav

res_mww_wav_eval <- mww_wav_eval(d,xwav,index,LU)
res_mww_wav_cov_eval <- mww_wav_cov_eval(d,xwav,index,psih,grid,LU)

## Using directly the time series

res_mww <- mww(x,filter,LU)
res_d_mww<-res_mww$d
res_rho_mww<-res_mww$cov[1,2]

### Eval MWW_wav

res_mww_eval <- mww_eval(d,x,filter,LU)
res_mww_cov_eval <- mww_cov_eval(d,x,filter,LU)

Time series obtained by an fMRI experiment on the brain

Description

Time series for each region of interest in the brain. These series are obtained by SPM preprocessing.

Usage

data(brainHCP)

Format

A data frame with 1200 observations on the following 89 variables.

Source

contact S. Achard ([email protected])

References

M. Termenon, A. Jaillard, C. Delon-Martin, S. Achard (2016) Reliability of graph analysis of resting state fMRI using test-retest dataset from the Human Connectome Project, Neuroimage, Vol 142, pages 172-187.

Examples

data(brainHCP)
## maybe str(brainHCP) ; plot(brainHCP) ...

Wavelets coefficients utilities

Description

Computes the number of wavelet coefficients at each scale.

Usage

compute_nj(n, N)

Arguments

n

sample size.

N

filter length.

Value

nj

number of coefficients at each scale.

J

Number of scales.

Author(s)

S. Achard and I. Gannaz

References

G. Fay, E. Moulines, F. Roueff, M. S. Taqqu (2009) Estimators of long-memory: Fourier versus wavelets. Journal of Econometrics, vol. 151, N. 2, pages 159-177.

S. Achard, I. Gannaz (2016) Multivariate wavelet Whittle estimation in long-range dependence. Journal of Time Series Analysis, Vol 37, N. 4, pages 476-512. http://arxiv.org/abs/1412.0391.

S. Achard, I Gannaz (2019) Wavelet-Based and Fourier-Based Multivariate Whittle Estimation: multiwave. Journal of Statistical Software, Vol 89, N. 6, pages 1-31.

See Also

DWTexact, scaling_filter

Examples

res_filter <- scaling_filter('Daubechies',8);
filter <- res_filter$h
n <- 5^10
N <- length(filter)
compute_nj(n,N)

Exact discrete wavelet decomposition

Description

Computes the discrete wavelet transform of the data using the pyramidal algorithm.

Usage

DWTexact(x, filter)

Arguments

x

vector of raw data

filter

Quadrature mirror filter (also called scaling filter, as returned by the scaling_filter function)

Value

dwt

computable Wavelet coefficients without taking into account the border effect.

indmaxband

vector containing the largest index of each band, i.e. for j>1j > 1 the wavelet coefficients of scale jj are dwt(k)\code{dwt}(k) for k[indmaxband(j1)+1,indmaxband(j)]k \in [\code{indmaxband}(j-1)+1,\code{indmaxband}(j)] and for j=1j=1, dwt(k)\code{dwt}(k) for k[1,indmaxband(1)]k \in [1,\code{indmaxband}(1)].

Jmax

largest available scale index (=length of indmaxband).

Note

This function was rewritten from an original matlab version by Fay et al. (2009)

Author(s)

S. Achard and I. Gannaz

References

G. Fay, E. Moulines, F. Roueff, M. S. Taqqu (2009) Estimators of long-memory: Fourier versus wavelets. Journal of Econometrics, vol. 151, N. 2, pages 159-177.

S. Achard, I. Gannaz (2016) Multivariate wavelet Whittle estimation in long-range dependence. Journal of Time Series Analysis, Vol 37, N. 4, pages 476-512. http://arxiv.org/abs/1412.0391.

S. Achard, I Gannaz (2019) Wavelet-Based and Fourier-Based Multivariate Whittle Estimation: multiwave. Journal of Statistical Software, Vol 89, N. 6, pages 1-31.

See Also

scaling_filter

Examples

res_filter <- scaling_filter('Daubechies',8);
filter <- res_filter$h
u <- rnorm(2^10,0,1)
x <- vfracdiff(u,d=0.2)

	resw <- DWTexact(x,filter)
		xwav <- resw$dwt
		index <- resw$indmaxband
		Jmax <- resw$Jmax

## Wavelet scale 1
ws_1 <- xwav[1:index[1]]
## Wavelet scale 2
ws_2 <- xwav[(index[1]+1):index[2]]
## Wavelet scale 3
ws_3 <- xwav[(index[2]+1):index[3]]
### upto Jmax

simulation of FIVARMA process

Description

Generates N observations of a realisation of a multivariate FIVARMA process X.

Usage

fivarma(N, d = 0, cov_matrix = diag(length(d)), VAR = NULL,
            VMA = NULL,skip = 2000)

Arguments

N

number of time points.

d

vector of parameters of long-memory.

cov_matrix

matrix of correlation between the innovations (optional, default is identity).

VAR

array of VAR coefficient matrices (optional).

VMA

array of VMA coefficient matrices (optional).

skip

number of initial observations omitted, after applying the ARMA operator and the fractional integration (optional, the default is 2000).

Details

Let (e(t))t(e(t))_t be a multivariate gaussian process with a covariance matrix cov_matrix. The values of the process X are given by the equations:

VAR(L)U(t)=VMA(L)e(t),VAR(L)U(t) = VMA(L)e(t),

and

diag((1L)d)X(t)=U(t)diag((1-L)^d)X(t) = U(t)

where L is the lag-operator.

Value

x

vector containing the N observations of the vector ARFIMA(arlags, d, malags) process.

long_run_cov

matrix of covariance of the spectral density of x around the zero frequency.

d

vector of parameters of long-range dependence, modified in case of cointegration.

Author(s)

S. Achard and I. Gannaz

References

R. J. Sela and C. M. Hurvich (2009) Computationaly efficient methods for two multivariate fractionnaly integrated models. Journal of Time Series Analysis, Vol 30, N. 6, pages 631-651.

S. Achard, I. Gannaz (2016) Multivariate wavelet Whittle estimation in long-range dependence. Journal of Time Series Analysis, Vol 37, N. 4, pages 476-512. http://arxiv.org/abs/1412.0391.

S. Achard, I Gannaz (2019) Wavelet-Based and Fourier-Based Multivariate Whittle Estimation: multiwave. Journal of Statistical Software, Vol 89, N. 6, pages 1-31.

See Also

varma, vfracdiff

Examples

rho1 <- 0.3
rho2 <- 0.8
cov <- matrix(c(1,rho1,rho2,rho1,1,rho1,rho2,rho1,1),3,3)
d <- c(0.2,0.3,0.4)

J <- 9
N <- 2^J
VMA <- diag(c(0.4,0.1,0))
### or another example VAR <- array(c(0.8,0,0,0,0.6,0,0,0,0.2,0,0,0,0,0.4,0,0,0,0.5),dim=c(3,3,2))
VAR <- diag(c(0.8,0.6,0))
resp <- fivarma(N, d, cov_matrix=cov, VAR=VAR, VMA=VMA)
x <- resp$x
long_run_cov <- resp$long_run_cov
d <- resp$d

Evaluation of function KK

Description

Computes the function KK as defined in (Achard and Gannaz 2014).

Usage

K_eval(psi_hat,u,d)

Arguments

psi_hat

Fourier transform of the wavelet mother at values u

u

grid for the approximation of the integral

d

vector of long-memory parameters.

Details

K_eval computes the matrix KK with elements

K(dl,dm)=u(dl+dm)psi_hat(u)2duK(d_l,d_m)=\int u^{(d_l+d_m)} |\code{psi}\_\code{hat}(u)|^2 du

Value

value of function K as a matrix.

Author(s)

S. Achard and I. Gannaz

References

S. Achard, I. Gannaz (2016) Multivariate wavelet Whittle estimation in long-range dependence. Journal of Time Series Analysis, Vol 37, N. 4, pages 476-512. http://arxiv.org/abs/1412.0391.

S. Achard, I Gannaz (2019) Wavelet-Based and Fourier-Based Multivariate Whittle Estimation: multiwave. Journal of Statistical Software, Vol 89, N. 6, pages 1-31.

See Also

psi_hat_exact

Examples

res_filter <- scaling_filter('Daubechies',8);
filter <- res_filter$h
M <- res_filter$M
alpha <- res_filter$alpha
res_psi <- psi_hat_exact(filter,J=10)
K_eval(res_psi$psih,res_psi$grid,d=c(0.2,0.2))

multivariate Fourier Whittle estimators

Description

Computes the multivariate Fourier Whittle estimators of the long-memory parameters and the long-run covariance matrix also called fractal connectivity.

Usage

mfw(x, m)

Arguments

x

data (matrix with time in rows and variables in columns).

m

truncation number used for the estimation of the periodogram.

Details

The choice of m determines the range of frequencies used in the computation of the periodogram, λj=2πj/N\lambda_j = 2\pi j/N, jj = 1,... , m. The optimal value depends on the spectral properties of the time series such as the presence of short range dependence. In Shimotsu (2007), m is chosen to be equal to N0.65N^{0.65}.

Value

d

estimation of the vector of long-memory parameters.

cov

estimation of the long-run covariance matrix.

Author(s)

S. Achard and I. Gannaz

References

K. Shimotsu (2007) Gaussian semiparametric estimation of multivariate fractionally integrated processes Journal of Econometrics Vol. 137, N. 2, pages 277-310.

S. Achard, I. Gannaz (2016) Multivariate wavelet Whittle estimation in long-range dependence. Journal of Time Series Analysis, Vol 37, N. 4, pages 476-512. http://arxiv.org/abs/1412.0391.

S. Achard, I Gannaz (2019) Wavelet-Based and Fourier-Based Multivariate Whittle Estimation: multiwave. Journal of Statistical Software, Vol 89, N. 6, pages 1-31.

See Also

mfw_eval, mfw_cov_eval

Examples

### Simulation of ARFIMA(0,d,0)
rho <- 0.4
cov <- matrix(c(1,rho,rho,1),2,2)
d <- c(0.4,0.2)
J <- 9
N <- 2^J

resp <- fivarma(N, d, cov_matrix=cov)
x <- resp$x
long_run_cov <- resp$long_run_cov

m <- 57 ## default value of Shimotsu 2007
res_mfw <- mfw(x,m)

multivariate Fourier Whittle estimators

Description

Computes the multivariate Fourier Whittle estimator of the long-run covariance matrix (also called fractal connectivity) for a given value of long-memory parameters d.

Usage

mfw_cov_eval(d, x, m)

Arguments

d

vector of long-memory parameters (dimension should match dimension of x)

x

data (matrix with time in rows and variables in columns)

m

truncation number used for the estimation of the periodogram

Details

The choice of m determines the range of frequencies used in the computation of the periodogram, λj=2πj/N\lambda_j = 2\pi j/N, jj = 1,... , m. The optimal value depends on the spectral properties of the time series such as the presence of short range dependence. In Shimotsu (2007), m is chosen to be equal to N0.65N^{0.65}.

Value

long-run covariance matrix estimation.

Author(s)

S. Achard and I. Gannaz

References

K. Shimotsu (2007) Gaussian semiparametric estimation of multivariate fractionally integrated processes Journal of Econometrics Vol. 137, N. 2, pages 277-310.

S. Achard, I. Gannaz (2016) Multivariate wavelet Whittle estimation in long-range dependence. Journal of Time Series Analysis, Vol 37, N. 4, pages 476-512. http://arxiv.org/abs/1412.0391.

S. Achard, I Gannaz (2019) Wavelet-Based and Fourier-Based Multivariate Whittle Estimation: multiwave. Journal of Statistical Software, Vol 89, N. 6, pages 1-31.

See Also

mfw_eval, mfw

Examples

### Simulation of ARFIMA(0,\code{d},0)
rho <- 0.4
cov <- matrix(c(1,rho,rho,1),2,2)
d <- c(0.4,0.2)
J <- 9
N <- 2^J

resp <- fivarma(N, d, cov_matrix=cov)
x <- resp$x
long_run_cov <- resp$long_run_cov

m <- 57 ## default value of Shimotsu
G <- mfw_cov_eval(d,x,m) # estimation of the covariance matrix when d is known

evaluation of multivariate Fourier Whittle estimator

Description

Evaluates the multivariate Fourier Whittle criterion at a given long-memory parameter value d.

Usage

mfw_eval(d, x, m)

Arguments

d

vector of long-memory parameters (dimension should match dimension of x).

x

data (matrix with time in rows and variables in columns).

m

truncation number used for the estimation of the periodogram.

Details

The choice of m determines the range of frequencies used in the computation of the periodogram, λj=2πj/N\lambda_j = 2\pi j/N, jj = 1,... , m. The optimal value depends on the spectral properties of the time series such as the presence of short range dependence. In Shimotsu (2007), m is chosen to be equal to N0.65N^{0.65}.

Value

multivariate Fourier Whittle estimator computed at point d.

Author(s)

S. Achard and I. Gannaz

References

K. Shimotsu (2007) Gaussian semiparametric estimation of multivariate fractionally integrated processes Journal of Econometrics Vol. 137, N. 2, pages 277-310.

S. Achard, I. Gannaz (2016) Multivariate wavelet Whittle estimation in long-range dependence. Journal of Time Series Analysis, Vol 37, N. 4, pages 476-512. http://arxiv.org/abs/1412.0391.

S. Achard, I Gannaz (2019) Wavelet-Based and Fourier-Based Multivariate Whittle Estimation: multiwave. Journal of Statistical Software, Vol 89, N. 6, pages 1-31.

See Also

mfw_cov_eval, mfw

Examples

### Simulation of ARFIMA(0,d,0)
rho <- 0.4
cov <- matrix(c(1,rho,rho,1),2,2)
d <- c(0.4,0.2)
J <- 9
N <- 2^J

resp <- fivarma(N, d, cov_matrix=cov)
x <- resp$x
long_run_cov <- resp$long_run_cov

m <- 57 ## default value of Shimotsu
res_mfw <- mfw(x,m)
d <- res_mfw$d
G <- mfw_eval(d,x,m)
k <- length(d)
res_d <- optim(rep(0,k),mfw_eval,x=x,m=m,method='Nelder-Mead',lower=-Inf,upper=Inf)$par

multivariate wavelet Whittle estimation

Description

Computes the multivariate wavelet Whittle estimation for the long-memory parameter vector d and the long-run covariance matrix, using DWTexact for the wavelet decomposition.

Usage

mww(x, filter, LU = NULL)

Arguments

x

data (matrix with time in rows and variables in columns).

filter

wavelet filter as obtain with scaling_filter.

LU

bivariate vector (optional) containing L, the lowest resolution in wavelet decomposition U, the maximal resolution in wavelet decomposition. (Default values are set to L=1, and U=Jmax.)

Details

L is fixing the lower limit of wavelet scales. L can be increased to avoid finest frequencies that can be corrupted by the presence of high frequency phenomena.

U is fixing the upper limit of wavelet scales. U can be decreased when highest frequencies have to be discarded.

Value

d

estimation of vector of long-memory parameters.

cov

estimation of long-run covariance matrix.

Author(s)

S. Achard and I. Gannaz

References

S. Achard, I. Gannaz (2016) Multivariate wavelet Whittle estimation in long-range dependence. Journal of Time Series Analysis, Vol 37, N. 4, pages 476-512. http://arxiv.org/abs/1412.0391.

S. Achard, I Gannaz (2019) Wavelet-Based and Fourier-Based Multivariate Whittle Estimation: multiwave. Journal of Statistical Software, Vol 89, N. 6, pages 1-31.

See Also

mww_eval, mww_cov_eval,mww_wav,mww_wav_eval,mww_wav_cov_eval

Examples

### Simulation of ARFIMA(0,d,0)
rho <- 0.4
cov <- matrix(c(1,rho,rho,1),2,2)
d <- c(0.4,0.2)
J <- 9
N <- 2^J

resp <- fivarma(N, d, cov_matrix=cov)
x <- resp$x
long_run_cov <- resp$long_run_cov

## wavelet coefficients definition
res_filter <- scaling_filter('Daubechies',8);
filter <- res_filter$h
M <- res_filter$M
alpha <- res_filter$alpha

LU <- c(2,11)

res_mww <- mww(x,filter,LU)

multivariate wavelet Whittle estimation of the long-run covariance matrix

Description

Computes the multivariate wavelet Whittle estimation of the long-run covariance matrix given the long-memory parameter vector d, using DWTexact for the wavelet decomposition.

Usage

mww_cov_eval(d, x, filter, LU)

Arguments

d

vector of long-memory parameters (dimension should match dimension of x).

x

data (matrix with time in rows and variables in columns).

filter

wavelet filter as obtain with scaling_filter.

LU

bivariate vector (optional) containing L, the lowest resolution in wavelet decomposition U, the maximal resolution in wavelet decomposition.

Details

L is fixing the lower limit of wavelet scales. L can be increased to avoid finest frequencies that can be corrupted by the presence of high frequency phenomena.

U is fixing the upper limit of wavelet scales. U can be decreased when highest frequencies have to be discarded.

Value

long-run covariance matrix estimation.

Author(s)

S. Achard and I. Gannaz

References

S. Achard, I. Gannaz (2016) Multivariate wavelet Whittle estimation in long-range dependence. Journal of Time Series Analysis, Vol 37, N. 4, pages 476-512. http://arxiv.org/abs/1412.0391.

S. Achard, I Gannaz (2019) Wavelet-Based and Fourier-Based Multivariate Whittle Estimation: multiwave. Journal of Statistical Software, Vol 89, N. 6, pages 1-31.

See Also

mww, mww_eval,mww_wav,mww_wav_eval,mww_wav_cov_eval

Examples

### Simulation of ARFIMA(0,d,0)
rho <- 0.4
cov <- matrix(c(1,rho,rho,1),2,2)
d <- c(0.4,0.2)
J <- 9
N <- 2^J

resp <- fivarma(N, d, cov_matrix=cov)
x <- resp$x
long_run_cov <- resp$long_run_cov

## wavelet coefficients definition
res_filter <- scaling_filter('Daubechies',8);
filter <- res_filter$h
M <- res_filter$M
alpha <- res_filter$alpha

LU <- c(2,11)

res_mww <- mww_cov_eval(d,x,filter,LU)

evaluation of multivariate wavelet Whittle estimation

Description

Evaluates the multivariate wavelet Whittle criterion at a given long-memory parameter vector d using DWTexact for the wavelet decomposition.

Usage

mww_eval(d, x, filter, LU = NULL)

Arguments

d

vector of long-memory parameters (dimension should match dimension of x).

x

data (matrix with time in rows and variables in columns).

filter

wavelet filter as obtain with scaling_filter.

LU

bivariate vector (optional) containing L, the lowest resolution in wavelet decomposition U, the maximal resolution in wavelet decomposition. (Default values are set to L=1, and U=Jmax.)

Details

L is fixing the lower limit of wavelet scales. L can be increased to avoid finest frequencies that can be corrupted by the presence of high frequency phenomena.

U is fixing the upper limit of wavelet scales. U can be decreased when highest frequencies have to be discarded.

Value

multivariate wavelet Whittle criterion.

Author(s)

S. Achard and I. Gannaz

References

E. Moulines, F. Roueff, M. S. Taqqu (2009) A wavelet whittle estimator of the memory parameter of a nonstationary Gaussian time series. Annals of statistics, vol. 36, N. 4, pages 1925-1956

S. Achard, I. Gannaz (2016) Multivariate wavelet Whittle estimation in long-range dependence. Journal of Time Series Analysis, Vol 37, N. 4, pages 476-512. http://arxiv.org/abs/1412.0391.

S. Achard, I Gannaz (2019) Wavelet-Based and Fourier-Based Multivariate Whittle Estimation: multiwave. Journal of Statistical Software, Vol 89, N. 6, pages 1-31.

See Also

mww, mww_cov_eval,mww_wav,mww_wav_eval,mww_wav_cov_eval

Examples

### Simulation of ARFIMA(0,d,0)
rho <- 0.4
cov <- matrix(c(1,rho,rho,1),2,2)
d <- c(0.4,0.2)
J <- 9
N <- 2^J

resp <- fivarma(N, d, cov_matrix=cov)
x <- resp$x
long_run_cov <- resp$long_run_cov

## wavelet coefficients definition
res_filter <- scaling_filter('Daubechies',8);
filter <- res_filter$h
M <- res_filter$M
alpha <- res_filter$alpha

LU <- c(2,11)

res_mww <- mww_eval(d,x,filter,LU)
k <- length(d)
res_d <- optim(rep(0,k),mww_eval,x=x,filter=filter,
	  	        LU=LU,method='Nelder-Mead',lower=-Inf,upper=Inf)$par

multivariate wavelet Whittle estimation for data as wavelet coefficients

Description

Computes the multivariate wavelet Whittle estimation of the long-memory parameter vector d and the long-run covariance matrix for the already wavelet decomposed data.

Usage

mww_wav(xwav, index, psih, grid_K, LU = NULL)

Arguments

xwav

wavelet coefficients matrix (with scales in rows and variables in columns).

index

vector containing the largest index of each band, i.e. for j>1j>1 the wavelet coefficients of scale jj are dwt(k)\code{dwt}(k) for k[indmaxband(j1)+1,indmaxband(j)]k \in [\code{indmaxband}(j-1)+1,\code{indmaxband}(j)] and for j=1j=1, dwt(k)\code{dwt}(k) for k[1,indmaxband(1)]k \in [1,\code{indmaxband}(1)].

psih

the Fourier transform of the wavelet mother at values grid_K.

grid_K

the grid for the approximation of the integral in K.

LU

bivariate vector (optional) containing L, the lowest resolution in wavelet decomposition U, the maximal resolution in wavelet decomposition. (Default values are set to L=1, and U=Jmax.)

Details

L is fixing the lower limit of wavelet scales. L can be increased to avoid finest frequencies that can be corrupted by the presence of high frequency phenomena.

U is fixing the upper limit of wavelet scales. U can be decreased when highest frequencies have to be discarded.

Value

d

estimation of the vector of long-memory parameters.

cov

estimation of the long-run covariance matrix.

Author(s)

S. Achard and I. Gannaz

References

S. Achard, I. Gannaz (2016) Multivariate wavelet Whittle estimation in long-range dependence. Journal of Time Series Analysis, Vol 37, N. 4, pages 476-512. http://arxiv.org/abs/1412.0391.

S. Achard, I Gannaz (2019) Wavelet-Based and Fourier-Based Multivariate Whittle Estimation: multiwave. Journal of Statistical Software, Vol 89, N. 6, pages 1-31.

See Also

mww_eval, mww_cov_eval,mww,mww_wav_eval,mww_wav_cov_eval

Examples

### Simulation of ARFIMA(0,d,0)
rho <- 0.4
cov <- matrix(c(1,rho,rho,1),2,2)
d <- c(0.4,0.2)
J <- 9
N <- 2^J

resp <- fivarma(N, d, cov_matrix=cov)
x <- resp$x
long_run_cov <- resp$long_run_cov

## wavelet coefficients definition
res_filter <- scaling_filter('Daubechies',8);
filter <- res_filter$h

LU <- c(2,11)

### wavelet decomposition

if(is.matrix(x)){
     N <- dim(x)[1]
     k <- dim(x)[2]
}else{
     N <- length(x)
     k <- 1
}
x <- as.matrix(x,dim=c(N,k))

     ## Wavelet decomposition
     xwav <- matrix(0,N,k)
     for(j in 1:k){
          xx <- x[,j]
             
          resw <- DWTexact(xx,filter)
          xwav_temp <- resw$dwt
          index <- resw$indmaxband
          Jmax <- resw$Jmax
          xwav[1:index[Jmax],j] <- xwav_temp;
     }
     ## we free some memory
     new_xwav <- matrix(0,min(index[Jmax],N),k)
     if(index[Jmax]<N){
          new_xwav[(1:(index[Jmax])),] <- xwav[(1:(index[Jmax])),]
     }
     xwav <- new_xwav
     index <- c(0,index)

##### Compute the wavelet functions 
res_psi <- psi_hat_exact(filter,10)
psih <- res_psi$psih
grid <- res_psi$grid

res_mww <- mww_wav(xwav,index, psih, grid,LU)

multivariate wavelet Whittle estimation of the long-run covariance matrix

Description

Computes the multivariate wavelet Whittle estimation of the long-run covariance matrix given the long-memory parameter vector d for the already wavelet decomposed data.

Usage

mww_wav_cov_eval(d, xwav, index,psih,grid_K, LU)

Arguments

d

vector of long-memory parameters (dimension should match dimension of xwav).

xwav

wavelet coefficients matrix (with scales in rows and variables in columns).

index

vector containing the largest index of each band, i.e. for j>1j>1 the wavelet coefficients of scale jj are dwt(k)\code{dwt}(k) for k[indmaxband(j1)+1,indmaxband(j)]k \in [\code{indmaxband}(j-1)+1,\code{indmaxband}(j)] and for j=1j=1, dwt(k)\code{dwt}(k) for k[1,indmaxband(1)]k \in [1,\code{indmaxband}(1)].

psih

the Fourier transform of the wavelet mother at values grid_K

grid_K

the grid for the approximation of the integral in K

LU

bivariate vector (optional) containing L, the lowest resolution in wavelet decomposition U, the maximal resolution in wavelet decomposition.

Details

L is fixing the lower limit of wavelet scales. L can be increased to avoid finest frequencies that can be corrupted by the presence of high frequency phenomena.

U is fixing the upper limit of wavelet scales. U can be decreased when highest frequencies have to be discarded.

Value

Long-run covariance matrix estimation.

Author(s)

S. Achard and I. Gannaz

References

S. Achard, I. Gannaz (2016) Multivariate wavelet Whittle estimation in long-range dependence. Journal of Time Series Analysis, Vol 37, N. 4, pages 476-512. http://arxiv.org/abs/1412.0391.

S. Achard, I Gannaz (2019) Wavelet-Based and Fourier-Based Multivariate Whittle Estimation: multiwave. Journal of Statistical Software, Vol 89, N. 6, pages 1-31.

See Also

mww, mww_eval,mww_wav,mww_wav_eval,mww_cov_eval

Examples

### Simulation of ARFIMA(0,d,0)
rho<-0.4
cov <- matrix(c(1,rho,rho,1),2,2)
d<-c(0.4,0.2)
J <- 9
N <- 2^J

resp <- fivarma(N, d, cov_matrix=cov)
x <- resp$x
long_run_cov <- resp$long_run_cov

## wavelet coefficients definition
res_filter <- scaling_filter('Daubechies',8);
filter <- res_filter$h
M <- res_filter$M
alpha <- res_filter$alpha

LU <- c(2,11)

### wavelet decomposition

if(is.matrix(x)){
     N <- dim(x)[1]
     k <- dim(x)[2]
}else{
     N <- length(x)
     k <- 1
}
x <- as.matrix(x,dim=c(N,k))

     ## Wavelet decomposition
     xwav <- matrix(0,N,k)
     for(j in 1:k){
          xx <- x[,j]
             
          resw <- DWTexact(xx,filter)
          xwav_temp <- resw$dwt
          index <- resw$indmaxband
          Jmax <- resw$Jmax
          xwav[1:index[Jmax],j] <- xwav_temp;
     }
     ## we free some memory
     new_xwav <- matrix(0,min(index[Jmax],N),k)
     if(index[Jmax]<N){
          new_xwav[(1:(index[Jmax])),] <- xwav[(1:(index[Jmax])),]
     }
     xwav <- new_xwav
     index <- c(0,index)

##### Compute the wavelet functions 
res_psi <- psi_hat_exact(filter,10)
psih<-res_psi$psih
grid<-res_psi$grid

res_mww <- mww_wav_cov_eval(d,xwav,index, psih, grid,LU)

multivariate wavelet Whittle estimation for data as wavelet coefficients

Description

Evaluates the multivariate wavelet Whittle criterion at a given long-memory parameter vector d for the already wavelet decomposed data.

Usage

mww_wav_eval(d, xwav, index, LU = NULL)

Arguments

d

vector of long-memory parameters (dimension should match dimension of x).

xwav

wavelet coefficients matrix (with scales in rows and variables in columns).

index

vector containing the largest index of each band, i.e. for j>1j>1 the wavelet coefficients of scale jj are dwt(k)\code{dwt}(k) for k[indmaxband(j1)+1,indmaxband(j)]k \in [\code{indmaxband}(j-1)+1,\code{indmaxband}(j)] and for j=1j=1, dwt(k)\code{dwt}(k) for k[1,indmaxband(1)]k \in [1,\code{indmaxband}(1)].

LU

bivariate vector (optional) containing L, the lowest resolution in wavelet decomposition U, the maximal resolution in wavelet decomposition. (Default values are set to L=1, and U=Jmax.)

Details

L is fixing the lower limit of wavelet scales. L can be increased to avoid finest frequencies that can be corrupted by the presence of high frequency phenomena.

U is fixing the upper limit of wavelet scales. U can be decreased when highest frequencies have to be discarded.

Value

multivariate wavelet Whittle criterion.

Author(s)

S. Achard and I. Gannaz

References

E. Moulines, F. Roueff, M. S. Taqqu (2009) A wavelet whittle estimator of the memory parameter of a nonstationary Gaussian time series. Annals of statistics, vol. 36, N. 4, pages 1925-1956

S. Achard, I. Gannaz (2016) Multivariate wavelet Whittle estimation in long-range dependence. Journal of Time Series Analysis, Vol 37, N. 4, pages 476-512. http://arxiv.org/abs/1412.0391.

S. Achard, I Gannaz (2019) Wavelet-Based and Fourier-Based Multivariate Whittle Estimation: multiwave. Journal of Statistical Software, Vol 89, N. 6, pages 1-31.

See Also

mww, mww_cov_eval,mww_wav,mww_eval,mww_wav_cov_eval

Examples

### Simulation of ARFIMA(0,d,0)
rho <- 0.4
cov <- matrix(c(1,rho,rho,1),2,2)
d <- c(0.4,0.2)
J <- 9
N <- 2^J

resp <- fivarma(N, d, cov_matrix=cov)
x <- resp$x
long_run_cov <- resp$long_run_cov

## wavelet coefficients definition
res_filter <- scaling_filter('Daubechies',8);
filter <- res_filter$h

LU <- c(2,11)

### wavelet decomposition

if(is.matrix(x)){
     N <- dim(x)[1]
     k <- dim(x)[2]
}else{
     N <- length(x)
     k <- 1
}
x <- as.matrix(x,dim=c(N,k))

     ## Wavelet decomposition
     xwav <- matrix(0,N,k)
     for(j in 1:k){
          xx <- x[,j]
             
          resw <- DWTexact(xx,filter)
          xwav_temp <- resw$dwt
          index <- resw$indmaxband
          Jmax <- resw$Jmax
          xwav[1:index[Jmax],j] <- xwav_temp;
     }
     ## we free some memory
     new_xwav <- matrix(0,min(index[Jmax],N),k)
     if(index[Jmax]<N){
          new_xwav[(1:(index[Jmax])),] <- xwav[(1:(index[Jmax])),]
     }
     xwav <- new_xwav
     index <- c(0,index)

res_mww <- mww_wav_eval(d,xwav,index,LU)
res_d <- optim(rep(0,k),mww_wav_eval,xwav=xwav,index=index,LU=LU,
          method='Nelder-Mead',lower=-Inf,upper=Inf)$par

discrete Fourier transform of the wavelet

Description

Computes the discrete Fourier transform of the wavelet associated to the given filter using scaling_function. The length of the Fourier transform is equal to the length of the grid where the wavelet is evaluated.

Usage

psi_hat_exact(filter,J=10)

Arguments

filter

wavelet filter as obtained with scaling_filter.

J

2^J corresponds to the size of the grid for the discretisation of the wavelet. The default value is set to 10.

Value

psih

Values of the discrete Fourier transform of the wavelet.

grid

Frequencies where the Fourier transform is evaluated.

Author(s)

S. Achard and I. Gannaz

References

G. Fay, E. Moulines, F. Roueff, M. S. Taqqu (2009) Estimators of long-memory: Fourier versus wavelets. Journal of Econometrics, vol. 151, N. 2, pages 159-177.

S. Achard, I. Gannaz (2016) Multivariate wavelet Whittle estimation in long-range dependence. Journal of Time Series Analysis, Vol 37, N. 4, pages 476-512. http://arxiv.org/abs/1412.0391.

S. Achard, I Gannaz (2019) Wavelet-Based and Fourier-Based Multivariate Whittle Estimation: multiwave. Journal of Statistical Software, Vol 89, N. 6, pages 1-31.

See Also

DWTexact, scaling_filter

Examples

res_filter <- scaling_filter('Daubechies',8);
filter <- res_filter$h
psi_hat_exact(filter,J=6)

wavelet scaling filter coefficients

Description

Computes the filter coefficients of the Haar or Daubechies wavelet family with a specific order

Usage

scaling_filter(family, parameter)

Arguments

family

Wavelet family, 'Haar' or 'Daubechies'

parameter

Order of the Daubechies wavelet (equal to twice the number of vanishing moments). The value of parameter can be 2,4,8,10,12,14 and 16.

Value

h

Vector of scaling filter coefficients.

M

Number of vanishing moments.

alpha

Fourier decay exponent.

Author(s)

S. Achard and I. Gannaz

References

G. Fay, E. Moulines, F. Roueff, M. S. Taqqu (2009) Estimators of long-memory: Fourier versus wavelets. Journal of Econometrics, vol. 151, N. 2, pages 159-177.

S. Achard, I. Gannaz (2016) Multivariate wavelet Whittle estimation in long-range dependence. Journal of Time Series Analysis, Vol 37, N. 4, pages 476-512. http://arxiv.org/abs/1412.0391.

S. Achard, I Gannaz (2019) Wavelet-Based and Fourier-Based Multivariate Whittle Estimation: multiwave. Journal of Statistical Software, Vol 89, N. 6, pages 1-31.

See Also

DWTexact

Examples

res_filter <- scaling_filter('Daubechies',8);
filter <- res_filter$h
M <- res_filter$M
alpha <- res_filter$alpha

scaling function and the wavelet function

Description

Computes the scaling function and the wavelet function (for compactly supported wavelet) using the cascade algorithm on the grid of dyadic integer 2J2^{-J}

Usage

scaling_function(filter,J)

Arguments

filter

wavelet filter as obtained with scaling_filter.

J

value of the largest scale.

Value

phi

Scaling function.

psi

Wavelet function.

Note

This function was rewritten from an original matlab version by Fay et al. (2009)

Author(s)

S. Achard and I. Gannaz

References

G. Fay, E. Moulines, F. Roueff, M. S. Taqqu (2009) Estimators of long-memory: Fourier versus wavelets. Journal of Econometrics, vol. 151, N. 2, pages 159-177.

S. Achard, I. Gannaz (2016) Multivariate wavelet Whittle estimation in long-range dependence. Journal of Time Series Analysis, Vol 37, N. 4, pages 476-512. http://arxiv.org/abs/1412.0391.

S. Achard, I Gannaz (2019) Wavelet-Based and Fourier-Based Multivariate Whittle Estimation: multiwave. Journal of Statistical Software, Vol 89, N. 6, pages 1-31.

See Also

DWTexact, scaling_filter

Examples

res_filter <- scaling_filter('Daubechies',8);
filter <- res_filter$h
scaling_function(filter,J=6)

Transform a vector in a non symmetric Toeplitz matrix

Description

Transform a vector in a non symmetric Toeplitz matrix

Usage

toeplitz_nonsym(vec)

Arguments

vec

input vector.

Value

the corresponding matrix.

Author(s)

S. Achard and I. Gannaz

References

S. Achard, I. Gannaz (2016) Multivariate wavelet Whittle estimation in long-range dependence. Journal of Time Series Analysis, Vol 37, N. 4, pages 476-512. http://arxiv.org/abs/1412.0391.

S. Achard, I Gannaz (2019) Wavelet-Based and Fourier-Based Multivariate Whittle Estimation: multiwave. Journal of Statistical Software, Vol 89, N. 6, pages 1-31.

See Also

scaling_function

Examples

res_filter <- scaling_filter('Daubechies',8);
filter <- res_filter$h
Htmp <- toeplitz_nonsym(filter)

simulation of multivariate ARMA process

Description

generates N observations of a k-vector ARMA process

Usage

varma(N, k = 1, VAR = NULL, VMA = NULL, cov_matrix = diag(k), innov=NULL)

Arguments

N

number of time points.

k

dimension of the vector ARMA (optional, default is univariate)

VAR

array of VAR coefficient matrices (optional).

VMA

array of VMA coefficient matrices (optional).

cov_matrix

matrix of correlation between the innovations (optional, default is identity).

innov

matrix of the innovations (optional, default is a gaussian process).

Value

vector containing the N observations of the k-vector ARMA process.

Author(s)

S. Achard and I. Gannaz

References

S. Achard, I. Gannaz (2016) Multivariate wavelet Whittle estimation in long-range dependence. Journal of Time Series Analysis, Vol 37, N. 4, pages 476-512. http://arxiv.org/abs/1412.0391.

S. Achard, I Gannaz (2019) Wavelet-Based and Fourier-Based Multivariate Whittle Estimation: multiwave. Journal of Statistical Software, Vol 89, N. 6, pages 1-31.

See Also

fivarma, vfracdiff

Examples

rho1 <- 0.3
rho2 <- 0.8
cov <- matrix(c(1,rho1,rho2,rho1,1,rho1,rho2,rho1,1),3,3)

J <- 9
N <- 2^J
VMA <- diag(c(0.4,0.1,0))
### or another example VAR <- array(c(0.8,0,0,0,0.6,0,0,0,0.2,0,0,0,0,0.4,0,0,0,0.5),dim=c(3,3,2))
VAR <- diag(c(0.8,0.6,0))
x <- varma(N, k=3, cov_matrix=cov, VAR=VAR, VMA=VMA)

simulation of vector fractional differencing process

Description

Given a vector process x and a vector of long memory parameters d, this function is producing the corresponding fractional differencing process.

Usage

vfracdiff(x, d)

Arguments

x

initial process.

d

vector of long-memory parameters

Details

Given a process x, this function applied a fractional difference procedure using the formula:

diag((1L)d)x,diag((1-L)^d) x,

where L is the lag operator.

Value

vector fractional differencing of x.

Author(s)

S. Achard and I. Gannaz

References

S. Achard, I. Gannaz (2016) Multivariate wavelet Whittle estimation in long-range dependence. Journal of Time Series Analysis, Vol 37, N. 4, pages 476-512. http://arxiv.org/abs/1412.0391.

K. Shimotsu (2007) Gaussian semiparametric estimation of multivariate fractionally integrated processes Journal of Econometrics Vol. 137, N. 2, pages 277-310.

S. Achard, I Gannaz (2019) Wavelet-Based and Fourier-Based Multivariate Whittle Estimation: multiwave. Journal of Statistical Software, Vol 89, N. 6, pages 1-31.

See Also

varma, fivarma

Examples

rho1 <- 0.3
rho2 <- 0.8
cov <- matrix(c(1,rho1,rho2,rho1,1,rho1,rho2,rho1,1),3,3)
d <- c(0.2,0.3,0.4)


J <- 9
N <- 2^J
VMA <- diag(c(0.4,0.1,0))
### or another example VAR <- array(c(0.8,0,0,0,0.6,0,0,0,0.2,0,0,0,0,0.4,0,0,0,0.5),dim=c(3,3,2))
VAR <- diag(c(0.8,0.6,0))
x <- varma(N, k=3, cov_matrix=cov, VAR=VAR, VMA=VMA)
vx<-vfracdiff(x,d)