Package 'iAR'

Title: Irregularly Observed Autoregressive Models
Description: Data sets, functions and scripts with examples to implement autoregressive models for irregularly observed time series. The models available in this package are the irregular autoregressive model (Eyheramendy et al.(2018) <doi:10.1093/mnras/sty2487>), the complex irregular autoregressive model (Elorrieta et al.(2019) <doi:10.1051/0004-6361/201935560>) and the bivariate irregular autoregressive model (Elorrieta et al.(2021) <doi:10.1093/mnras/stab1216>).
Authors: Elorrieta Felipe [aut, cre], Ojeda Cesar [aut], Eyheramendy Susana [aut], Palma Wilfredo [aut]
Maintainer: Elorrieta Felipe <[email protected]>
License: GPL-2
Version: 1.2.0
Built: 2024-11-18 06:53:31 UTC
Source: CRAN

Help Index


Active Galactic Nuclei

Description

Time series of the AGN MCG-6-30-15 measured in the K-band between 2006 August and 2011 July with the ANDICAM camera mounted on the 1.3 m telescope at Cerro Tololo Inter-American Observatory (CTIO)

Usage

agn

Format

A data frame with 237 observations on the following 3 variables:

t

heliocentric Julian Day - 2450000

m

Flux $(10^(-15) ergs/s/cm^2 /A)$

merr

measurement error standard deviations.

References

Lira P, Arévalo P, Uttley P, McHardy IMM, Videla L (2015). “Long-term monitoring of the archetype Seyfert galaxy MCG-6-30-15: X-ray, optical and near-IR variability of the corona, disc and torus.” Monthly Notices of the Royal Astronomical Society, 454(1), 368-379. ISSN 0035-8711, doi:10.1093/mnras/stv1945.

Examples

data(agn)
plot(agn$t,agn$m,type="l",ylab="",xlab="")

Fitted Values of BIAR model

Description

Fit a BIAR model to a bivariate irregularly observed time series.

Usage

BIARfit(phiValues, y1, y2, t, yerr1, yerr2, zeroMean = TRUE)

Arguments

phiValues

An array with the parameters of the BIAR model. The elements of the array are, in order, the autocorrelation and the cross correlation parameter of the BIAR model.

y1

Array with the observations of the first time series of the BIAR process.

y2

Array with the observations of the second time series of the BIAR process.

t

Array with the irregular observational times.

yerr1

Array with the measurements error standard deviations of the first time series of the BIAR process.

yerr2

Array with the measurements error standard deviations of the second time series of the BIAR process.

zeroMean

logical; if TRUE, the array y has zero mean; if FALSE, y has a mean different from zero.

Value

A list with the following components:

  • rho Estimated value of the contemporary correlation coefficient.

  • innov.var Estimated value of the innovation variance.

  • fitted Fitted values of the BIAR model.

  • fitted.state Fitted state values of the BIAR model.

  • Lambda Lambda value estimated by the BIAR model at the last time point.

  • Theta Theta array estimated by the BIAR model at the last time point.

  • Sighat Covariance matrix estimated by the BIAR model at the last time point.

  • Qt Covariance matrix of the state equation estimated by the BIAR model at the last time point.

References

Elorrieta F, Eyheramendy S, Palma W, Ojeda C (2021). “A novel bivariate autoregressive model for predicting and forecasting irregularly observed time series.” Monthly Notices of the Royal Astronomical Society, 505(1), 1105-1116. ISSN 0035-8711, doi:10.1093/mnras/stab1216, https://academic.oup.com/mnras/article-pdf/505/1/1105/38391762/stab1216.pdf.

See Also

gentime, BIARsample, BIARphikalman, BIARkalman

Examples

n=80
set.seed(6714)
st<-gentime(n)
x=BIARsample(n=n,phiR=0.9,phiI=0.3,st=st,rho=0.9)
y=x$y
y1=y/apply(y,1,sd)
yerr1=rep(0,n)
yerr2=rep(0,n)
biar=BIARkalman(y1=y1[1,],y2=y1[2,],t=st,delta1 = yerr1,delta2=yerr2)
biar
predbiar=BIARfit(phiValues=c(biar$phiR,biar$phiI),y1=y1[1,],y2=y1[2,],t=st,yerr1
 = rep(0,length(y[1,])),yerr2=rep(0,length(y[1,])))
rho=predbiar$rho
print(rho)
yhat=predbiar$fitted

Forecast from BIAR model

Description

Forecast from models fitted by BIARkalman

Usage

BIARforecast(phiR, phiI, y1, y2, t, tAhead)

Arguments

phiR

Autocorrelation coefficient of BIAR model.

phiI

Cross-correlation coefficient of BIAR model.

y1

Array with the observations of the first time series of the BIAR process.

y2

Array with the observations of the second time series of the BIAR process.

t

Array with the observational times.

tAhead

The time ahead for which the forecast is required.

Value

A list with the following components:

  • fitted Fitted values by the BIAR model.

  • forecast Point forecast in the time ahead required.

  • Lambda Lambda value estimated by the BIAR model at the last time point.

  • Sighat Covariance matrix estimated by the BIAR model at the last time point.

References

Elorrieta F, Eyheramendy S, Palma W, Ojeda C (2021). “A novel bivariate autoregressive model for predicting and forecasting irregularly observed time series.” Monthly Notices of the Royal Astronomical Society, 505(1), 1105-1116. ISSN 0035-8711, doi:10.1093/mnras/stab1216, https://academic.oup.com/mnras/article-pdf/505/1/1105/38391762/stab1216.pdf.

See Also

BIARsample, BIARkalman, BIARfit

Examples

#Simulated Data
n=100
set.seed(6714)
st<-gentime(n)
x=BIARsample(n=n,phiR=0.9,phiI=0.3,st=st)
biar=iAR::BIARkalman(y1=x$y[1,],y2=x$y[2,],t=st)
forBIAR<-BIARforecast(phiR=biar$phiR,phiI=biar$phiI,y1=x$y[1,],y2=x$y[2,],t=st,tAhead=c(1.3))

Interpolation from BIAR model

Description

Interpolation of missing values from models fitted by BIARkalman

Usage

BIARinterpolation(
  x,
  y1,
  y2,
  t,
  delta1 = 0,
  delta2 = 0,
  yini1 = 0,
  yini2 = 0,
  zero.mean = TRUE,
  niter = 10,
  seed = 1234,
  nsmooth = 1
)

Arguments

x

An array with the parameters of the BIAR model. The elements of the array are, in order, the real (phiR) and the imaginary (phiI) part of the coefficient of BIAR model.

y1

Array with the observations of the first time series of the BIAR process.

y2

Array with the observations of the second time series of the BIAR process.

t

Array with the irregular observational times.

delta1

Array with the measurements error standard deviations of the first time series of the BIAR process.

delta2

Array with the measurements error standard deviations of the second time series of the BIAR process.

yini1

a single value, initial value of the estimation of the missing value of the first time series of the BIAR process.

yini2

a single value, initial value of the estimation of the missing value of the second time series of the BIAR process.

zero.mean

logical; if TRUE, the array y has zero mean; if FALSE, y has a mean different from zero.

niter

Number of iterations in which the function nlminb will be repeated.

seed

a single value, interpreted as the seed of the random process.

nsmooth

a single value; If 1, only one time series of the BIAR process has a missing value. If 2, both time series of the BIAR process have a missing value.

Value

A list with the following components:

  • fitted Estimation of the missing values of the BIAR process.

  • ll Value of the negative log likelihood evaluated in the fitted missing values.

References

Elorrieta F, Eyheramendy S, Palma W, Ojeda C (2021). “A novel bivariate autoregressive model for predicting and forecasting irregularly observed time series.” Monthly Notices of the Royal Astronomical Society, 505(1), 1105-1116. ISSN 0035-8711, doi:10.1093/mnras/stab1216, https://academic.oup.com/mnras/article-pdf/505/1/1105/38391762/stab1216.pdf.

See Also

gentime, BIARsample, BIARphikalman

Examples

set.seed(6713)
n=100
st<-gentime(n)
x=BIARsample(n=n,phiR=0.9,phiI=0.3,st=st,rho=0.9)
y=x$y
y1=y/apply(y,1,sd)
yerr1=rep(0,n)
yerr2=rep(0,n)
biar=BIARkalman(y1=y1[1,],y2=y1[2,],t=st,delta1 = yerr1,delta2=yerr2)
biar
napos=10
y0=y1
y1[1,napos]=NA
xest=c(biar$phiR,biar$phiI)
yest=BIARinterpolation(xest,y1=y1[1,],y2=y1[2,],t=st,delta1=yerr1,
delta2=yerr2,nsmooth=1)
yest$fitted
mse=(y0[1,napos]-yest$fitted)^2
print(mse)
par(mfrow=c(2,1))
plot(st,x$y[1,],type='l',xlim=c(st[napos-5],st[napos+5]))
points(st,x$y[1,],pch=20)
points(st[napos],yest$fitted*apply(y,1,sd)[1],col="red",pch=20)
plot(st,x$y[2,],type='l',xlim=c(st[napos-5],st[napos+5]))
points(st,x$y[2,],pch=20)

