bayesian: Estimating the time delay via the Bayesian method

Description Usage Arguments Details Value Author(s) References Examples

View source: R/bayesian.R

Description

bayesian produces posterior samples of all the model parameters one of which is the time delay. This function has options for two MCMC techniques, ancillarity-sufficiency interweaving strategy (ASIS) and adaptive MCMC, to improve the convergence rate of the MCMC.

Usage

1
2
3
4
5
6
7
bayesian(data.lcA, data.lcB, data.flux, theta.ini, 
   delta.ini, delta.uniform.range, delta.proposal.scale, 
   tau.proposal.scale, tau.prior.shape, tau.prior.scale,
   sigma.prior.shape, sigma.prior.scale, asis = TRUE, micro, multimodality = FALSE, 
   adaptive.frequency, adaptive.delta =  TRUE, adaptive.delta.factor,
   adaptive.tau = TRUE, adaptive.tau.factor,
   sample.size, warmingup.size)

Arguments

data.lcA

A (n by 3) matrix for image A (light curve A); the first column has n observation times, the second column contains n flux (or magnitude) values, the third column includes n measurement errors.

data.lcB

A (w by 3) matrix for image B (light curve B); the first column has w observation times, the second column contains w flux (or magnitude) values, the third column includes w measurement errors.

data.flux

"True" if data are recorded on flux scale or "FALSE" if data are on magnitude scale. If your observed time series can take on negative values, then data.flux = FALSE.

theta.ini

Initial values of theta = (mu, sigma, tau) for MCMC.

delta.ini

Initial values of the time delay for MCMC.

delta.uniform.range

The range of the Uniform prior distribution for the time delay. The feasible entire support is c(min(simple[, 1]) - max(simple[, 1]), max(simple[, 1]) - min(simple[, 1])).

delta.proposal.scale

The proposal scale of the Metropolis step for the time delay.

tau.proposal.scale

The proposal scale of the Metropolis-Hastings step for tau.

tau.prior.shape

The shape parameter of the Inverse-Gamma hyper-prior distribution for tau.

tau.prior.scale

The scale parameter of the Inverse-Gamma hyper-prior distribution for tau.

sigma.prior.shape

The shape parameter of the Inverse-Gamma hyper-prior distribution for sigma^2.

sigma.prior.scale

The scale parameter of the Inverse-Gamma hyper-prior distribution for sigma^2. If no prior information is available, we recommend using 2 * 10^(-7).

asis

(Optional) "TRUE" if we use the ancillarity-sufficiency interweaving strategy (ASIS) for c (always recommended). Default is "TRUE".

micro

It determines the order of a polynomial regression model that accounts for the difference between microlensing trends. Default is 3. When zero is assigned, the Bayesian model fits a curve-shifted model.

multimodality

(Optional) "TRUE" if we use a repelling-attracting Metropolis algorithm for sampling Delta. When it is "TRUE" it is recommended that adaptive.delta = FALSE so that the adaptive MCMC is not used in the presence of multimodality. Default is "FALSE".

adaptive.frequency

(If "adaptive.delta = TRUE" or "adaptive.tau = TRUE") The adaptive MCMC is applied for every specified frequency. If it is specified as 500, the adaptive MCMC is applied to every 500th iterstion.

adaptive.delta

(Optional) "TRUE" if we use the adaptive MCMC for the time delay. Default is "TRUE".

adaptive.delta.factor

(If "adaptive.delta = TRUE") The factor, exp(adaptive.delta.factor) or exp(-adaptive.delta.factor), multiplied to the proposal scale of the time delay for adaptive MCMC.

adaptive.tau

(Optional) "TRUE" if we use the adaptive MCMC for tau. Default is "TRUE".

adaptive.tau.factor

(If "adaptive.tau = TRUE") The factor, exp(adaptive.tau.factor) or exp(-adaptive.tau.factor), multiplied to the proposal scale of tau for adaptive MCMC.

sample.size

The number of the posterior samples of each model parameter.

warmingup.size

The number of burn-ins for MCMC.

Details

The function bayesian produces posterior samples of the model parameters one of which is the time delay.

Value

The outcome of bayesian comprises of:

delta

Posterior samples of the time delay

X

Posterior samples of the latent magnitudes

beta

Posterior samples of the polynomial regression coefficients, beta

