Introduction to BLiSS method

  library(bliss)

This vignette describes step by step how to use the BLiSS method. Below, you can find the following implemented features:

  • Simulate data to test the BLiSS model
  • Obtain a sample of the posterior distribution with a Gibbs Sampler
  • Plot the posterior distribution of the coefficient function and the posterior distribution of the support
  • Compute the different Bayesian estimators

One single functional covariate case

Simulate a data set

In order to simulate a proper dataset for Bliss application, some characteristics must be specified:

  • n (the number of observations),
  • p (number of instant of measure),
  • beta_types (the shape of the coefficient function), and
  • grids_lim, a 2-vector (to define the domain of curves xi(.)).

Based on these parameters, data can be simulated (curves xi(.) and real values yi) from the functional linear regression model by using the sim function, as suggested in the following chunck.

  set.seed(1)
  param <- list(                        # define the "param" to simulate data
                Q=1,                    # the number of functional covariate
                n=50,                  # n is the sample size and p is the
                p=c(20),                # number of time observations of the curves
                beta_types=c("smooth"), # define the shape of the "true" coefficient function
                grids_lim=list(c(0,1))) # Give the beginning and the end of the observation's domain of the functions.

  data <- sim(param) # Simulate the data

How to apply the Bliss method

In order to apply the Bliss method, the main function to use is fit_Bliss. This function provides the following outputs:

  • a posterior sample of the Bliss model,
  • an approximation of the posterior distribution of the coefficient function,
  • a piecewise constant estimate (stepfunction) of the coefficient function, which is computed thanks to an optimization algorithm,
  • an estimation of the support, which shoulb be useful if your purpose is to detect periods for which the functional coviarate has an (linear) impact on the dependent scalar variable,
  • the posterior densities of the posterior sample, which should be used to compute model choice criterion.