Maximum Likelihood Estimation of the BIAR Model via Kalman Recursions

Description

Maximum Likelihood Estimation of the BIAR model parameters phiR and phiI. The estimation procedure uses the Kalman Filter to find the maximum of the likelihood.

Usage

BIARkalman(
  y1,
  y2,
  t,
  delta1 = 0,
  delta2 = 0,
  zero.mean = "TRUE",
  niter = 10,
  seed = 1234
)

Arguments

y1

Array with the observations of the first time series of the BIAR process.

y2

Array with the observations of the second time series of the BIAR process.

t

Array with the irregular observational times.

delta1

Array with the measurements error standard deviations of the first time series of the BIAR process.

delta2

Array with the measurements error standard deviations of the second time series of the BIAR process.

zero.mean

logical; if true, the array y has zero mean; if false, y has a mean different from zero.

niter

Number of iterations in which the function nlminb will be repeated.

seed

a single value, interpreted as the seed of the random process.

Value

A list with the following components:

  • phiR MLE of the autocorrelation coefficient of BIAR model (phiR).

  • phiI MLE of the cross-correlation coefficient of the BIAR model (phiI).

  • ll Value of the negative log likelihood evaluated in phiR and phiI.

References

Elorrieta F, Eyheramendy S, Palma W, Ojeda C (2021). “A novel bivariate autoregressive model for predicting and forecasting irregularly observed time series.” Monthly Notices of the Royal Astronomical Society, 505(1), 1105-1116. ISSN 0035-8711, doi:10.1093/mnras/stab1216, https://academic.oup.com/mnras/article-pdf/505/1/1105/38391762/stab1216.pdf.

See Also

gentime, BIARsample, BIARphikalman

Examples

n=80
set.seed(6714)
st<-gentime(n)
x=BIARsample(n=n,phiR=0.9,phiI=0,st=st,rho=0)
y=x$y
y1=y/apply(y,1,sd)
biar=BIARkalman(y1=y1[1,],y2=y1[2,],t=st,delta1 = rep(0,length(y[1,])),
delta2=rep(0,length(y[1,])))
biar

Minus Log Likelihood of the BIAR Model

Description

This function return the negative log likelihood of the BIAR process given specific values of phiR and phiI

Usage

BIARphikalman(yest, phiValues, y1, y2, t, yerr1, yerr2, zeroMean = TRUE)

Arguments

yest

An array with the estimate of a missing value in one or both time series of the bivariate process. This function recognizes a missing value with a NA. If the bivariate time series does not have a missing value, this value does not affect the computation of the likelihood.

phiValues

An array with the parameters of the BIAR model. The elements of the array are, in order, the real (phiR) and the imaginary (phiI) part of the coefficient of BIAR model.

y1

Array with the observations of the first time series of the BIAR process.

y2

Array with the observations of the second time series of the BIAR process.

t

Array with the irregular observational times.

yerr1

Array with the measurements error standard deviations of the first time series of the BIAR process.

yerr2

Array with the measurements error standard deviations of the second time series of the BIAR process.

zeroMean

logical; if TRUE, the array y has zero mean; if FALSE, y has a mean different from zero.

Value

Value of the negative log likelihood evaluated in phiR and phiI.

References

Elorrieta F, Eyheramendy S, Palma W, Ojeda C (2021). “A novel bivariate autoregressive model for predicting and forecasting irregularly observed time series.” Monthly Notices of the Royal Astronomical Society, 505(1), 1105-1116. ISSN 0035-8711, doi:10.1093/mnras/stab1216, https://academic.oup.com/mnras/article-pdf/505/1/1105/38391762/stab1216.pdf.

See Also

gentime, BIARsample

Examples

n=300
set.seed(6714)
st<-gentime(n)
x=BIARsample(n=n,phiR=0.9,phiI=0.3,st=st)
y=x$y
y1=y[1,]
y2=y[2,]
yerr1=rep(0,n)
yerr2=rep(0,n)
BIARphikalman(phiValues=c(0.8,0.2),y1=y1,y2=y2,t=st,yerr1=yerr1,yerr2=yerr2,yest=c(0,0))

Simulate from a BIAR Model

Description

Simulates a BIAR Time Series Model

Usage

BIARsample(n, st, phiR, phiI, delta1 = 0, delta2 = 0, rho = 0)

Arguments

n

Length of the output bivariate time series. A strictly positive integer.

st

Array with observational times.

phiR

Autocorrelation coefficient of BIAR model. A value between -1 and 1.

phiI

Crosscorrelation coefficient of BIAR model. A value between -1 and 1.

delta1

Array with the measurements error standard deviations of the first time series of the bivariate process.

delta2

Array with the measurements error standard deviations of the second time series of the bivariate process.

rho

Contemporary correlation coefficient of BIAR model. A value between -1 and 1.

Details

The chosen phiR and phiI values must satisfy the condition $|phiR + i phiI| < 1$.

Value

A list with the following components:

  • y Matrix with the simulated BIAR process.

  • t Array with observation times.

  • Sigma Covariance matrix of the process.

References

Elorrieta F, Eyheramendy S, Palma W, Ojeda C (2021). “A novel bivariate autoregressive model for predicting and forecasting irregularly observed time series.” Monthly Notices of the Royal Astronomical Society, 505(1), 1105-1116. ISSN 0035-8711, doi:10.1093/mnras/stab1216, https://academic.oup.com/mnras/article-pdf/505/1/1105/38391762/stab1216.pdf.

See Also

gentime

Examples

n=300
set.seed(6714)
st<-gentime(n)
x=BIARsample(n=n,phiR=0.9,phiI=0.3,st=st)
plot(st,x$y[1,],type='l')
plot(st,x$y[2,],type='l')
x=BIARsample(n=n,phiR=-0.9,phiI=-0.3,st=st)
plot(st,x$y[1,],type='l')
plot(st,x$y[2,],type='l')

Fitted Values of CIAR model

Description

Fit a CIAR model to an irregularly observed time series.

Usage

CIARfit(phiValues, y, t, standardized = TRUE, c = 1)

Arguments

phiValues

An array with the parameters of the CIAR model. The elements of the array are, in order, the real and the imaginary part of the phi parameter of the CIAR model.

y

Array with the time series observations.

t

Array with the irregular observational times.

standardized

logical; if TRUE, the array y is standardized; if FALSE, y contains the raw time series

c

Nuisance parameter corresponding to the variance of the imaginary part.

Value

A list with the following components:

  • yhat Fitted values of the observable part of CIAR model.

  • xhat Fitted values of both observable part and imaginary part of CIAR model.

  • Lambda Lambda value estimated by the CIAR model at the last time point.

  • Theta Theta array estimated by the CIAR model at the last time point.

  • Sighat Covariance matrix estimated by the CIAR model at the last time point.

  • Qt Covariance matrix of the state equation estimated by the CIAR model at the last time point.

References

Elorrieta, F, Eyheramendy, S, Palma, W (2019). “Discrete-time autoregressive model for unequally spaced time-series observations.” A&A, 627, A120. doi:10.1051/0004-6361/201935560.

See Also

gentime, CIARsample, CIARphikalman,CIARkalman

Examples

n=100
set.seed(6714)
st<-gentime(n)
x=CIARsample(n=n,phiR=0.9,phiI=0,st=st,c=1)
y=x$y
y1=y/sd(y)
ciar=CIARkalman(y=y1,t=st)
ciar
yhat=CIARfit(phiValues=c(ciar$phiR,ciar$phiI),y=y1,t=st)

Forecast from CIAR model

Description

Forecast from models fitted by CIARkalman

Usage

CIARforecast(phiR, phiI, y1, st, tAhead)

Arguments

phiR

Real part of the phi coefficient of CIAR model.

phiI

Imaginary part of the phi coefficient of CIAR model.

y1

Array with the time series observations.

st

Array with the observational times.

tAhead

The time ahead for which the forecast is required.

Value

A list with the following components:

  • fitted Fitted values by the CIAR model.

  • forecast Point forecast in the time ahead required.

  • Lambda Lambda value estimated by the CIAR model at the last time point.

  • Sighat Covariance matrix estimated by the CIAR model at the last time point.

References

Elorrieta, F, Eyheramendy, S, Palma, W (2019). “Discrete-time autoregressive model for unequally spaced time-series observations.” A&A, 627, A120. doi:10.1051/0004-6361/201935560.

See Also

CIARsample, CIARkalman, CIARfit

Examples

#Simulated Data
n=100
set.seed(6714)
st<-gentime(n)
x=CIARsample(n=n,phiR=0.9,phiI=0,st=st,c=1)
y=x$y
y1=y/sd(y)
n=length(y1)
p=trunc(n*0.99)
ytr=y1[1:p]
yte=y1[(p+1):n]
str=st[1:p]
ste=st[(p+1):n]
tahead=ste-str[p]