mu

Posterior samples of the mean parameter of the O-U process, mu

sigma

Posterior samples of the short term variability of the O-U process, sigma

tau

Posterior samples of the mean reversion time of the O-U process, tau

tau.accept.rate

The acceptance rate of the MCMC for tau

delta.accept.rate

The acceptance rate of the MCMC for the time delay

Author(s)

Hyungsuk Tak

References

Hyungsuk Tak, Kaisey Mandel, David A. van Dyk, Vinay L. Kashyap, Xiao-Li Meng, and Aneta Siemiginowska (2017+). "Bayesian Estimates of Astronomical Time Delays between Gravitationally Lensed Stochastic Light Curves," tentatively accepted in Annals of Applied Statistics (ArXiv 1602.01462). Hyungsuk Tak, Xiao-Li Meng, David A. van Dyk (2017+). "A Repelling-Attracting Metropolis Algorithm for Multimodality," submitted to Journal of Computational and Graphical Statistics (ArXiv 1601.05633).

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
  # Loading datasets
  data(simple)

  # A typical quasar data set  
  head(simple)

  # Subset (data for image A) of the typical quasar data set
  lcA <- simple[, 1 : 3]

  # Another subset (data for image B) of the typical quasar data set
  # The observation times for image B are not necessarily the same as those for image A
  lcB <- simple[, c(1, 4, 5)]

  # The two subsets do not need to have the same number of observations
  # For example, here we add one more observation time for image B
  lcB <- rbind(lcB, c(290, 1.86, 0.006))

  dim(lcA)
  dim(lcB)

  ###############################################
  # Time delay estimation via Bayesian approach #
  ###############################################

  # Cubic microlensing model (m = 3)
  output = bayesian(data.lcA = lcA, data.lcB = lcB, 
                    data.flux = FALSE, theta.ini = c(0, 0.01, 200), 
                    delta.ini = 50, delta.uniform.range = c(0, 100), 
                    delta.proposal.scale = 1, 
                    tau.proposal.scale = 3, 
                    tau.prior.shape = 1, tau.prior.scale = 1,
                    sigma.prior.shape = 1, sigma.prior.scale = 2 / 10^7, 
                    asis = TRUE, micro = 3,
                    sample.size = 100, warmingup.size = 50)

  names(output)

  # hist(output$delta, 20)
  # plot(output$delta, type = "l")
  # acf(output$delta)
  # output$delta.accept
  # output$tau.accept

  # If multimodality exists, then please add two arguments:
  # multimodality = TRUE, adaptive.delta = FALSE
  # This is to prevent the Markov chain from adapting to a local mode.
  # If we know the distance between the most distant modes,
  # try the smallest value of delta.proposal.scale that still makes 
  # the Markov chain jump between those modes frequently.
  # For example, when the distance is d, then try 
  # delta.proposal.scale = d * 0.7, decreasing or increasing 0.7.

  # Graphical model checking 
  # beta.hat <- colMeans(output$beta)
  # delta.hat <- mean(output$delta)
  # time.x <- lcB[, 1] - delta.hat
  # time.covariate <- cbind(rep(1, length(time.x)), time.x, time.x^2, time.x^3)
  ##### This time.covariate is when micro = 3, a third-order polynomial regression.
  ##### With micro = 0, "time.covariate <- rep(1, length(time.x))" is used.
  # predicted <- time.covariate %*% beta.hat
  # plot(lcA[, 1], lcA[, 2], col = 2)     
  ##### Adjust the range of the x-axis by the argument of plot.
  ##### For example "xlim = c(-1, 2)" 
  # points(lcB[, 1] - delta.hat, lcB[, 2] - predicted, col = 4, pch = 1)
  # for (i in 1 : length(output$delta)) {
  #   temp.time <- c(lcA[, 1], lcB[, 1] - output$delta[i])
  #   points(sort(temp.time), output$X[i, ], 
  #          col = "gray", cex = 0.5, pch = 1, lwd = 0.01)
  # }
  # points(lcA[, 1], lcA[, 2], pch = 0, col = 2, lwd = 1)
  # points(lcB[, 1] - delta.hat, lcB[, 2] - predicted, 
  #        col = 4, pch = 1, lwd = 1)

timedelay documentation built on May 16, 2018, 5:06 p.m.