An important required argument of the previous function is param, which is a list containing:

  • iter, the number of iterations for the Gibbs algorithm,
  • burnin the number of iteration to drop at the beginning of the Gibbs Sampler,
  • K, hyperparameter K of the Bliss model,
  • grids, the grid of instant for which the curves xi(.) are measured,
  • prior_beta, an argument specifying a prior distribution of the slope coefficient β, (only the Ridge_Zellner case is considered in this vignette),
  • and phi_l, an argument specifying a prior distribution for the half-width of the intervals (only the Gamma` case is considered in this vignette),

Find below, an example of use of this function and a sketch of the structure of the returned object.

  param <- list(            # define the required values of the Bliss method.
                iter=5e2,   # The number of iteration of the main numerical algorithm of Bliss.
                burnin=2e2, # The number of burnin iteration for the Gibbs Sampler
                K=c(3))     # The number of intervals of the beta

   
  res_bliss<-fit_Bliss(data=data,param=param)
#> ===========================================================
#> - Pretreatment.
#> - Model fiting.
#>    Gibbs Sampler: Sample from the posterior distribution.
#>   Initialization.
#>   Determine the starting point.
#>   Start the Gibbs Sampler loop.
#>    Compute beta functions related to the posterior sample.
#> - Derive estimates.
#>    Simulated Annealing (to get the Bliss and smooth estimates):
#> Functional Dimension 1/1 Functional Dimension 1/1 Repeat 1/5 Functional
#> Dimension 1/1 Repeat 2/5 Functional Dimension 1/1 Repeat 3/5 Functional
#> Dimension 1/1 Repeat 4/5 Functional Dimension 1/1 Repeat 5/5
#> 
#>    Compute the approximation of the posterior distribution.
#>    Support estimation.
#> - Posterior quantities.
#>    Compute the (log) densities of the posterior sample. 
#>    Fitting posterior y values. 
#>    Usual quantities. 
#> ===========================================================
#>     nb_param       loglik          BIC          MSE 
#> 11.000000000  0.007160281 43.017932497  0.468118502 
#> ===========================================================
#> * Useful fields 
#>    $Bliss_estime, $smooth_estimate, $support_estimate, $alpha,
#>    $beta_posterior_density, $posterior_sample, $beta_sample 
#> * Useful post-treatment functions 
#>    image_Bliss(), interpretation_plot()
  
  # Structure of a Bliss object
  # str(res_bliss)

Graphical results

This section presents how to obtain main graphical results (posterior quantities) derived from the Bliss method.

Coefficient function

Considering Functional Linear Regression model (FLR), and the specific scalar-on-function case, the major model parameter to infer is the coefficient function β(.). The following chunck shows how to plot the posterior distribution of the coefficient function:

  library(ggplot2)
  image_Bliss(res_bliss$beta_posterior_density,param,q=1) 

Additionnaly to this plot, one could usually want to display a point estimate of the coefficient function (which is a function). By using the following code, you can access to: * Bliss estimate, a piecewise constant version of the coefficient function, and * the smooth estimate, the standard bayesian estimate of the coefficient function (standard means that it minimizes the posterior L2-loss).

  image_Bliss(res_bliss$beta_posterior_density,param,q=1) + 
    lines_bliss(res_bliss$data$grids[[1]],res_bliss$Bliss_estimate[[1]]) + 
    lines_bliss(res_bliss$data$grids[[1]],res_bliss$smooth_estimate[[1]],lty = "dashed")+ 
    lines_bliss(res_bliss$data$grids[[1]],data$betas[[1]],col="purple")

The solid black line is Bliss estimate, the dashed black line is the smooth estimate and the solid purple line is the true coefficien function.

Support of coefficient function

According to the scientific problematic, one could aim to infer the coefficient function, but it is possible to alternatively focus only on the support of the coefficient function. In this case, the sign and the magniture of the coefficient function could be considered as nuisance parameters. Therefore, the Bliss method provides a specific estimation procedure for the support of the coefficient function (which relies on the posterior distribution of the coefficient function). It consists in deriving the posterior probabilities α(t|D), for each t in the domain 𝒯 of the functional data, which correspond to the probabilities (conditionnaly to the observed data) that the support of the coefficient function covers the time t.

To plot the posterior probabilities, you have to use the following code :

  plot(res_bliss$alpha[[1]],type="o",xlab="time",ylab="posterior probabilities")

From these posterior probabilities, the support estimate is derived by thresholding the probabilities. Without prior information guiding the estimation procedure, the default threshold is 0.5. The estimate support is then defined as the collection of time t for which the posterior probability α(t|D) > 0.5.

  plot(res_bliss$alpha[[1]],type="o",xlab="time",ylab="posterior probabilities")
  abline(h=0.5,col=2,lty=2)
  
  for(i in 1:nrow(res_bliss$support_estimate[[1]])){
  segments(res_bliss$support_estimate[[1]]$begin[i],0.05,
           res_bliss$support_estimate[[1]]$end[i],0.05,col="red"
           )
  points(res_bliss$support_estimate[[1]]$begin[i],0.05,col="red",pch="|",lwd=2)
  points(res_bliss$support_estimate[[1]]$end[i],0.05,col="red",pch="|",lwd=2)
  }

A resume of the support estimate is provided with:

res_bliss$support_estimate[[1]]
#>   begin end
#> 2     2   6
#> 4     8  11
#> 6    17  19

Multiple functional covariates

To avoid unnecesseray computational time, this section is not executed. You could figure out that the functions, objects and procedures are mostly similar to the previous one (single functional covariate case). The main differences are that:

  • number of function covariates have to be specified in the param object, and
  • posterior quantities are available in elements of output lists.

Simulate a data set

   param <- list(Q=2,
                 n=50,
                 p=c(40,10),
                 beta_shapes=c("simple","smooth"),
                 grids_lim=list(c(0,1),c(0,2)))

  data <- sim(param)

How to apply the Bliss method

  param <- list(       # define the required values of the Bliss method.
     iter=1e3,         # The number of iteration of the main numerical algorithm of Bliss.
     burnin=2e2,       # The number of burnin iteration for the Gibbs Sampler
     K=c(3,3))         # The number of intervals of the beta

  res_Bliss_mult <- fit_Bliss(data=data,param=param)

Graphical results

  image_Bliss(res_Bliss_mult$beta_posterior_density,param,q=1) + 
    lines_bliss(res_Bliss_mult$data$grids[[1]],res_Bliss_mult$Bliss_estimate[[1]]) + 
    lines_bliss(res_Bliss_mult$data$grids[[1]],res_Bliss_mult$smooth_estimate[[1]],lty = "dashed")+ 
    lines_bliss(res_Bliss_mult$data$grids[[1]],data$betas[[1]],col="purple")

  image_Bliss(res_Bliss_mult$beta_posterior_density,param,q=2) + 
    lines_bliss(res_Bliss_mult$data$grids[[2]],res_Bliss_mult$Bliss_estimate[[2]]) + 
    lines_bliss(res_Bliss_mult$data$grids[[2]],res_Bliss_mult$smooth_estimate[[2]],lty = "dashed")+ 
    lines_bliss(res_Bliss_mult$data$grids[[2]],data$betas[[2]],col="purple")

Session informations

#> R version 4.4.2 (2024-10-31)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.1 LTS
#> 
#> Matrix products: default
#> BLAS:   /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3 
#> LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so;  LAPACK version 3.12.0
#> 
#> locale:
#>  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
#>  [3] LC_TIME=en_US.UTF-8        LC_COLLATE=C              
#>  [5] LC_MONETARY=en_US.UTF-8    LC_MESSAGES=en_US.UTF-8   
#>  [7] LC_PAPER=en_US.UTF-8       LC_NAME=C                 
#>  [9] LC_ADDRESS=C               LC_TELEPHONE=C            
#> [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C       
#> 
#> time zone: Etc/UTC
#> tzcode source: system (glibc)
#> 
#> attached base packages:
#> [1] stats     graphics  grDevices utils     datasets  methods   base     
#> 
#> other attached packages:
#> [1] ggplot2_3.5.1 bliss_1.1.1  
#> 
#> loaded via a namespace (and not attached):
#>  [1] gtable_0.3.6           jsonlite_1.8.9         compiler_4.4.2        
#>  [4] highr_0.11             Rcpp_1.0.13-1          jquerylib_0.1.4       
#>  [7] scales_1.3.0           yaml_2.3.10            fastmap_1.2.0         
#> [10] R6_2.5.1               labeling_0.4.3         knitr_1.48            
#> [13] MASS_7.3-61            tibble_3.2.1           maketools_1.3.1       
#> [16] munsell_0.5.1          bslib_0.8.0            pillar_1.9.0          
#> [19] rlang_1.1.4            utf8_1.2.4             cachem_1.1.0          
#> [22] xfun_0.49              sass_0.4.9             sys_3.4.3             
#> [25] cli_3.6.3              withr_3.0.2            magrittr_2.0.3        
#> [28] digest_0.6.37          grid_4.4.2             lifecycle_1.0.4       
#> [31] RcppArmadillo_14.0.2-1 vctrs_0.6.5            evaluate_1.0.1        
#> [34] glue_1.8.0             farver_2.1.2           buildtools_1.0.0      
#> [37] fansi_1.0.6            colorspace_2.1-1       rmarkdown_2.29        
#> [40] tools_4.4.2            pkgconfig_2.0.3        htmltools_0.5.8.1