ciar=CIARkalman(y=ytr,t=str)
forCIAR<-CIARforecast(ciar$phiR,ciar$phiI,ytr,str,tAhead=tahead)

Interpolation from CIAR model

Description

Interpolation of missing values from models fitted by CIARkalman

Usage

CIARinterpolation(
  x,
  y,
  t,
  delta = 0,
  yini = 0,
  zero.mean = TRUE,
  standardized = TRUE,
  c = 1,
  seed = 1234
)

Arguments

x

An array with the parameters of the CIAR model. The elements of the array are, in order, the real (phiR) and the imaginary (phiI) part of the coefficient of CIAR model.

y

Array with the time series observations.

t

Array with the irregular observational times.

delta

Array with the measurements error standard deviations.

yini

a single value, initial value for the estimation of the missing value of the time series.

zero.mean

logical; if TRUE, the array y has zero mean; if FALSE, y has a mean different from zero.

standardized

logical; if TRUE, the array y is standardized; if FALSE, y contains the raw time series.

c

Nuisance parameter corresponding to the variance of the imaginary part.

seed

a single value, interpreted as the seed of the random process.

Value

A list with the following components:

  • fitted Estimation of a missing value of the CIAR process.

  • ll Value of the negative log likelihood evaluated in the fitted missing values.

References

Elorrieta, F, Eyheramendy, S, Palma, W (2019). “Discrete-time autoregressive model for unequally spaced time-series observations.” A&A, 627, A120. doi:10.1051/0004-6361/201935560.

See Also

gentime, CIARsample, CIARkalman

Examples

n=100
set.seed(6714)
st<-gentime(n)
x=CIARsample(n=n,phiR=0.9,phiI=0,st=st,c=1)
y=x$y
y1=y/sd(y)
ciar=CIARkalman(y=y1,t=st)
ciar
napos=10
y0=y1
y1[napos]=NA
xest=c(ciar$phiR,ciar$phiI)
yest=CIARinterpolation(xest,y=y1,t=st)
yest$fitted
mse=(y0[napos]-yest$fitted)^2
print(mse)
plot(st,y,type='l',xlim=c(st[napos-5],st[napos+5]))
points(st,y,pch=20)
points(st[napos],yest$fitted*sd(y),col="red",pch=20)

Maximum Likelihood Estimation of the CIAR Model via Kalman Recursions

Description

Maximum Likelihood Estimation of the CIAR model parameters phiR and phiI. The estimation procedure uses the Kalman Filter to find the maximum of the likelihood.

Usage

CIARkalman(
  y,
  t,
  delta = 0,
  zero.mean = TRUE,
  standardized = TRUE,
  c = 1,
  niter = 10,
  seed = 1234
)

Arguments

y

Array with the time series observations.

t

Array with the irregular observational times.

delta

Array with the measurements error standard deviations.

zero.mean

logical; if TRUE, the array y has zero mean; if FALSE, y has a mean different from zero.

standardized

logical; if TRUE, the array y is standardized; if FALSE, y contains the raw time series.

c

Nuisance parameter corresponding to the variance of the imaginary part.

niter

Number of iterations in which the function nlminb will be repeated.

seed

a single value, interpreted as the seed of the random process.

Value

A list with the following components:

  • phiR MLE of the Real part of the coefficient of CIAR model (phiR).

  • phiI MLE of the Imaginary part of the coefficient of the CIAR model (phiI).

  • ll Value of the negative log likelihood evaluated in phiR and phiI.

References

Elorrieta, F, Eyheramendy, S, Palma, W (2019). “Discrete-time autoregressive model for unequally spaced time-series observations.” A&A, 627, A120. doi:10.1051/0004-6361/201935560.

See Also

gentime, CIARsample, CIARphikalman

Examples

n=100
set.seed(6714)
st<-gentime(n)
x=CIARsample(n=n,phiR=0.9,phiI=0,st=st,c=1)
y=x$y
y1=y/sd(y)
ciar=CIARkalman(y=y1,t=st)
ciar
Mod(complex(real=ciar$phiR,imaginary=ciar$phiI))

Minus Log Likelihood of the CIAR Model

Description

This function return the negative log likelihood of the CIAR process given specific values of phiR and phiI

Usage

CIARphikalman(yest, x, y, t, yerr, zeroMean = TRUE, standardized = TRUE, c = 1)

Arguments

yest

The estimate of a missing value in the time series. This function recognizes a missing value with a NA. If the time series does not have a missing value, this value does not affect the computation of the likelihood.

x

An array with the parameters of the CIAR model. The elements of the array are, in order, the real (phiR) and the imaginary (phiI) part of the coefficient of CIAR model.

y

Array with the time series observations.

t

Array with the irregular observational times.

yerr

Array with the measurements error standard deviations.

zeroMean

logical; if TRUE, the array y has zero mean; if FALSE, y has a mean different from zero.

standardized

logical; if TRUE, the array y is standardized; if FALSE, y contains the raw time series.

c

Nuisance parameter corresponding to the variance of the imaginary part.

Value

Value of the negative log likelihood evaluated in phiR and phiI.

References

Elorrieta, F, Eyheramendy, S, Palma, W (2019). “Discrete-time autoregressive model for unequally spaced time-series observations.” A&A, 627, A120. doi:10.1051/0004-6361/201935560.

See Also

gentime, CIARsample

Examples

n=300
set.seed(6714)
st<-gentime(n)
x=CIARsample(n=n,phiR=0.9,phiI=0,st=st,c=1)
y=x$y
yerr=rep(0,n)
CIARphikalman(x=c(0.8,0),y=y,t=st,yerr=yerr,yest=0)

Simulate from a CIAR Model

Description

Simulates a CIAR Time Series Model

Usage

CIARsample(n, phiR, phiI, st, rho = 0L, c = 1L)

Arguments

n

Length of the output time series. A strictly positive integer.

phiR

Real part of the coefficient of CIAR model. A value between -1 and 1.

phiI

Imaginary part of the coefficient of CIAR model. A value between -1 and 1.

st

Array with observational times.

rho

Correlation between the real and the imaginary part of the process. A value between -1 and 1.

c

Nuisance parameter corresponding to the variance of the imaginary part.

Details

The chosen phiR and phiI values must satisfy the condition $|phiR + i phiI| < 1$.

Value

A list with the following components:

  • yArray with the simulated real part of the CIAR process.

  • t Array with observation times.

  • Sigma Covariance matrix of the process.

References

Elorrieta, F, Eyheramendy, S, Palma, W (2019). “Discrete-time autoregressive model for unequally spaced time-series observations.” A&A, 627, A120. doi:10.1051/0004-6361/201935560.

See Also

gentime

Examples

n=300
set.seed(6714)
st<-gentime(n)
x=CIARsample(n=n,phiR=0.9,phiI=0,st=st,c=1)
plot(st,x$y,type='l')
x=CIARsample(n=n,phiR=-0.9,phiI=0,st=st,c=1)
plot(st,x$y,type='l')

Classical Cepheid

Description

Time series of a classical cepheid variable star obtained from HIPPARCOS.

Usage

clcep

Format

A data frame with 109 observations on the following 3 variables:

t

heliocentric Julian Day

m

magnitude

merr

measurement error of the magnitude (in mag).

Details

The frequency computed by GLS for this light curve is 0.060033386. Catalogs and designations of this star: HD 1989: HD 305996 TYCHO-2 2000:TYC 8958-2333-1 USNO-A2.0:USNO-A2 0225-10347916 HIP: HIP-54101

Examples

data(clcep)
f1=0.060033386
foldlc(clcep,f1)

ZTF g-band Cataclysmic Variable/Nova

Description

Time series of a cataclysmic variable/nova object observed in the g-band of the ZTF survey and processed by the ALeRCE broker.ZTF Object code: ZTF18aayzpbr

Usage

cvnovag

Format

A data frame with 67 observations on the following 3 variables:

t

heliocentric Julian Day - 2400000

m

magnitude

merr

measurement error standard deviations.

References

Förster F, Cabrera-Vives G, Castillo-Navarrete E, Estévez PA, Sánchez-Sáez P, Arredondo J, Bauer FE, Carrasco-Davis R, Catelan M, Elorrieta F, Eyheramendy S, Huijse P, Pignata G, Reyes E, Reyes I, Rodríguez-Mancini D, Ruz-Mieres D, Valenzuela C, Álvarez-Maldonado I, Astorga N, Borissova J, Clocchiatti A, Cicco DD, Donoso-Oliva C, Hernández-García L, Graham MJ, Jordán A, Kurtev R, Mahabal A, Maureira JC, Muñoz-Arancibia A, Molina-Ferreiro R, Moya A, Palma W, Pérez-Carrasco M, Protopapas P, Romero M, Sabatini-Gacitua L, Sánchez A, Martín JS, Sepúlveda-Cobo C, Vera E, Vergara JR (2021). “The Automatic Learning for the Rapid Classification of Events (ALeRCE) Alert Broker.” The Astronomical Journal, 161(5), 242. doi:10.3847/1538-3881/abe9bc.

Examples

data(cvnovag)
plot(cvnovag$t,cvnovag$m,type="l",ylab="",xlab="",col="green")

ZTF r-band Cataclysmic Variable/Nova

Description

Time series of a cataclysmic variable/nova object observed in the r-band of the ZTF survey and processed by the ALeRCE broker.ZTF Object code: ZTF18aayzpbr

Usage

cvnovar

Format

A data frame with 65 observations on the following 3 variables:

t

heliocentric Julian Day - 2400000

m

magnitude

merr

measurement error standard deviations.

References

Förster F, Cabrera-Vives G, Castillo-Navarrete E, Estévez PA, Sánchez-Sáez P, Arredondo J, Bauer FE, Carrasco-Davis R, Catelan M, Elorrieta F, Eyheramendy S, Huijse P, Pignata G, Reyes E, Reyes I, Rodríguez-Mancini D, Ruz-Mieres D, Valenzuela C, Álvarez-Maldonado I, Astorga N, Borissova J, Clocchiatti A, Cicco DD, Donoso-Oliva C, Hernández-García L, Graham MJ, Jordán A, Kurtev R, Mahabal A, Maureira JC, Muñoz-Arancibia A, Molina-Ferreiro R, Moya A, Palma W, Pérez-Carrasco M, Protopapas P, Romero M, Sabatini-Gacitua L, Sánchez A, Martín JS, Sepúlveda-Cobo C, Vera E, Vergara JR (2021). “The Automatic Learning for the Rapid Classification of Events (ALeRCE) Alert Broker.” The Astronomical Journal, 161(5), 242. doi:10.3847/1538-3881/abe9bc.

Examples

data(cvnovar)
plot(cvnovar$t,cvnovar$m,type="l",ylab="",xlab="",col="red")

Double Mode Cepheid.

Description

Time series of a double mode cepheid variable star obtained from OGLE.

Usage

dmcep

Format

A data frame with 191 observations on the following 3 variables:

t

heliocentric Julian Day

m

magnitude

merr

measurement error of the magnitude (in mag).

Details

The dominant frequency computed by GLS for this light curve is 0.7410152. The second frequency computed by GLS for this light curve is 0.5433353. OGLE-ID:175210

Examples

data(dmcep)
f1=0.7410152
foldlc(dmcep,f1)
fit=harmonicfit(dmcep,f1)
f2=0.5433353
foldlc(cbind(dmcep$t,fit$res,dmcep$merr),f2)

Delta Scuti

Description

Time series of a Delta Scuti variable star obtained from HIPPARCOS.

Usage

dscut

Format

A data frame with 116 observations on the following 3 variables:

t

heliocentric Julian Day

m

magnitude

merr

measurement error of the magnitude (in mag).

Details

The frequency computed by GLS for this light curve is 14.88558646. Catalogs and designations of this star: HD 1989: HD 199757 TYCHO-2 2000: TYC 7973-401-1 USNO-A2.0: USNO-A2 0450-39390397 HIP: HIP 103684

Examples

data(dscut)
f1=14.88558646
foldlc(dscut,f1)

Eclipsing Binaries (Beta Lyrae)

Description

Time series of a Beta Lyrae variable star obtained from OGLE.

Usage

eb

Format

A data frame with 470 observations on the following 3 variables:

t

heliocentric Julian Day

m

magnitude

merr

measurement error of the magnitude (in mag).

Details

The frequency computed by GLS for this light curve is 1.510571586. Catalogs and designations of this star:OGLE051951.22-694002.7

Examples

data(eb)
f1=1.510571586
foldlc(eb,f1)

Plotting folded light curves

Description

This function plots a time series folded on its period.

Usage

foldlc(file, f1, plot = TRUE)

Arguments

file

Matrix with the light curve observations. The first column must have the irregular times, the second column must have the brightness magnitudes and the third column must have the measurement errors.

f1

Frequency (1/Period) of the light curve.

plot

logical; if TRUE, the function returns the plot of folded time series.

Value

A matrix whose first column has the folded (phased) observational times.

Examples

data(clcep)
f1=0.060033386
foldlc(clcep,f1)

Forecast from iAR package model's

Description

Forecast with any of the models available in the iAR package

Usage

Forecast_iARModels(
  phi,
  y,
  st,
  tAhead,
  model = "iAR",
  mu = NULL,
  phiI = NULL,
  nu = NULL,
  level = 95
)

Arguments

phi

Autocorrelation coefficient estimated by the method specified.

y

Array with the time series observations.

st

Array with the observational times.

tAhead

The time ahead for which the forecast is required.

model

model to be used for the forecast. The default is to use the iAR model. Other models available are "iAR-T", "iAR-Gamma", "CiAR" and "BiAR".

mu

Level parameter of the IAR-Gamma process. A positive value.

phiI

Imaginary parameter of CIAR model or Cross-correlation parameter of BIAR model.

nu

degrees of freedom parameter of iAR-T model.

level

significance level for the confidence interval. The default value is 95.

Value

A dataframe with the following columns:

  • tAhead The time ahead used for the forecast.

  • forecast Point forecast in the time ahead required.

  • stderror Standard error of the forecast.

  • lowerCI Lower limit of the confidence interval.

  • upperCI Upper limit of the confidence interval.

References

Eyheramendy S, Elorrieta F, Palma W (2018). “An irregular discrete time series model to identify residuals with autocorrelation in astronomical light curves.” Monthly Notices of the Royal Astronomical Society, 481(4), 4311-4322. ISSN 0035-8711, doi:10.1093/mnras/sty2487, https://academic.oup.com/mnras/article-pdf/481/4/4311/25906473/sty2487.pdf.

See Also

gentime, IARforecast, IARgforecast, IARforecast, BIARforecast

Examples

st <- gentime(n=200,lambda1=15,lambda2=2)
y  <- IARsample(phi=0.9,n=200,st=st)
model<-IARloglik(y=y$series,st=st)
phi=model$phi
forIAR<-IARforecast(phi=phi,y$series,st=st,tAhead=c(1.3),standardized=FALSE,zero.mean=FALSE)
forIAR
forIAR<-Forecast_iARModels(phi=phi,y=y$series,st=st,tAhead=c(1.3,2.6))
forIAR

Generating Irregularly spaced times

Description

Function to generate irregularly spaced times from a mixture of exponential distributions.

Usage

gentime(
  n,
  distribution = "expmixture",
  lambda1 = 130,
  lambda2 = 6.5,
  p1 = 0.15,
  p2 = 0.85,
  a = 0,
  b = 1
)

Arguments

n

A positive integer. Length of observation times.

distribution

Distribution of the observation times that will be generated. Default value is "expmixture" for a mixture of exponential distributions. Alternative distributions are "uniform", "exponential" and "gamma".

lambda1

Mean (1/rate) of the exponential distribution or the first exponential distribution in a mixture of exponential distributions.

lambda2

Mean (1/rate) of the second exponential distribution in a mixture of exponential distributions.

p1

Weight of the first exponential distribution in a mixture of exponential distributions.

p2

Weight of the second exponential distribution in a mixture of exponential distributions.

a

Shape parameter of a gamma distribution or lower limit of the uniform distribution.

b

Scale parameter of a gamma distribution or upper limit of the uniform distribution.

Value

Array with irregularly spaced observations times

References

Eyheramendy S, Elorrieta F, Palma W (2018). “An irregular discrete time series model to identify residuals with autocorrelation in astronomical light curves.” Monthly Notices of the Royal Astronomical Society, 481(4), 4311-4322. ISSN 0035-8711, doi:10.1093/mnras/sty2487, https://academic.oup.com/mnras/article-pdf/481/4/4311/25906473/sty2487.pdf.

See Also

IARsample

Examples

st<-gentime(n=100)
st<-gentime(n=100,distribution="uniform")
st<-gentime(n=100,distribution="gamma",a=1,b=1)
st<-gentime(n=100,distribution="exponential",lambda1=1)

Harmonic Fit to Time Series

Description

This function fit an k-harmonic function to time series data.

Usage

harmonicfit(
  file,
  f1,
  nham = 4,
  weights = NULL,
  print = FALSE,
  remove_trend = TRUE
)

Arguments

file

A matrix with two columns. The first column corresponds to the observations times, and the second column corresponds to the measures.

f1

Frequency (1/Period) of the time series

nham

Number of harmonic components in the model

weights

An array with the weights of each observation

print

logical; if true, the summary of the harmonic fitted model will be printed. The default value is false.

remove_trend

logical; if true, the linear trend of time series will be removed before the the harmonic model is fitted.

Value

A list with the following components:

  • res Residuals to the harmonic fit of the time series.

  • t Observations times.

  • R2 Adjusted R-Squared.

  • MSE Mean Squared Error.

  • coef Summary of the coefficients estimated by the harmonic model.

Examples

data(clcep)
f1=0.060033386
results=harmonicfit(file=clcep[,1:2],f1=f1)
results$R2
results$MSE
results=harmonicfit(file=clcep[,1:2],f1=f1,nham=3)
results$R2
results$MSE
results=harmonicfit(file=clcep[,1:2],f1=f1,weights=clcep[,3])
results$R2
results$MSE

iAR: Irregularly Observed Autoregressive Models

Description

Description: Data sets, functions and scripts with examples to implement autoregressive models for irregularly observed time series. The models available in this package are the irregular autoregressive model (Eyheramendy et al.(2018) <doi:10.1093/mnras/sty2487>), the complex irregular autoregressive model (Elorrieta et al.(2019) <doi:10.1051/0004-6361/201935560>) and the bivariate irregular autoregressive model (Elorrieta et al.(2021) <doi:10.1093/mnras/stab1216>)

BIAR functions

The foo functions ...

CIAR functions

The foo functions ...

IAR functions

heloo


Fitted Values of IAR model

Description

Fit an IAR model to an irregularly observed time series.

Usage

IARfit(phi, y, st, standardized = TRUE, zero.mean = TRUE)

Arguments

phi

Estimated phi parameter by the iAR model.

y

Array with the time series observations.

st

Array with the irregular observational times.

standardized

logical; if TRUE, the array y is standardized; if FALSE, y contains the raw time series

zero.mean

logical; if TRUE, the array y has zero mean; if FALSE, y has a mean different from zero.

Value

Fitted values of the iAR model

References

Eyheramendy S, Elorrieta F, Palma W (2018). “An irregular discrete time series model to identify residuals with autocorrelation in astronomical light curves.” Monthly Notices of the Royal Astronomical Society, 481(4), 4311-4322. ISSN 0035-8711, doi:10.1093/mnras/sty2487, https://academic.oup.com/mnras/article-pdf/481/4/4311/25906473/sty2487.pdf.

See Also

gentime, IARsample, IARloglik, IARkalman

Examples

set.seed(6714)
st<-gentime(n=100)
y<-IARsample(phi=0.99,st=st,n=100)
y<-y$series
phi=IARloglik(y=y,st=st)$phi
fit=IARfit(phi=phi,y=y,st=st)

Forecast from IAR model

Description

Forecast from models fitted by IARloglik

Usage

IARforecast(phi, y, st, standardized = TRUE, zero.mean = TRUE, tAhead)

Arguments

phi

Estimated phi parameter by the iAR model.

y

Array with the time series observations.

st

Array with the irregular observational times.

standardized

logical; if TRUE, the array y is standardized; if FALSE, y contains the raw time series

zero.mean

logical; if TRUE, the array y has zero mean; if FALSE, y has a mean different from zero.

tAhead

The time ahead for forecast is required.

Value

Forecasted value from the iAR model

References

Eyheramendy S, Elorrieta F, Palma W (2018). “An irregular discrete time series model to identify residuals with autocorrelation in astronomical light curves.” Monthly Notices of the Royal Astronomical Society, 481(4), 4311-4322. ISSN 0035-8711, doi:10.1093/mnras/sty2487, https://academic.oup.com/mnras/article-pdf/481/4/4311/25906473/sty2487.pdf.

See Also

gentime, IARsample, IARloglik, IARkalman, IARfit

Examples

set.seed(6714)
st<-gentime(n=100)
y<-IARsample(phi=0.99,st=st,n=100)
y<-y$series
n=length(y)
p=trunc(n*0.99)
ytr=y[1:p]
yte=y[(p+1):n]
str=st[1:p]
ste=st[(p+1):n]
tahead=ste-str[p]
phi=IARloglik(y=ytr,st=str)$phi
forIAR=IARforecast(phi=phi,y=ytr,st=str,tAhead=tahead)

Maximum Likelihood Estimation of the IAR-Gamma model

Description

Maximum Likelihood Estimation of the IAR-Gamma model.

Usage

IARgamma(y, st)

Arguments

y

Array with the time series observations

st

Array with the irregular observational times

Value

A list with the following components:

  • phi MLE of the phi parameter of the IAR-Gamma model.

  • mu MLE of the mu parameter of the IAR-Gamma model.

  • sigma MLE of the sigma parameter of the IAR-Gamma model.

  • ll Value of the negative log likelihood evaluated in phi, mu and sigma.

References

Eyheramendy S, Elorrieta F, Palma W (2018). “An irregular discrete time series model to identify residuals with autocorrelation in astronomical light curves.” Monthly Notices of the Royal Astronomical Society, 481(4), 4311-4322. ISSN 0035-8711, doi:10.1093/mnras/sty2487, https://academic.oup.com/mnras/article-pdf/481/4/4311/25906473/sty2487.pdf.

See Also

gentime, IARgsample, IARphigamma

Examples

n=300
set.seed(6714)
st<-gentime(n)
y<-IARgsample(phi=0.9,st=st,n=n,sigma2=1,mu=1)
model<-IARgamma(y$y, st=st)
phi=model$phi
muest=model$mu
sigmaest=model$sigma

Fitted Values of IAR-Gamma model

Description

Fit an IAR-Gamma model to an irregularly observed time series.

Usage

IARgfit(phi, mu, y, st)

Arguments

phi

Estimated phi parameter by the iAR-Gamma model.

mu

Estimated mu parameter by the iAR-Gamma model.

y

Array with the time series observations.

st

Array with the irregular observational times.

Value

Fitted values of the iAR-Gamma model

References

Eyheramendy S, Elorrieta F, Palma W (2018). “An irregular discrete time series model to identify residuals with autocorrelation in astronomical light curves.” Monthly Notices of the Royal Astronomical Society, 481(4), 4311-4322. ISSN 0035-8711, doi:10.1093/mnras/sty2487, https://academic.oup.com/mnras/article-pdf/481/4/4311/25906473/sty2487.pdf.

See Also

gentime, IARgsample, IARgamma

Examples

n=300
set.seed(6714)
st<-gentime(n)
y<-IARgsample(phi=0.9,st=st,n=n,sigma2=1,mu=1)
model<-IARgamma(y$y, st=st)
phi=model$phi
muest=model$mu
sigmaest=model$sigma
fit=IARgfit(phi=phi,mu=muest,y=y$y,st=st)

Forecast from IAR-Gamma model

Description

Forecast from models fitted by IARgamma

Usage

IARgforecast(phi, mu, y, st, tAhead)

Arguments

phi

Estimated phi parameter by the iAR-Gamma model.

mu

Estimated mu parameter by the iAR-Gamma model.

y

Array with the time series observations.

st

Array with the irregular observational times.

tAhead

The time ahead for forecast is required.

Value

Forecasted value from the iAR-Gamma model

References

Eyheramendy S, Elorrieta F, Palma W (2018). “An irregular discrete time series model to identify residuals with autocorrelation in astronomical light curves.” Monthly Notices of the Royal Astronomical Society, 481(4), 4311-4322. ISSN 0035-8711, doi:10.1093/mnras/sty2487, https://academic.oup.com/mnras/article-pdf/481/4/4311/25906473/sty2487.pdf.

See Also

gentime, IARgsample, IARgamma, IARgfit

Examples

n=100
set.seed(6714)
st<-gentime(n)
y<-IARgsample(phi=0.9,st=st,n=n,sigma2=1,mu=1)
y<-y$y
n=length(y)
p=trunc(n*0.99)
ytr=y[1:p]
yte=y[(p+1):n]
str=st[1:p]
ste=st[(p+1):n]
tahead=ste-str[p]
model<-IARgamma(ytr, st=str)
phi=model$phi
muest=model$mu
sigmaest=model$sigma
fit=IARgforecast(phi=phi,mu=muest,y=ytr,st=str,tAhead=tahead)

Interpolation from IAR-Gamma model

Description

Interpolation of missing values from models fitted by IARgamma

Usage

IARginterpolation(x, y, st, yini = 1)

Arguments

x

A given array with the parameters of the IAR-Gamma model. The first element of the array corresponding to the phi parameter, the second to the level parameter mu, and the last one to the scale parameter sigma.

y

Array with the time series observations.

st

Array with the irregular observational times.

yini

a single value, initial value for the estimation of the missing value of the time series.

Value

A list with the following components:

  • fitted Estimation of a missing value of the IAR-Gamma process.

  • ll Value of the negative log likelihood evaluated in the fitted missing values.

References

Eyheramendy S, Elorrieta F, Palma W (2018). “An irregular discrete time series model to identify residuals with autocorrelation in astronomical light curves.” Monthly Notices of the Royal Astronomical Society, 481(4), 4311-4322. ISSN 0035-8711, doi:10.1093/mnras/sty2487, https://academic.oup.com/mnras/article-pdf/481/4/4311/25906473/sty2487.pdf.

See Also

gentime, IARgsample, IARgamma

Examples

set.seed(6714)
n<-100
st<-gentime(n)
y<-IARgsample(phi=0.9,st=st,n=n,sigma2=1,mu=1)
model<-IARgamma(y$y, st=st)
y<-y$y
napos=10
y0=y
y[napos]=NA
xest=c(model$phi,model$mu,model$sigma)
yest=IARginterpolation(x=xest,y=y,st=st)
yest$fitted
mse=(y0[napos]-yest$fitted)^2
print(mse)
plot(st,y,type='l',xlim=c(st[napos-5],st[napos+5]))
points(st,y,pch=20)
points(st[napos],yest$fitted,col="red",pch=20)

Simulate from an IAR-Gamma Model

Description

Simulates an IAR-Gamma Time Series Model.

Usage

IARgsample(phi, st, n = 100L, sigma2 = 1L, mu = 1L)

Arguments

phi

A coefficient of IAR-Gamma model. A value between 0 and 1.

st

Array with observational times.

n

Length of the output time series. A strictly positive integer.

sigma2

Scale parameter of the IAR-Gamma process. A positive value.

mu

Level parameter of the IAR-Gamma process. A positive value.

Value

A list with the following components:

  • y Array with simulated IAR-Gamma process.

  • st Array with observation times.

References

Eyheramendy S, Elorrieta F, Palma W (2018). “An irregular discrete time series model to identify residuals with autocorrelation in astronomical light curves.” Monthly Notices of the Royal Astronomical Society, 481(4), 4311-4322. ISSN 0035-8711, doi:10.1093/mnras/sty2487, https://academic.oup.com/mnras/article-pdf/481/4/4311/25906473/sty2487.pdf.

See Also

gentime

Examples

n=100
set.seed(6714)
st<-gentime(n)
y<-IARgsample(phi=0.9,st=st,n=n,sigma2=1,mu=1)
plot(st,y$y,type='l')
hist(y$y,breaks=20)

Interpolation from IAR model

Description

Interpolation of missing values from models fitted by IARkalman

Usage

IARinterpolation(
  x,
  y,
  st,
  delta = 0,
  yini = 0,
  zero.mean = TRUE,
  standardized = TRUE
)

Arguments

x

A given phi coefficient of the IAR model.

y

Array with the time series observations.

st

Array with the irregular observational times.

delta

Array with the measurements error standard deviations.

yini

a single value, initial value for the estimation of the missing value of the time series.

zero.mean

logical; if TRUE, the array y has zero mean; if FALSE, y has a mean different from zero.

standardized

logical; if TRUE, the array y is standardized; if FALSE, y contains the raw time series.

Value

A list with the following components:

  • fitted Estimation of a missing value of the IAR process.

  • ll Value of the negative log likelihood evaluated in the fitted missing values.

References

Eyheramendy S, Elorrieta F, Palma W (2018). “An irregular discrete time series model to identify residuals with autocorrelation in astronomical light curves.” Monthly Notices of the Royal Astronomical Society, 481(4), 4311-4322. ISSN 0035-8711, doi:10.1093/mnras/sty2487, https://academic.oup.com/mnras/article-pdf/481/4/4311/25906473/sty2487.pdf.

See Also

gentime, IARsample, IARkalman

Examples

set.seed(6714)
st<-gentime(n=100)
y<-IARsample(phi=0.99,st=st,n=100)
y<-y$series
phi=IARkalman(y=y,st=st)$phi
print(phi)
napos=10
y0=y
y[napos]=NA
xest=phi
yest=IARinterpolation(xest,y=y,st=st)
yest$fitted
mse=(y0[napos]-yest$fitted)^2
print(mse)
plot(st,y,type='l',xlim=c(st[napos-5],st[napos+5]))
points(st,y,pch=20)
points(st[napos],yest$fitted,col="red",pch=20)

Maximum Likelihood Estimation of the IAR Model via Kalman Recursions

Description

Maximum Likelihood Estimation of the IAR model parameter phi. The estimation procedure uses the Kalman Filter to find the maximum of the likelihood.

Usage

IARkalman(y, st, delta = 0, zero.mean = TRUE, standardized = TRUE)

Arguments

y

Array with the time series observations.

st

Array with the irregular observational times.

delta

Array with the measurements error standard deviations.

zero.mean

logical; if TRUE, the array y has zero mean; if FALSE, y has a mean different from zero.

standardized

logical; if TRUE, the array y is standardized; if FALSE, y contains the raw time series.

Value

A list with the following components:

  • phi MLE of the phi parameter of the IAR model.

  • ll Value of the negative log likelihood evaluated in phi.

References

Eyheramendy S, Elorrieta F, Palma W (2018). “An irregular discrete time series model to identify residuals with autocorrelation in astronomical light curves.” Monthly Notices of the Royal Astronomical Society, 481(4), 4311-4322. ISSN 0035-8711, doi:10.1093/mnras/sty2487, https://academic.oup.com/mnras/article-pdf/481/4/4311/25906473/sty2487.pdf.

See Also

gentime, IARsample, arima,IARphikalman

Examples

set.seed(6714)
st<-gentime(n=100)
y<-IARsample(phi=0.99,st=st,n=100)
y<-y$series
phi=IARkalman(y=y,st=st)$phi
print(phi)

Maximum Likelihood Estimation of the IAR Model

Description

Maximum Likelihood Estimation of the IAR Model.

Usage

IARloglik(y, st, delta = 0, zero.mean = TRUE, standardized = TRUE)

Arguments

y

Array with the time series observations.

st

Array with the irregular observational times.

delta

Array with the measurements error standard deviations.

zero.mean

logical; if TRUE, the array y has zero mean; if FALSE, y has a mean different from zero.

standardized

logical; if TRUE, the array y is standardized; if FALSE, y contains the raw time series.

Value

A list with the following components:

  • phi MLE of the phi parameter of the IAR model.

  • ll Value of the negative log likelihood evaluated in phi.

References

Eyheramendy S, Elorrieta F, Palma W (2018). “An irregular discrete time series model to identify residuals with autocorrelation in astronomical light curves.” Monthly Notices of the Royal Astronomical Society, 481(4), 4311-4322. ISSN 0035-8711, doi:10.1093/mnras/sty2487, https://academic.oup.com/mnras/article-pdf/481/4/4311/25906473/sty2487.pdf.

See Also

gentime, IARsample, arima, IARphiloglik

Examples

#Generating IAR sample
set.seed(6714)
st<-gentime(n=100)
y<-IARsample(phi=0.99,st=st,n=100)
y<-y$series
#Compute Phi
phi=IARloglik(y=y,st=st)$phi
print(phi)
#Compute the standard deviation of innovations
n=length(y)
d=c(0,diff(st))
phi1=phi**d
yhat=phi1*as.vector(c(0,y[1:(n-1)]))
plot(st,y,type='l')
lines(st,yhat,col='red')
sigma=var(y)
nu=c(sigma,sigma*(1-phi1**(2))[-1])
tau<-nu/sigma
sigmahat<-mean(c((y-yhat)**2/tau))
nuhat<-sigmahat*(1-phi1**(2))
nuhat2<-sqrt(nuhat)
#Equally spaced models
require(arfima)
fit2<-arfima(y,order=c(1,0,0))
fit<-arima(y,order=c(1,0,0),include.mean=FALSE)
syarf<-tacvfARFIMA(phi=fit2$modes[[1]]$phi,dfrac=fit2$modes[[1]]$dfrac,
sigma2=fit2$modes[[1]]$sigma,maxlag=20)[1]
syar<-fit$sigma/(1-fit$coef[1]**2)
print(sigmahat)
print(syar)
print(syarf)
carf<-fit2$modes[[1]]$sigma/syarf
car<-(1-fit$coef[1]**2)
ciar<-(1-phi1**(2))
#Compute the standard deviation of innovations (regular case)
sigma=var(y)
nuhat3=sqrt(sigma*ciar)
searf<-sqrt(sigma*carf)
sear<-sqrt(sigma*car)
#Plot the standard deviation of innovations
plot(st[-1], nuhat3[-1], t="n", axes=FALSE,xlab='Time',ylab='Standard Deviation of Innovations')
axis(1)
axis(2)
segments(x0=st[-1], y0=nuhat3[-1], y1=0, col=8)
points(st, nuhat3, pch=20, col=1, bg=1)
abline(h=sd(y),col='red',lwd=2)
abline(h=sear,col='blue',lwd=2)
abline(h=searf,col='green',lwd=2)
abline(h=mean(nuhat3[-1]),col='black',lwd=2)

Test for the significance of the autocorrelation estimated by the iAR package models

Description

This function perform a test for the significance of the autocorrelation estimated by the iAR package models. This test is based in to take N disordered samples of the original data.

Usage

IARPermutation(
  y,
  st,
  merr = 0,
  iter = 100,
  phi,
  model = "iAR",
  plot = TRUE,
  xlim = c(-1, 0),
  nu = 3
)

Arguments

y

Array with the time series observations.

st

Array with the irregular observational times.

merr

Array with the variance of the measurement errors.

iter

Number of disordered samples of the original data (N).

phi

autocorrelation estimated by one of the iAR package models.

model

model used to estimate the autocorrelation parameter ("iAR", "iAR-Gamma", "iAR-T", "CiAR" or "BiAR").

plot

logical; if true, the function return a density plot of the distribution of the bad fitted examples; if false, this function does not return a plot.

xlim

The x-axis limits (x1, x2) of the plot. Only works if plot='TRUE'. See plot.default for more details.

nu

degrees of freedom parameter of the iAR-T model.

Details

The null hypothesis of the test is: The autocorrelation coefficient estimated for the time series belongs to the distribution of the coefficients estimated on the disordered data, which are assumed to be uncorrelated. Therefore, if the hypothesis is accepted, it can be concluded that the observations of the time series are uncorrelated.The statistic of the test is log(phi) which was contrasted with a normal distribution with parameters corresponding to the log of the mean and the variance of the phi computed for the N samples of the disordered data. This test differs for IARTest in that to perform this test it is not necessary to know the period of the time series.

Value

A list with the following components:

  • phi MLE of the autocorrelation parameter of the model.

  • bad MLEs of the autocorrelation parameters of the models that has been fitted to the disordered samples.

  • norm Mean and variance of the normal distribution of the disordered data.

  • z0 Statistic of the test (log(abs(phi))).

  • pvalue P-value computed for the test.

References

Eyheramendy S, Elorrieta F, Palma W (2018). “An irregular discrete time series model to identify residuals with autocorrelation in astronomical light curves.” Monthly Notices of the Royal Astronomical Society, 481(4), 4311-4322. ISSN 0035-8711, doi:10.1093/mnras/sty2487, https://academic.oup.com/mnras/article-pdf/481/4/4311/25906473/sty2487.pdf.

See Also

Planets,IARloglik, IARTest, CIARkalman

Examples

data(Planets)
t<-Planets[,1]
res<-Planets[,2]
y=res/sqrt(var(res))
res3=IARloglik(y,t,standardized=TRUE)[1]
res3$phi
set.seed(6713)
require(ggplot2)
test<-IARPermutation(y=y,st=t,phi=res3$phi,model="iAR",plot=TRUE,xlim=c(-9.6,-9.45))

Minus Log Likelihood IAR-Gamma Model

Description

This function return the negative log likelihood of the IAR-Gamma given specific values of phi, mu and sigma.

Usage

IARphigamma(yest, x_input, y, st)

Arguments

yest

The estimate of a missing value in the time series. This function recognizes a missing value with a NA. If the time series does not have a missing value, this value does not affect the computation of the likelihood.

x_input

An array with the parameters of the IAR-Gamma model. The first element of the array corresponding to the phi parameter, the second to the level parameter mu, and the last one to the scale parameter sigma.

y

Array with the time series observations.

st

Array with the irregular observational times.

Value

Value of the negative log likelihood evaluated in phi, mu and sigma.

References

Eyheramendy S, Elorrieta F, Palma W (2018). “An irregular discrete time series model to identify residuals with autocorrelation in astronomical light curves.” Monthly Notices of the Royal Astronomical Society, 481(4), 4311-4322. ISSN 0035-8711, doi:10.1093/mnras/sty2487, https://academic.oup.com/mnras/article-pdf/481/4/4311/25906473/sty2487.pdf.

See Also

gentime, IARgsample

Examples

n=100
set.seed(6714)
st<-gentime(n)
y<-IARgsample(phi=0.9,st=st,n=n,sigma2=1,mu=1)
IARphigamma(x_input=c(0.9,1,1),y=y$y,st=st,yest=0)

Minus Log Likelihood of the IAR Model estimated via Kalman Recursions

Description

This function return the negative log likelihood of the IAR process given a specific value of phi.

Usage

IARphikalman(yest, x, y, yerr, st, zeroMean = TRUE, standardized = TRUE)

Arguments

yest

The estimate of a missing value in the time series. This function recognizes a missing value with a NA. If the time series does not have a missing value, this value does not affect the computation of the likelihood.

x

A given phi coefficient of the IAR model.

y

Array with the time series observations.

yerr

Array with the measurements error standard deviations.

st

Array with the irregular observational times.

zeroMean

logical; if TRUE, the array y has zero mean; if FALSE, y has a mean different from zero.

standardized

logical; if TRUE, the array y is standardized; if FALSE, y contains the raw time series.

Value

Value of the negative log likelihood evaluated in phi.

References

Eyheramendy S, Elorrieta F, Palma W (2018). “An irregular discrete time series model to identify residuals with autocorrelation in astronomical light curves.” Monthly Notices of the Royal Astronomical Society, 481(4), 4311-4322. ISSN 0035-8711, doi:10.1093/mnras/sty2487, https://academic.oup.com/mnras/article-pdf/481/4/4311/25906473/sty2487.pdf.

See Also

gentime, IARsample

Examples

set.seed(6714)
st<-gentime(n=100)
y<-IARsample(phi=0.99,st=st,n=100)
y<-y$series
yerr=rep(0,100)
IARphikalman(x=0.8,y=y,yerr=yerr,st=st,yest=0)

Minus Log Likelihood of the IAR Model

Description

This function return the negative log likelihood of the IAR Model for a specific value of phi.

Usage

IARphiloglik(x, y, st, delta_input, zeroMean = TRUE, standardized = TRUE)

Arguments

x

A given phi coefficient of the IAR model.

y

Array with the time series observations.

st

Array with the irregular observational times.

delta_input

Array with the measurements error standard deviations.

zeroMean

logical; if TRUE, the array y has zero mean; if FALSE, y has a mean different from zero.

standardized

logical; if TRUE, the array y was standardized; if FALSE, y contains the raw data

Value

Value of the negative log likelihood evaluated in phi.

References

Eyheramendy S, Elorrieta F, Palma W (2018). “An irregular discrete time series model to identify residuals with autocorrelation in astronomical light curves.” Monthly Notices of the Royal Astronomical Society, 481(4), 4311-4322. ISSN 0035-8711, doi:10.1093/mnras/sty2487, https://academic.oup.com/mnras/article-pdf/481/4/4311/25906473/sty2487.pdf.

See Also

gentime, IARsample

Examples

set.seed(6714)
st<-gentime(n=100)
y<-IARsample(phi=0.99,st=st,n=100)
y<-y$series
IARphiloglik(x=0.8,y=y,st=st,delta_input=c(0))

Minus Log Likelihood IAR-T Model

Description

This function return the negative log likelihood of the IAR-T given specific values of phi and sigma.

Usage

IARphit(yest, x, y, st, nu = 3)

Arguments

yest

The estimate of a missing value in the time series. This function recognizes a missing value with a NA. If the time series does not have a missing value, this value does not affect the computation of the likelihood.

x

An array with the parameters of the IAR-T model. The first element of the array corresponding to the phi parameter and the second element to the scale parameter sigma

y

Array with the time series observations

st

Array with the irregular observational times

nu

degrees of freedom

Value

Value of the negative log likelihood evaluated in phi,sigma and nu.

References

Eyheramendy S, Elorrieta F, Palma W (2018). “An irregular discrete time series model to identify residuals with autocorrelation in astronomical light curves.” Monthly Notices of the Royal Astronomical Society, 481(4), 4311-4322. ISSN 0035-8711, doi:10.1093/mnras/sty2487, https://academic.oup.com/mnras/article-pdf/481/4/4311/25906473/sty2487.pdf.

See Also

gentime, IARtsample

Examples

n=300
set.seed(6714)
st<-gentime(n)
y<-IARtsample(n,0.9,st,sigma2=1,nu=3)
IARphit(x=c(0.9,1),y=y$y,st=st,yest=0)

Simulate from an IAR Model

Description

Simulates an IAR Time Series Model.

Usage

IARsample(phi, st, n = 100L)

Arguments

phi

A coefficient of IAR model. A value between 0 and 1

st

Array with observational times.

n

Length of the output time series. A strictly positive integer.

Value

A list with the following components:

  • times Array with observation times.

  • series Array with simulated IAR data.

References

Eyheramendy S, Elorrieta F, Palma W (2018). “An irregular discrete time series model to identify residuals with autocorrelation in astronomical light curves.” Monthly Notices of the Royal Astronomical Society, 481(4), 4311-4322. ISSN 0035-8711, doi:10.1093/mnras/sty2487, https://academic.oup.com/mnras/article-pdf/481/4/4311/25906473/sty2487.pdf.

See Also

gentime

Examples

set.seed(6714)
st<-gentime(n=100)
y<-IARsample(phi=0.99,st=st, n=100)
y<-y$series
plot(st,y,type='l')

Maximum Likelihood Estimation of the IAR-T model

Description

Maximum Likelihood Estimation of the IAR-T model.

Usage

IARt(y, st, nu = 3)

Arguments

y

Array with the time series observations

st

Array with the irregular observational times

nu

degrees of freedom

Value

A list with the following components:

  • phi MLE of the phi parameter of the IAR-T model.

  • sigma MLE of the sigma parameter of the IAR-T model.

  • ll Value of the negative log likelihood evaluated in phi and sigma.

References

Eyheramendy S, Elorrieta F, Palma W (2018). “An irregular discrete time series model to identify residuals with autocorrelation in astronomical light curves.” Monthly Notices of the Royal Astronomical Society, 481(4), 4311-4322. ISSN 0035-8711, doi:10.1093/mnras/sty2487, https://academic.oup.com/mnras/article-pdf/481/4/4311/25906473/sty2487.pdf.

See Also

gentime, IARtsample, IARphit

Examples

n=300
set.seed(6714)
st<-gentime(n)
y<-IARtsample(n,0.9,st,sigma2=1,nu=3)
model<-IARt(y$y, st=st)
phi=model$phi
sigmaest=model$sigma

Test for the significance of the autocorrelation estimated by the iAR package models in periodic irregularly observed time series

Description

This function perform a test for the significance of the autocorrelation estimated by the iAR package models. This test is based on the residuals of the periodical time series fitted with an harmonic model using an incorrect period.

Usage

IARTest(
  y,
  st,
  merr = 0,
  f,
  phi,
  model = "iAR",
  plot = TRUE,
  xlim = c(-1, 0),
  nu = 3
)

Arguments

y

Array with the time series observations.

st

Array with the irregular observational times.

merr

Array with the variance of the measurement errors.

f

Frequency (1/Period) of the raw time series.

phi

autocorrelation estimated by one of the iAR package models.

model

model used to estimate the autocorrelation parameter ("iAR", "iAR-Gamma", "iAR-T", "CiAR" or "BiAR").

plot

logical; if true, the function return a density plot of the distribution of the bad fitted examples; if false, this function does not return a plot.

xlim

The x-axis limits (x1, x2) of the plot. Only works if plot='TRUE'. See plot.default for more details.

nu

degrees of freedom parameter of the iAR-T model.

Details

The null hypothesis of the test is: The autocorrelation estimated in the time series belongs to the distribution of the coefficients estimated for the residuals of the data fitted using wrong periods. Therefore, if the hypothesis is rejected, it can be concluded that the residuals of the harmonic model do not remain a time dependency structure.The statistic of the test is log(phi) which was contrasted with a normal distribution with parameters corresponding to the log of the mean and the variance of the phi computed for the residuals of the bad fitted light curves.

Value

A list with the following components:

  • phi MLE of the autocorrelation parameter of the IAR/CIAR model.

  • bad A matrix with two columns. The first column contains the incorrect frequencies used to fit each harmonic model. The second column has the MLEs of the autocorrelation parameters of the IAR/CIAR model that has been fitted to the residuals of the harmonic model fitted using the frequencies of the first column.

  • norm Mean and variance of the normal distribution of the bad fitted examples.

  • z0 Statistic of the test (log(abs(phi))).

  • pvalue P-value computed for the test.

References

Eyheramendy S, Elorrieta F, Palma W (2018). “An irregular discrete time series model to identify residuals with autocorrelation in astronomical light curves.” Monthly Notices of the Royal Astronomical Society, 481(4), 4311-4322. ISSN 0035-8711, doi:10.1093/mnras/sty2487, https://academic.oup.com/mnras/article-pdf/481/4/4311/25906473/sty2487.pdf.

See Also

clcep, harmonicfit, IARloglik, CIARkalman,IARPermutation

Examples

data(clcep)
f1=0.060033386
results=harmonicfit(file=clcep,f1=f1)
y=results$res/sqrt(var(results$res))
st=results$t
res3=IARloglik(y,st,standardized=TRUE)[1]
res3$phi
require(ggplot2)
test<-IARTest(y=clcep[,2],st=clcep[,1],f=f1,phi=res3$phi,model="iAR",plot=TRUE,xlim=c(-10,0.5))
test

Interpolation from IAR-T model

Description

Interpolation of missing values from models fitted by IARt

Usage

IARtinterpolation(x, y, st, nu = 3, yini = 0)

Arguments

x

A given array with the parameters of the IAR-T model. The first element of the array corresponding to the phi parameter and the second element to the scale parameter sigma

y

Array with the time series observations.

st

Array with the irregular observational times.

nu

degrees of freedom

yini

a single value, initial value for the estimation of the missing value of the time series.

Value

A list with the following components:

  • fitted Estimation of a missing value of the IAR-T process.

  • ll Value of the negative log likelihood evaluated in the fitted missing values.

References

Eyheramendy S, Elorrieta F, Palma W (2018). “An irregular discrete time series model to identify residuals with autocorrelation in astronomical light curves.” Monthly Notices of the Royal Astronomical Society, 481(4), 4311-4322. ISSN 0035-8711, doi:10.1093/mnras/sty2487, https://academic.oup.com/mnras/article-pdf/481/4/4311/25906473/sty2487.pdf.

See Also

gentime, IARtsample, IARt

Examples

set.seed(6714)
n<-100
st<-gentime(n)
y<-IARtsample(n,0.9,st,sigma2=1,nu=3)
model<-IARt(y$y, st=st)
napos=10
y0=y$y
y=y$y
y[napos]=NA
xest=c(model$phi,model$sigma)
yest=IARtinterpolation(x=xest,y=y,st=st)
yest$fitted
mse=(y0[napos]-yest$fitted)^2
print(mse)
plot(st,y,type='l',xlim=c(st[napos-5],st[napos+5]))
points(st,y,pch=20)
points(st[napos],yest$fitted,col="red",pch=20)

Simulate from an IAR-T Model

Description

Simulates an IAR-T Time Series Model.

Usage

IARtsample(n, phi, st, sigma2 = 1, nu = 3)

Arguments

n

Length of the output time series. A strictly positive integer.

phi

A coefficient of IAR-T model. A value between 0 and 1.

st

Array with observational times.

sigma2

Scale parameter of the IAR-T process. A positive value.

nu

degrees of freedom.

Value

A list with the following components:

  • y Array with simulated IAR-t process.

  • st Array with observation times.

References

Eyheramendy S, Elorrieta F, Palma W (2018). “An irregular discrete time series model to identify residuals with autocorrelation in astronomical light curves.” Monthly Notices of the Royal Astronomical Society, 481(4), 4311-4322. ISSN 0035-8711, doi:10.1093/mnras/sty2487, https://academic.oup.com/mnras/article-pdf/481/4/4311/25906473/sty2487.pdf.

See Also

gentime

Examples

n=300
set.seed(6714)
st<-gentime(n)
y<-IARtsample(n,0.9,st,sigma2=1,nu=3)
plot(st,y$y,type='l')
hist(y$y,breaks=20)

Pairing two irregularly observed time series

Description

Pairing the observational times of two irregularly observed time series

Usage

pairingits(lc1, lc2, tol = 0.1)

Arguments

lc1

data frame with three columns corresponding to the first irregularly observed time series. The columns must be ordered as follow: First the observational times, second the measures of each time, and third the measurement errors.

lc2

data frame with three columns corresponding to the second irregularly observed time series. The columns must be ordered as follow: First the observational times, second the measures of each time, and third the measurement errors.

tol

tolerance parameter. Minimum time gap to consider that two observations have measured at different times.

Value

A list with the following components:

  • n Number of observations paired by their observational times.

  • parData Frame with the paired datasets.

References

Elorrieta F, Eyheramendy S, Palma W, Ojeda C (2021). “A novel bivariate autoregressive model for predicting and forecasting irregularly observed time series.” Monthly Notices of the Royal Astronomical Society, 505(1), 1105-1116. ISSN 0035-8711, doi:10.1093/mnras/stab1216, https://academic.oup.com/mnras/article-pdf/505/1/1105/38391762/stab1216.pdf.

See Also

cvnovag, cvnovar, BIARkalman

Examples

data(cvnovag)
data(cvnovar)
pargr=pairingits(cvnovag,cvnovar,tol=0.1)

Transit of an extrasolar planet

Description

Time series corresponding to the residuals of the parametric model fitted by Jordan et al (2013) for a transit of an extrasolar planet.

Usage

Planets

Format

A data frame with 91 observations on the following 2 variables:

t

Time from mid-transit (hours).

r

Residuals of the parametric model fitted by Jordan et al (2013).

References

Jordán A, Espinoza N, Rabus M, Eyheramendy S, Sing D~K, Désert J, Bakos G~Á, Fortney J~J, López-Morales M, Maxted P~F~L, Triaud A~H~M~J, Szentgyorgyi A (2013). “A Ground-based Optical Transmission Spectrum of WASP-6b.” The Astrophysical Journal, 778, 184. doi:10.1088/0004-637X/778/2/184, 1310.6048.

Examples

data(Planets)
plot(Planets[,1],Planets[,2],xlab='Time from mid-transit (hours)',ylab='Noise',pch=20)