| Title: | Bayesian Kernelized Tensor Regression |
|---|---|
| Description: | Facilitates scalable spatiotemporally varying coefficient modelling with Bayesian kernelized tensor regression. The important features of this package are: (a) Enabling local temporal and spatial modeling of the relationship between the response variable and covariates. (b) Implementing the model described by Lei et al. (2023) <doi:10.48550/arXiv.2109.00046>. (c) Using a Bayesian Markov Chain Monte Carlo (MCMC) algorithm to sample from the posterior distribution of the model parameters. (d) Employing a tensor decomposition to reduce the number of estimated parameters. (e) Accelerating tensor operations and enabling graphics processing unit (GPU) acceleration with the 'torch' package. |
| Authors: | Julien Lanthier [aut, cre, cph] (ORCID: <https://orcid.org/0009-0008-8728-4996>), Mengying Lei [aut] (ORCID: <https://orcid.org/0000-0001-7343-3323>), Aurélie Labbe [aut] (ORCID: <https://orcid.org/0000-0002-4207-8143>), Lijun Sun [aut] (ORCID: <https://orcid.org/0000-0001-9488-0712>) |
| Maintainer: | Julien Lanthier <[email protected]> |
| License: | MIT + file LICENSE |
| Version: | 0.2.0.9000 |
| Built: | 2026-05-26 08:11:47 UTC |
| Source: | https://github.com/julien-hec/bktr |
Operator overloading for kernel multiplication
## S3 method for class 'Kernel' k1 * k2## S3 method for class 'Kernel' k1 * k2
k1 |
Kernel: The left kernel to use for composition |
k2 |
Kernel: The right kernel to use for composition |
A new KernelMulComposed object.
# Create a new locally periodic kernel k_loc_per <- KernelSE$new() * KernelPeriodic$new() # Set the kernel's positions positions_df <- data.frame(x=c(-4, 0, 3), y=c(-2, 0, 2)) k_loc_per$set_positions(positions_df) # Generate the kernel's covariance matrix k_loc_per$kernel_gen()# Create a new locally periodic kernel k_loc_per <- KernelSE$new() * KernelPeriodic$new() # Set the kernel's positions positions_df <- data.frame(x=c(-4, 0, 3), y=c(-2, 0, 2)) k_loc_per$set_positions(positions_df) # Generate the kernel's covariance matrix k_loc_per$kernel_gen()
Operator overloading for kernel addition
## S3 method for class 'Kernel' k1 + k2## S3 method for class 'Kernel' k1 + k2
k1 |
Kernel: The left kernel to use for composition |
k2 |
Kernel: The right kernel to use for composition |
A new KernelAddComposed object.
# Create a new additive kernel k_rq_plus_per <- KernelRQ$new() + KernelPeriodic$new() # Set the kernel's positions positions_df <- data.frame(x=c(-4, 0, 3), y=c(-2, 0, 2)) k_rq_plus_per$set_positions(positions_df) # Generate the kernel's covariance matrix k_rq_plus_per$kernel_gen()# Create a new additive kernel k_rq_plus_per <- KernelRQ$new() + KernelPeriodic$new() # Set the kernel's positions positions_df <- data.frame(x=c(-4, 0, 3), y=c(-2, 0, 2)) k_rq_plus_per$set_positions(positions_df) # Generate the kernel's covariance matrix k_rq_plus_per$kernel_gen()
These data represent 14 spatial features (columns) for 587 bike sharing stations (rows) located at different geographical coordinates (longitude, latitude) in Montreal. The Montreal based bike sharing company is named BIXI. The first column contains the descriptive label affected to each station and the other columns contain information about the infrastructure, points of interests, walkscore and population surrounding each station for 2019.
data("bixi_spatial_features")data("bixi_spatial_features")
A data frame with 587 observations on the following 14 variables.
locationa character vector
area_parka numeric vector
len_cycle_patha numeric vector
len_major_roada numeric vector
len_minor_roada numeric vector
num_metro_stationsa numeric vector
num_other_commerciala numeric vector
num_restaurantsa numeric vector
num_universitya numeric vector
num_popa numeric vector
num_bus_stationsa numeric vector
num_bus_routesa numeric vector
walkscorea numeric vector
capacitya numeric vector
Wang, X., Cheng, Z., Trépanier, M., & Sun, L. (2021). Modeling bike-sharing demand using a regression model with spatially varying coefficients. Journal of Transport Geography, 93, 103059.
Reference for the BIXI station informations: BIXI Montréal (2023). “Open Data.” Accessed: 2023-07-11, URL https://bixi.com/en/ open-data.
Reference for point of interests and infrastructure informations: DMTI Spatial Inc (2019). “Enhanced Point of Interest (DMTI).” URL https://www.dmtispatial.com.
Reference for Walkscore: Walk Score (2023). “Walk Score Methodology.” Accessed: 2023-07-11, URL https://www.walkscore.com/methodology.shtml.
The population information comes from the 2016 Canada census data at a dissemination block level.
Data points representing the spatial locations of 587 bike sharing stations for the Montreal based bike sharing company named BIXI. The dataframe contains a label to identify each station and its associated longitude and latitude coordinates.
data("bixi_spatial_locations")data("bixi_spatial_locations")
A data frame with 587 observations on the following 3 variables.
locationa character vector
latitudea numeric vector
longitudea numeric vector
Wang, X., Cheng, Z., Trépanier, M., & Sun, L. (2021). Modeling bike-sharing demand using a regression model with spatially varying coefficients. Journal of Transport Geography, 93, 103059.
BIXI Montréal (2023). “Open Data.” Accessed: 2023-07-11, URL https://bixi.com/en/ open-data.
These data capture the number of daily departure for 587 bike sharing stations (rows) through 197 days (columns). The data is limited to the 2019 season of a Montreal based bike sharing company named BIXI.
data("bixi_station_departures")data("bixi_station_departures")
A data frame with 587 rows and 197 columns.
Wang, X., Cheng, Z., Trépanier, M., & Sun, L. (2021). Modeling bike-sharing demand using a regression model with spatially varying coefficients. Journal of Transport Geography, 93, 103059.
BIXI Montréal (2023). “Open Data.” Accessed: 2023-07-11, URL https://bixi.com/en/ open-data.
These data represent the temporal features in Montreal applicable to a Montreal based bike sharing company named BIXI. The data include six features (columns) for 196 days (rows). The time column represent the label associated to each captured time for the 2019 season of BIXI. The other columns contain information about Montral weather and applicable holidays for each day.
data("bixi_temporal_features")data("bixi_temporal_features")
A data frame with 196 observations on the following 6 variables.
timea IDate
humiditya numeric vector
max_temp_fa numeric vector
mean_temp_ca numeric vector
total_precip_mma numeric vector
holidaya numeric vector
Lei, M., Labbe, A., & Sun, L. (2021). Scalable Spatiotemporally Varying Coefficient Modelling with Bayesian Kernelized Tensor Regression. arXiv preprint arXiv:2109.00046.
The weather data is sourced from the Environment and Climate Change Canada Historical Climate Data website.
The holiday column is specifying if a date is a holiday or not, according to the Quebec government.
These data represent 196 temporal indices (rows) related to each day of the 2019 season of Montreal based bike sharing company named BIXI. The time column represent the label associated to each day and the time_index column represent the location in time space of each day when compared to each other. Since no days are missing and they are all spaced by exactly one day, the time_index is simply a range from 0 to 195.
data("bixi_temporal_locations")data("bixi_temporal_locations")
A data frame with 196 observations on the following 2 variables.
timea IDate
time_indexa numeric vector
Lei, M., Labbe, A., & Sun, L. (2021). Scalable Spatiotemporally Varying Coefficient Modelling with Bayesian Kernelized Tensor Regression. arXiv preprint arXiv:2109.00046.
R6 class encapsulating all BIXI dataframes. It is also
possible to use a light version of the dataset by using the is_light
parameter. In this case, the dataset is reduced to its first 25 stations
and first 50 days. The light version is only used for testing and short examples.
departure_dfThe departure dataframe
spatial_features_dfThe spatial features dataframe
temporal_features_dfThe temporal features dataframe
spatial_positions_dfThe spatial positions dataframe
temporal_positions_dfThe temporal positions dataframe
data_dfThe data dataframe
is_lightWhether the light version of the dataset is used
new()
Initialize the BIXI data class
BixiData$new(is_light = FALSE)
is_lightWhether the light version of the dataset is used, defaults to FALSE.
A new BIXI data instance
clone()
The objects of this class are cloneable with this method.
BixiData$clone(deep = FALSE)
deepWhether to make a deep clone.
# Create a light BIXI data collection instance containing multiple dataframes # This only uses the first 25 stations and 50 days of the full dataset bixi_data <- BixiData$new(is_light = TRUE) # Dataframe containing the position (latitude and longitude) of M stations bixi_data$spatial_positions_df # Dataframe containing the time position of N days (O to N-1) bixi_data$temporal_positions_df # Dataframe with spatial and temporal features for each day and station (M x N rows) bixi_data$data_df# Create a light BIXI data collection instance containing multiple dataframes # This only uses the first 25 stations and 50 days of the full dataset bixi_data <- BixiData$new(is_light = TRUE) # Dataframe containing the position (latitude and longitude) of M stations bixi_data$spatial_positions_df # Dataframe containing the time position of N days (O to N-1) bixi_data$temporal_positions_df # Dataframe with spatial and temporal features for each day and station (M x N rows) bixi_data$data_df
A BKTRRegressor holds all the key elements to accomplish the MCMC sampling algorithm (Algorithm 1 of the paper).
data_dfThe dataframe containing all the covariates through time and space (including the response variable)
yThe response variable tensor
omegaThe tensor indicating which response values are not missing
covariatesThe tensor containing all the covariates
covariates_dimThe dimensions of the covariates tensor
logged_params_tensorThe tensor containing all the sampled hyperparameters
tauThe precision hyperparameter
spatial_decompThe spatial covariate decomposition
temporal_decompThe temporal covariate decomposition
covs_decompThe feature covariate decomposition
result_loggerThe result logger instance used to store the results of the MCMC sampling
has_completed_samplingBoolean showing wheter the MCMC sampling has been completed
spatial_kernelThe spatial kernel used
temporal_kernelThe temporal kernel used
spatial_positions_dfThe dataframe containing the spatial positions
temporal_positions_dfThe dataframe containing the temporal positions
spatial_params_samplerThe spatial kernel hyperparameter sampler
temporal_params_samplerThe temporal kernel hyperparameter sampler
tau_samplerThe tau hyperparameter sampler
precision_matrix_samplerThe precision matrix sampler
spatial_ll_evaluatorThe spatial likelihood evaluator
temporal_ll_evaluatorThe temporal likelihood evaluator
rank_decompThe rank of the CP decomposition
burn_in_iterThe number of burn in iterations
sampling_iterThe number of sampling iterations
max_iterThe total number of iterations
a_0The initial value for the shape in the gamma function generating tau
b_0The initial value for the rate in the gamma function generating tau
formulaThe formula used to specify the relation between the response variable and the covariates
spatial_labelsThe spatial labels
temporal_labelsThe temporal labels
feature_labelsThe feature labels
geo_coords_projectorThe geographic coordinates projector
summaryA summary of the BKTRRegressor instance
beta_covariates_summaryA dataframe containing the summary of the beta covariates
y_estimatesA dataframe containing the y estimates
imputed_y_estimatesA dataframe containing the imputed y estimates
beta_estimatesA dataframe containing the beta estimates
hyperparameters_per_iter_dfA dataframe containing the beta estimates per iteration
decomposition_tensorsList of all used decomposition tensors
new()
Create a new BKTRRegressor object.
BKTRRegressor$new( data_df, spatial_positions_df, temporal_positions_df, rank_decomp = 10, burn_in_iter = 500, sampling_iter = 500, formula = NULL, spatial_kernel = KernelMatern$new(smoothness_factor = 3), temporal_kernel = KernelSE$new(), sigma_r = 0.01, a_0 = 1e-06, b_0 = 1e-06, has_geo_coords = TRUE, geo_coords_scale = 10 )
data_dfdata.table: A dataframe containing all the covariates through time and space. It is important that the dataframe has a two indexes named 'location' and 'time' respectively. The dataframe should also contain every possible combinations of 'location' and 'time' (i.e. even missing rows should be filled present but filled with NaN). So if the dataframe has 10 locations and 5 time points, it should have 50 rows (10 x 5). If formula is None, the dataframe should contain the response variable 'Y' as the first column. Note that the covariate columns cannot contain NaN values, but the response variable can.
spatial_positions_dfdata.table: Spatial kernel input tensor used to calculate covariates' distance. Vector of length equal to the number of location points.
temporal_positions_dfdata.table: Temporal kernel input tensor used to calculate covariate distance. Vector of length equal to the number of time points.
rank_decompInteger: Rank of the CP decomposition (Paper – ). Defaults to 10.
burn_in_iterInteger: Number of iteration before sampling (Paper – ). Defaults to 500.
sampling_iterInteger: Number of sampling iterations (Paper – ). Defaults to 500.
formulaA Wilkinson R formula to specify the relation between the response variable 'Y' and the covariates. If Null, the first column of the data frame will be used as the response variable and all the other columns will be used as the covariates. Defaults to Null.
spatial_kernelKernel: Spatial kernel Used. Defaults to a KernelMatern(smoothness_factor=3).
temporal_kernelKernel: Temporal kernel used. Defaults to KernelSE().
sigma_rNumeric: Variance of the white noise process ()
defaults to 1E-2.
a_0Numeric: Initial value for the shape () in the gamma function
generating tau defaults to 1E-6.
b_0Numeric: Initial value for the rate () in the gamma function
generating tau defaults to 1E-6.
has_geo_coordsBoolean: Whether the spatial positions df use geographic coordinates (latitude, longitude). Defaults to TRUE.
geo_coords_scaleNumeric: Scale factor to convert geographic coordinates to euclidean 2D space via Mercator projection using x & y domains of [-scale/2, +scale/2]. Only used if has_geo_coords is TRUE. Defaults to 10.
A new BKTRRegressor object.
mcmc_sampling()
Launch the MCMC sampling process.
For a predefined number of iterations:
Sample spatial kernel hyperparameters
Sample temporal kernel hyperparameters
Sample the precision matrix from a wishart distribution
Sample a new spatial covariate decomposition
Sample a new feature covariate decomposition
Sample a new temporal covariate decomposition
Calculate respective errors for the iterations
Sample a new tau value
Collect all the important data for the iteration
BKTRRegressor$mcmc_sampling()
NULL Results are stored and can be accessed via summary()
predict()
Use interpolation to predict betas and response values for new data.
BKTRRegressor$predict( new_data_df, new_spatial_positions_df = NULL, new_temporal_positions_df = NULL, jitter = 1e-05 )
new_data_dfdata.table: New covariates. Must have the same columns as the covariates used to fit the model. The index should contain the combination of all old spatial coordinates with all new temporal coordinates, the combination of all new spatial coordinates with all old temporal coordinates, and the combination of all new spatial coordinates with all new temporal coordinates.
new_spatial_positions_dfdata.table or NULL: A data frame containing the new spatial positions. Defaults to NULL.
new_temporal_positions_dfdata.table or NULL: A data frame containing the new temporal positions. Defaults to NULL.
jitterNumeric or NULL: A small value to add to the diagonal of the precision matrix. Defaults to NULL.
List: A list of two dataframes. The first represents the beta forecasted for all new spatial locations or temporal points. The second represents the forecasted response for all new spatial locations or temporal points.
get_iterations_betas()
Return all sampled betas through sampling iterations for a given set of spatial, temporal and feature labels. Useful for plotting the distribution of sampled beta values.
BKTRRegressor$get_iterations_betas( spatial_label, temporal_label, feature_label )
spatial_labelString: The spatial label for which we want to get the betas
temporal_labelString: The temporal label for which we want to get the betas
feature_labelString: The feature label for which we want to get the betas
A list containing the sampled betas through iteration for the given labels
get_beta_summary_df()
Get a summary of estimated beta values. If no labels are given, then the summary is for all the betas. If labels are given, then the summary is for the given labels.
BKTRRegressor$get_beta_summary_df( spatial_labels = NULL, temporal_labels = NULL, feature_labels = NULL )
spatial_labelsvector: The spatial labels used in summary. If NULL, then all spatial labels are used. Defaults to NULL.
temporal_labelsvector: The temporal labels used in summary. If NULL, then all temporal labels are used. Defaults to NULL.
feature_labelsvector: The feature labels used in summary. If NULL, then all feature labels are used. Defaults to NULL.
A new data.table with the beta summary for the given labels.
clone()
The objects of this class are cloneable with this method.
BKTRRegressor$clone(deep = FALSE)
deepWhether to make a deep clone.
# Create a BIXI data collection instance containing multiple dataframes bixi_data <- BixiData$new(is_light = TRUE) # Use light version for example # Create a BKTRRegressor instance bktr_regressor <- BKTRRegressor$new( formula = nb_departure ~ 1 + mean_temp_c + area_park, data_df <- bixi_data$data_df, spatial_positions_df = bixi_data$spatial_positions_df, temporal_positions_df = bixi_data$temporal_positions_df, burn_in_iter = 5, sampling_iter = 10) # For example only (too few iterations) # Launch the MCMC sampling bktr_regressor$mcmc_sampling() # Get the summary of the bktr regressor summary(bktr_regressor) # Get estimated response variables for missing values bktr_regressor$imputed_y_estimates # Get the list of sampled betas for given spatial, temporal and feature labels bktr_regressor$get_iterations_betas( spatial_label = bixi_data$spatial_positions_df$location[1], temporal_label = bixi_data$temporal_positions_df$time[1], feature_label = 'mean_temp_c') # Get the summary of all betas for the 'mean_temp_c' feature bktr_regressor$get_beta_summary_df(feature_labels = 'mean_temp_c') ## PREDICTION EXAMPLE ## # Create a light version of the BIXI data collection instance bixi_data <- BixiData$new(is_light = TRUE) # Simplify variable names data_df <- bixi_data$data_df spa_pos_df <- bixi_data$spatial_positions_df temp_pos_df <- bixi_data$temporal_positions_df # Keep some data aside for prediction new_spa_pos_df <- spa_pos_df[1:2, ] new_temp_pos_df <- temp_pos_df[1:5, ] reg_spa_pos_df <- spa_pos_df[-(1:2), ] reg_temp_pos_df <- temp_pos_df[-(1:5), ] reg_data_df_mask <- data_df$location %in% reg_spa_pos_df$location & data_df$time %in% reg_temp_pos_df$time reg_data_df <- data_df[reg_data_df_mask, ] new_data_df <- data_df[!reg_data_df_mask, ] # Launch mcmc sampling on regression data bktr_regressor <- BKTRRegressor$new( formula = nb_departure ~ 1 + mean_temp_c + area_park, data_df = reg_data_df, spatial_positions_df = reg_spa_pos_df, temporal_positions_df = reg_temp_pos_df, burn_in_iter = 5, sampling_iter = 10) # For example only (too few iterations) bktr_regressor$mcmc_sampling() # Predict response values for new data bktr_regressor$predict( new_data_df = new_data_df, new_spatial_positions_df = new_spa_pos_df, new_temporal_positions_df = new_temp_pos_df)# Create a BIXI data collection instance containing multiple dataframes bixi_data <- BixiData$new(is_light = TRUE) # Use light version for example # Create a BKTRRegressor instance bktr_regressor <- BKTRRegressor$new( formula = nb_departure ~ 1 + mean_temp_c + area_park, data_df <- bixi_data$data_df, spatial_positions_df = bixi_data$spatial_positions_df, temporal_positions_df = bixi_data$temporal_positions_df, burn_in_iter = 5, sampling_iter = 10) # For example only (too few iterations) # Launch the MCMC sampling bktr_regressor$mcmc_sampling() # Get the summary of the bktr regressor summary(bktr_regressor) # Get estimated response variables for missing values bktr_regressor$imputed_y_estimates # Get the list of sampled betas for given spatial, temporal and feature labels bktr_regressor$get_iterations_betas( spatial_label = bixi_data$spatial_positions_df$location[1], temporal_label = bixi_data$temporal_positions_df$time[1], feature_label = 'mean_temp_c') # Get the summary of all betas for the 'mean_temp_c' feature bktr_regressor$get_beta_summary_df(feature_labels = 'mean_temp_c') ## PREDICTION EXAMPLE ## # Create a light version of the BIXI data collection instance bixi_data <- BixiData$new(is_light = TRUE) # Simplify variable names data_df <- bixi_data$data_df spa_pos_df <- bixi_data$spatial_positions_df temp_pos_df <- bixi_data$temporal_positions_df # Keep some data aside for prediction new_spa_pos_df <- spa_pos_df[1:2, ] new_temp_pos_df <- temp_pos_df[1:5, ] reg_spa_pos_df <- spa_pos_df[-(1:2), ] reg_temp_pos_df <- temp_pos_df[-(1:5), ] reg_data_df_mask <- data_df$location %in% reg_spa_pos_df$location & data_df$time %in% reg_temp_pos_df$time reg_data_df <- data_df[reg_data_df_mask, ] new_data_df <- data_df[!reg_data_df_mask, ] # Launch mcmc sampling on regression data bktr_regressor <- BKTRRegressor$new( formula = nb_departure ~ 1 + mean_temp_c + area_park, data_df = reg_data_df, spatial_positions_df = reg_spa_pos_df, temporal_positions_df = reg_temp_pos_df, burn_in_iter = 5, sampling_iter = 10) # For example only (too few iterations) bktr_regressor$mcmc_sampling() # Predict response values for new data bktr_regressor$predict( new_data_df = new_data_df, new_spatial_positions_df = new_spa_pos_df, new_temporal_positions_df = new_temp_pos_df)
Kernel Composition Operations Enum. Possibilities of operation between
two kernels to generate a new composed kernel. The values are: MUL and ADD.
CompositionOpsCompositionOps
An object of class list of length 2.
Abstract base class for kernels (Should not be instantiated)
kernel_varianceThe variance of the kernel
jitter_valueThe jitter value to add to the kernel matrix
distance_matrixThe distance matrix between points in a tensor format
nameThe kernel's name
parametersThe parameters of the kernel (list of KernelParameter)
covariance_matrixThe covariance matrix of the kernel in a tensor format
positions_dfThe positions of the points in a dataframe format
has_dist_matrixIdentify if the kernel has a distance matrix or not
new()
Kernel abstract base constructor
Kernel$new(kernel_variance, jitter_value)
kernel_varianceNumeric: The variance of the kernel
jitter_valueNumeric: The jitter value to add to the kernel matrix
A new Kernel object.
core_kernel_fn()
Abstract method to compute the core kernel's covariance matrix
Kernel$core_kernel_fn()
add_jitter_to_kernel()
Method to add jitter to the kernel's covariance matrix
Kernel$add_jitter_to_kernel()
kernel_gen()
Method to compute the kernel's covariance matrix
Kernel$kernel_gen()
set_positions()
Method to set the kernel's positions and compute the distance matrix
Kernel$set_positions(positions_df)
positions_dfDataframe: The positions of the points in a dataframe format
plot()
Method to plot the kernel's covariance matrix
Kernel$plot(show_figure = TRUE)
show_figureBoolean: If TRUE, the figure is shown, otherwise it is returned
If show_figure is TRUE, the figure is shown, otherwise it is returned
clone()
The objects of this class are cloneable with this method.
Kernel$clone(deep = FALSE)
deepWhether to make a deep clone.
R6 class automatically generated when adding two kernels together.
BKTR::Kernel -> BKTR::KernelComposed -> KernelAddComposed
new()
Create a new KernelAddComposed object.
KernelAddComposed$new(left_kernel, right_kernel, new_name)
left_kernelKernel: The left kernel to use for composition
right_kernelKernel: The right kernel to use for composition
new_nameString: The name of the composed kernel
A new KernelAddComposed object.
clone()
The objects of this class are cloneable with this method.
KernelAddComposed$clone(deep = FALSE)
deepWhether to make a deep clone.
# Create a new additive kernel k_rq_plus_per <- KernelAddComposed$new( left_kernel = KernelRQ$new(), right_kernel = KernelPeriodic$new(), new_name = 'SE + Periodic Kernel' ) # Set the kernel's positions positions_df <- data.frame(x=c(-4, 0, 3), y=c(-2, 0, 2)) k_rq_plus_per$set_positions(positions_df) # Generate the kernel's covariance matrix k_rq_plus_per$kernel_gen()# Create a new additive kernel k_rq_plus_per <- KernelAddComposed$new( left_kernel = KernelRQ$new(), right_kernel = KernelPeriodic$new(), new_name = 'SE + Periodic Kernel' ) # Set the kernel's positions positions_df <- data.frame(x=c(-4, 0, 3), y=c(-2, 0, 2)) k_rq_plus_per$set_positions(positions_df) # Generate the kernel's covariance matrix k_rq_plus_per$kernel_gen()
R6 class for Composed Kernels
BKTR::Kernel -> KernelComposed
nameThe kernel's name
parametersThe parameters of the kernel (list of KernelParameter)
left_kernelThe left kernel to use for composition
right_kernelThe right kernel to use for composition
composition_operationThe operation to use for composition
has_dist_matrixIdentify if the kernel has a distance matrix or not
new()
Create a new KernelComposed object.
KernelComposed$new(left_kernel, right_kernel, new_name, composition_operation)
left_kernelKernel: The left kernel to use for composition
right_kernelKernel: The right kernel to use for composition
new_nameString: The name of the composed kernel
composition_operationCompositionOps: The operation to use for composition
core_kernel_fn()
Method to compute the core kernel's covariance matrix
KernelComposed$core_kernel_fn()
The core kernel's covariance matrix
set_positions()
Method to set the kernel's positions and compute the distance matrix
KernelComposed$set_positions(positions_df)
positions_dfDataframe: The positions of the points in a dataframe format
NULL, set the kernel's positions and compute the distance matrix
clone()
The objects of this class are cloneable with this method.
KernelComposed$clone(deep = FALSE)
deepWhether to make a deep clone.
# Create a new locally periodic kernel k_loc_per <- KernelComposed$new( left_kernel = KernelSE$new(), right_kernel = KernelPeriodic$new(), new_name = 'Locally Periodic Kernel', composition_operation = CompositionOps$MUL ) # Set the kernel's positions positions_df <- data.frame(x=c(-4, 0, 3), y=c(-2, 0, 2)) k_loc_per$set_positions(positions_df) # Generate the kernel's covariance matrix k_loc_per$kernel_gen()# Create a new locally periodic kernel k_loc_per <- KernelComposed$new( left_kernel = KernelSE$new(), right_kernel = KernelPeriodic$new(), new_name = 'Locally Periodic Kernel', composition_operation = CompositionOps$MUL ) # Set the kernel's positions positions_df <- data.frame(x=c(-4, 0, 3), y=c(-2, 0, 2)) k_loc_per$set_positions(positions_df) # Generate the kernel's covariance matrix k_loc_per$kernel_gen()
R6 class for Matern Kernels
BKTR::Kernel -> KernelMatern
lengthscaleThe lengthscale parameter instance of the kernel
smoothness_factorThe smoothness factor of the kernel
has_dist_matrixIdentify if the kernel has a distance matrix or not
new()
Create a new KernelMatern object.
KernelMatern$new( smoothness_factor = 5, lengthscale = KernelParameter$new(2), kernel_variance = 1, jitter_value = NULL )
smoothness_factorNumeric: The smoothness factor of the kernel (1, 3 or 5)
lengthscaleKernelParameter: The lengthscale parameter instance of the kernel
kernel_varianceNumeric: The variance of the kernel
jitter_valueNumeric: The jitter value to add to the kernel matrix
get_smoothness_kernel_fn()
Method to the get the smoothness kernel function for a given integer smoothness factor
KernelMatern$get_smoothness_kernel_fn()
The smoothness kernel function
core_kernel_fn()
Method to compute the core kernel's covariance matrix
KernelMatern$core_kernel_fn()
The core kernel's covariance matrix
clone()
The objects of this class are cloneable with this method.
KernelMatern$clone(deep = FALSE)
deepWhether to make a deep clone.
# Create a new Matern 3/2 kernel k_matern <- KernelMatern$new(smoothness_factor = 3) # Set the kernel's positions positions_df <- data.frame(x=c(-4, 0, 3), y=c(-2, 0, 2)) k_matern$set_positions(positions_df) # Generate the kernel's covariance matrix k_matern$kernel_gen()# Create a new Matern 3/2 kernel k_matern <- KernelMatern$new(smoothness_factor = 3) # Set the kernel's positions positions_df <- data.frame(x=c(-4, 0, 3), y=c(-2, 0, 2)) k_matern$set_positions(positions_df) # Generate the kernel's covariance matrix k_matern$kernel_gen()
R6 class automatically generated when multiplying two kernels together.
BKTR::Kernel -> BKTR::KernelComposed -> KernelMulComposed
new()
Create a new KernelMulComposed object.
KernelMulComposed$new(left_kernel, right_kernel, new_name)
left_kernelKernel: The left kernel to use for composition
right_kernelKernel: The right kernel to use for composition
new_nameString: The name of the composed kernel
A new KernelMulComposed object.
clone()
The objects of this class are cloneable with this method.
KernelMulComposed$clone(deep = FALSE)
deepWhether to make a deep clone.
# Create a new locally periodic kernel k_loc_per <- KernelMulComposed$new( left_kernel = KernelSE$new(), right_kernel = KernelPeriodic$new(), new_name = 'Locally Periodic Kernel' ) # Set the kernel's positions positions_df <- data.frame(x=c(-4, 0, 3), y=c(-2, 0, 2)) k_loc_per$set_positions(positions_df) # Generate the kernel's covariance matrix k_loc_per$kernel_gen()# Create a new locally periodic kernel k_loc_per <- KernelMulComposed$new( left_kernel = KernelSE$new(), right_kernel = KernelPeriodic$new(), new_name = 'Locally Periodic Kernel' ) # Set the kernel's positions positions_df <- data.frame(x=c(-4, 0, 3), y=c(-2, 0, 2)) k_loc_per$set_positions(positions_df) # Generate the kernel's covariance matrix k_loc_per$kernel_gen()
KernelParameter contains all information and behaviour related to a kernel parameters.
valueThe hyperparameter mean's prior value or its constant value
is_fixedSays if the kernel parameter is fixed or not (if fixed, there is no sampling)
lower_boundThe hyperparameter's minimal value during sampling
upper_boundThe hyperparameter's maximal value during sampling
slice_sampling_scaleThe sampling range's amplitude
hparam_precisionPrecision of the hyperparameter
kernelThe kernel associated with the parameter (it is set at kernel instanciation)
nameThe kernel parameter's name
full_nameThe kernel parameter's full name
new()
Create a new KernelParameter object.
KernelParameter$new( value, is_fixed = FALSE, lower_bound = DEFAULT_LBOUND, upper_bound = DEFAULT_UBOUND, slice_sampling_scale = log(10), hparam_precision = 1 )
valueNumeric: The hyperparameter mean's prior value (Paper - ) or its constant value
is_fixedBoolean: Says if the kernel parameter is fixed or not (if fixed, there is no sampling)
lower_boundNumeric: Hyperparameter's minimal value during sampling (Paper - )
upper_boundNumeric: Hyperparameter's maximal value during sampling (Paper - )
slice_sampling_scaleNumeric: The sampling range's amplitude (Paper - )
hparam_precisionNumeric: The hyperparameter's precision
A new KernelParameter object.
set_kernel()
Set Kernel for a given KernelParameter object.
KernelParameter$set_kernel(kernel, param_name)
kernelKernel: The kernel to associate with the parameter
param_nameString: The parameter's name
NULL, set a new kernel for the parameter
clone()
The objects of this class are cloneable with this method.
KernelParameter$clone(deep = FALSE)
deepWhether to make a deep clone.
# A kernel parameter can be a constant value const_param <- KernelParameter$new(7, is_fixed = TRUE) # It can otherwise be sampled and have its value updated through sampling samp_param <- KernelParameter$new(1, lower_bound = 0.1, upper_bound = 10, slice_sampling_scale = 4) # A kernel parameter can be associated with any type of kernel KernelPeriodic$new(period_length = const_param, lengthscale = samp_param)# A kernel parameter can be a constant value const_param <- KernelParameter$new(7, is_fixed = TRUE) # It can otherwise be sampled and have its value updated through sampling samp_param <- KernelParameter$new(1, lower_bound = 0.1, upper_bound = 10, slice_sampling_scale = 4) # A kernel parameter can be associated with any type of kernel KernelPeriodic$new(period_length = const_param, lengthscale = samp_param)
R6 class for Periodic Kernels
BKTR::Kernel -> KernelPeriodic
lengthscaleThe lengthscale parameter instance of the kernel
period_lengthThe period length parameter instance of the kernel
has_dist_matrixIdentify if the kernel has a distance matrix or not
nameThe kernel's name
new()
Create a new KernelPeriodic object.
KernelPeriodic$new( lengthscale = KernelParameter$new(2), period_length = KernelParameter$new(2), kernel_variance = 1, jitter_value = NULL )
lengthscaleKernelParameter: The lengthscale parameter instance of the kernel
period_lengthKernelParameter: The period length parameter instance of the kernel
kernel_varianceNumeric: The variance of the kernel
jitter_valueNumeric: The jitter value to add to the kernel matrix
A new KernelPeriodic object.
core_kernel_fn()
Method to compute the core kernel's covariance matrix
KernelPeriodic$core_kernel_fn()
The core kernel's covariance matrix
clone()
The objects of this class are cloneable with this method.
KernelPeriodic$clone(deep = FALSE)
deepWhether to make a deep clone.
# Create a new Periodic kernel k_periodic <- KernelPeriodic$new() # Set the kernel's positions positions_df <- data.frame(x=c(-4, 0, 3), y=c(-2, 0, 2)) k_periodic$set_positions(positions_df) # Generate the kernel's covariance matrix k_periodic$kernel_gen()# Create a new Periodic kernel k_periodic <- KernelPeriodic$new() # Set the kernel's positions positions_df <- data.frame(x=c(-4, 0, 3), y=c(-2, 0, 2)) k_periodic$set_positions(positions_df) # Generate the kernel's covariance matrix k_periodic$kernel_gen()
R6 class for Rational Quadratic Kernels
BKTR::Kernel -> KernelRQ
lengthscaleThe lengthscale parameter instance of the kernel
alphaThe alpha parameter instance of the kernel
has_dist_matrixThe distance matrix between points in a tensor format
nameThe kernel's name
new()
Create a new KernelRQ object.
KernelRQ$new( lengthscale = KernelParameter$new(2), alpha = KernelParameter$new(2), kernel_variance = 1, jitter_value = NULL )
lengthscaleKernelParameter: The lengthscale parameter instance of the kernel
alphaKernelParameter: The alpha parameter instance of the kernel
kernel_varianceNumeric: The variance of the kernel
jitter_valueNumeric: The jitter value to add to the kernel matrix
A new KernelRQ object.
core_kernel_fn()
Method to compute the core kernel's covariance matrix
KernelRQ$core_kernel_fn()
The core kernel's covariance matrix
clone()
The objects of this class are cloneable with this method.
KernelRQ$clone(deep = FALSE)
deepWhether to make a deep clone.
# Create a new RQ kernel k_rq <- KernelRQ$new() # Set the kernel's positions positions_df <- data.frame(x=c(-4, 0, 3), y=c(-2, 0, 2)) k_rq$set_positions(positions_df) # Generate the kernel's covariance matrix k_rq$kernel_gen()# Create a new RQ kernel k_rq <- KernelRQ$new() # Set the kernel's positions positions_df <- data.frame(x=c(-4, 0, 3), y=c(-2, 0, 2)) k_rq$set_positions(positions_df) # Generate the kernel's covariance matrix k_rq$kernel_gen()
R6 class for Square Exponential Kernels
BKTR::Kernel -> KernelSE
lengthscaleThe lengthscale parameter instance of the kernel
has_dist_matrixIdentify if the kernel has a distance matrix or not
nameThe kernel's name
new()
Create a new KernelSE object.
KernelSE$new( lengthscale = KernelParameter$new(2), kernel_variance = 1, jitter_value = NULL )
lengthscaleKernelParameter: The lengthscale parameter instance of the kernel
kernel_varianceNumeric: The variance of the kernel
jitter_valueNumeric: The jitter value to add to the kernel matrix
A new KernelSE object.
core_kernel_fn()
Method to compute the core kernel's covariance matrix
KernelSE$core_kernel_fn()
The core kernel's covariance matrix
clone()
The objects of this class are cloneable with this method.
KernelSE$clone(deep = FALSE)
deepWhether to make a deep clone.
# Create a new SE kernel k_se <- KernelSE$new() # Set the kernel's positions positions_df <- data.frame(x=c(-4, 0, 3), y=c(-2, 0, 2)) k_se$set_positions(positions_df) # Generate the kernel's covariance matrix k_se$kernel_gen()# Create a new SE kernel k_se <- KernelSE$new() # Set the kernel's positions positions_df <- data.frame(x=c(-4, 0, 3), y=c(-2, 0, 2)) k_se$set_positions(positions_df) # Generate the kernel's covariance matrix k_se$kernel_gen()
R6 class for White Noise Kernels
BKTR::Kernel -> KernelWhiteNoise
has_dist_matrixIdentify if the kernel has a distance matrix or not
nameThe kernel's name
new()
KernelWhiteNoise$new(kernel_variance = 1, jitter_value = NULL)
kernel_varianceNumeric: The variance of the kernel
jitter_valueNumeric: The jitter value to add to the kernel matrix
A new KernelWhiteNoise object.
core_kernel_fn()
Method to compute the core kernel's covariance matrix
KernelWhiteNoise$core_kernel_fn()
The core kernel's covariance matrix
clone()
The objects of this class are cloneable with this method.
KernelWhiteNoise$clone(deep = FALSE)
deepWhether to make a deep clone.
# Create a new white noise kernel k_white_noise <- KernelWhiteNoise$new() # Set the kernel's positions positions_df <- data.frame(x=c(-4, 0, 3), y=c(-2, 0, 2)) k_white_noise$set_positions(positions_df) # Generate the kernel's covariance matrix k_white_noise$kernel_gen()# Create a new white noise kernel k_white_noise <- KernelWhiteNoise$new() # Set the kernel's positions positions_df <- data.frame(x=c(-4, 0, 3), y=c(-2, 0, 2)) k_white_noise$set_positions(positions_df) # Generate the kernel's covariance matrix k_white_noise$kernel_gen()
Plot the distribution of beta values for a given list of labels.
plot_beta_dists( bktr_reg, labels_list, show_figure = TRUE, fig_width = 9, fig_height = 6, fig_resolution = 200 )plot_beta_dists( bktr_reg, labels_list, show_figure = TRUE, fig_width = 9, fig_height = 6, fig_resolution = 200 )
bktr_reg |
BKTRRegressor: BKTRRegressor object. |
labels_list |
List: List of labels tuple (spatial, temporal, feature) for which to plot the beta distribution through iterations |
show_figure |
Boolean: Whether to show the figure. Defaults to True. |
fig_width |
Integer: Figure width. Defaults to 9. |
fig_height |
Integer: Figure height. Defaults to 6. |
fig_resolution |
Numeric: Figure resolution PPI. Defaults to 200. |
ggplot or NULL: ggplot object or NULL if show_figure is set to FALSE.
# Launch MCMC sampling on a light version of the BIXI dataset bixi_data <- BixiData$new(is_light = TRUE) bktr_regressor <- BKTRRegressor$new( data_df <- bixi_data$data_df, spatial_positions_df = bixi_data$spatial_positions_df, temporal_positions_df = bixi_data$temporal_positions_df, burn_in_iter = 5, sampling_iter = 10) # For example only (too few iterations) bktr_regressor$mcmc_sampling() # Plot temporal beta coefficients for the first station and the first feature spa_lab <- bixi_data$spatial_positions_df$location[3] plot_beta_dists( bktr_regressor, labels_list = list( c(spa_lab, '2019-04-15', 'area_park'), c(spa_lab, '2019-04-16', 'area_park'), c(spa_lab, '2019-04-16', 'mean_temp_c') ), )# Launch MCMC sampling on a light version of the BIXI dataset bixi_data <- BixiData$new(is_light = TRUE) bktr_regressor <- BKTRRegressor$new( data_df <- bixi_data$data_df, spatial_positions_df = bixi_data$spatial_positions_df, temporal_positions_df = bixi_data$temporal_positions_df, burn_in_iter = 5, sampling_iter = 10) # For example only (too few iterations) bktr_regressor$mcmc_sampling() # Plot temporal beta coefficients for the first station and the first feature spa_lab <- bixi_data$spatial_positions_df$location[3] plot_beta_dists( bktr_regressor, labels_list = list( c(spa_lab, '2019-04-15', 'area_park'), c(spa_lab, '2019-04-16', 'area_park'), c(spa_lab, '2019-04-16', 'mean_temp_c') ), )
Plot the distribution of beta estimates regrouped by covariates.
plot_covariates_beta_dists( bktr_reg, feature_labels = NULL, show_figure = TRUE, fig_width = 9, fig_height = 6, fig_resolution = 200 )plot_covariates_beta_dists( bktr_reg, feature_labels = NULL, show_figure = TRUE, fig_width = 9, fig_height = 6, fig_resolution = 200 )
bktr_reg |
BKTRRegressor: BKTRRegressor object. |
feature_labels |
Array or NULL: Array of feature labels for which to plot the beta estimates distribution. If NULL plot for all features. |
show_figure |
Boolean: Whether to show the figure. Defaults to True. |
fig_width |
Integer: Figure width. Defaults to 9. |
fig_height |
Integer: Figure height. Defaults to 6. |
fig_resolution |
Numeric: Figure resolution PPI. Defaults to 200. |
ggplot or NULL: ggplot object or NULL if show_figure is set to FALSE.
# Launch MCMC sampling on a light version of the BIXI dataset bixi_data <- BixiData$new(is_light = TRUE) bktr_regressor <- BKTRRegressor$new( formula = 'nb_departure ~ 1 + area_park + mean_temp_c + total_precip_mm', data_df <- bixi_data$data_df, spatial_positions_df = bixi_data$spatial_positions_df, temporal_positions_df = bixi_data$temporal_positions_df, burn_in_iter = 5, sampling_iter = 10) # For example only (too few iterations) bktr_regressor$mcmc_sampling() # Plot beta estimates distribution for all features plot_covariates_beta_dists(bktr_regressor) # Or plot for a subset of features plot_covariates_beta_dists(bktr_regressor, c('area_park', 'mean_temp_c'))# Launch MCMC sampling on a light version of the BIXI dataset bixi_data <- BixiData$new(is_light = TRUE) bktr_regressor <- BKTRRegressor$new( formula = 'nb_departure ~ 1 + area_park + mean_temp_c + total_precip_mm', data_df <- bixi_data$data_df, spatial_positions_df = bixi_data$spatial_positions_df, temporal_positions_df = bixi_data$temporal_positions_df, burn_in_iter = 5, sampling_iter = 10) # For example only (too few iterations) bktr_regressor$mcmc_sampling() # Plot beta estimates distribution for all features plot_covariates_beta_dists(bktr_regressor) # Or plot for a subset of features plot_covariates_beta_dists(bktr_regressor, c('area_park', 'mean_temp_c'))
Plot the distribution of hyperparameters through iterations
plot_hyperparams_dists( bktr_reg, hyperparameters = NULL, show_figure = TRUE, fig_width = 9, fig_height = 6, fig_resolution = 200 )plot_hyperparams_dists( bktr_reg, hyperparameters = NULL, show_figure = TRUE, fig_width = 9, fig_height = 6, fig_resolution = 200 )
bktr_reg |
BKTRRegressor: BKTRRegressor object. |
hyperparameters |
Array or NULL: Array of hyperparameters to plot. If NULL, plot all hyperparameters. Defaults to NULL. |
show_figure |
Boolean: Whether to show the figure. Defaults to True. |
fig_width |
Integer: Figure width. Defaults to 9. |
fig_height |
Integer: Figure height. Defaults to 6. |
fig_resolution |
Numeric: Figure resolution PPI. Defaults to 200. |
ggplot or NULL: ggplot object or NULL if show_figure is set to FALSE.
# Launch MCMC sampling on a light version of the BIXI dataset bixi_data <- BixiData$new(is_light = TRUE) k_matern <- KernelMatern$new() k_periodic <- KernelPeriodic$new() bktr_regressor <- BKTRRegressor$new( data_df <- bixi_data$data_df, spatial_kernel = k_matern, temporal_kernel = k_periodic, spatial_positions_df = bixi_data$spatial_positions_df, temporal_positions_df = bixi_data$temporal_positions_df, burn_in_iter = 5, sampling_iter = 10) # For example only (too few iterations) bktr_regressor$mcmc_sampling() # Plot the distribution of all hyperparameters plot_hyperparams_dists(bktr_regressor) # Plot the distribution of the spatial kernel hyperparameters spa_par_name <- paste0('Spatial - ', k_matern$parameters[[1]]$full_name) plot_hyperparams_dists(bktr_regressor, spa_par_name) # Plot the distribution of the temporal kernel hyperparameters temp_par_names <- sapply(k_periodic$parameters, function(x) x$full_name) temp_par_names <- paste0('Temporal - ', temp_par_names) plot_hyperparams_dists(bktr_regressor, temp_par_names)# Launch MCMC sampling on a light version of the BIXI dataset bixi_data <- BixiData$new(is_light = TRUE) k_matern <- KernelMatern$new() k_periodic <- KernelPeriodic$new() bktr_regressor <- BKTRRegressor$new( data_df <- bixi_data$data_df, spatial_kernel = k_matern, temporal_kernel = k_periodic, spatial_positions_df = bixi_data$spatial_positions_df, temporal_positions_df = bixi_data$temporal_positions_df, burn_in_iter = 5, sampling_iter = 10) # For example only (too few iterations) bktr_regressor$mcmc_sampling() # Plot the distribution of all hyperparameters plot_hyperparams_dists(bktr_regressor) # Plot the distribution of the spatial kernel hyperparameters spa_par_name <- paste0('Spatial - ', k_matern$parameters[[1]]$full_name) plot_hyperparams_dists(bktr_regressor, spa_par_name) # Plot the distribution of the temporal kernel hyperparameters temp_par_names <- sapply(k_periodic$parameters, function(x) x$full_name) temp_par_names <- paste0('Temporal - ', temp_par_names) plot_hyperparams_dists(bktr_regressor, temp_par_names)
Plot the evolution of hyperparameters through iterations. (Traceplot)
plot_hyperparams_traceplot( bktr_reg, hyperparameters = NULL, show_figure = TRUE, fig_width = 9, fig_height = 5.5, fig_resolution = 200 )plot_hyperparams_traceplot( bktr_reg, hyperparameters = NULL, show_figure = TRUE, fig_width = 9, fig_height = 5.5, fig_resolution = 200 )
bktr_reg |
BKTRRegressor: BKTRRegressor object. |
hyperparameters |
Array or NULL: Array of hyperparameters to plot. If NULL, plot all hyperparameters. Defaults to NULL. |
show_figure |
Boolean: Whether to show the figure. Defaults to True. |
fig_width |
Integer: Figure width. Defaults to 9. |
fig_height |
Integer: Figure height. Defaults to 5.5. |
fig_resolution |
Numeric: Figure resolution PPI. Defaults to 200. |
ggplot or NULL: ggplot object or NULL if show_figure is set to FALSE.
# Launch MCMC sampling on a light version of the BIXI dataset bixi_data <- BixiData$new(is_light = TRUE) k_matern <- KernelMatern$new() k_periodic <- KernelPeriodic$new() bktr_regressor <- BKTRRegressor$new( data_df <- bixi_data$data_df, spatial_kernel = k_matern, temporal_kernel = k_periodic, spatial_positions_df = bixi_data$spatial_positions_df, temporal_positions_df = bixi_data$temporal_positions_df, burn_in_iter = 5, sampling_iter = 10) # For example only (too few iterations) bktr_regressor$mcmc_sampling() # Plot the traceplot of all hyperparameters plot_hyperparams_traceplot(bktr_regressor) # Plot the traceplot of the spatial kernel hyperparameters spa_par_name <- paste0('Spatial - ', k_matern$parameters[[1]]$full_name) plot_hyperparams_traceplot(bktr_regressor, spa_par_name) # Plot the traceplot of the temporal kernel hyperparameters temp_par_names <- sapply(k_periodic$parameters, function(x) x$full_name) temp_par_names <- paste0('Temporal - ', temp_par_names) plot_hyperparams_traceplot(bktr_regressor, temp_par_names)# Launch MCMC sampling on a light version of the BIXI dataset bixi_data <- BixiData$new(is_light = TRUE) k_matern <- KernelMatern$new() k_periodic <- KernelPeriodic$new() bktr_regressor <- BKTRRegressor$new( data_df <- bixi_data$data_df, spatial_kernel = k_matern, temporal_kernel = k_periodic, spatial_positions_df = bixi_data$spatial_positions_df, temporal_positions_df = bixi_data$temporal_positions_df, burn_in_iter = 5, sampling_iter = 10) # For example only (too few iterations) bktr_regressor$mcmc_sampling() # Plot the traceplot of all hyperparameters plot_hyperparams_traceplot(bktr_regressor) # Plot the traceplot of the spatial kernel hyperparameters spa_par_name <- paste0('Spatial - ', k_matern$parameters[[1]]$full_name) plot_hyperparams_traceplot(bktr_regressor, spa_par_name) # Plot the traceplot of the temporal kernel hyperparameters temp_par_names <- sapply(k_periodic$parameters, function(x) x$full_name) temp_par_names <- paste0('Temporal - ', temp_par_names) plot_hyperparams_traceplot(bktr_regressor, temp_par_names)
Create a plot of beta values through space for a given temporal point and a set of feature labels. We use ggmap under the hood, so you need to provide a Google or Stadia API token to plot on a map. See: https://cran.r-project.org/web/packages/ggmap/readme/README.html for more details on how to get an API token.
plot_spatial_betas( bktr_reg, plot_feature_labels, temporal_point_label, nb_cols = 1, use_dark_mode = TRUE, show_figure = TRUE, zoom = 11, google_token = NULL, stadia_token = NULL, fig_width = 8.5, fig_height = 5.5, fig_resolution = 200 )plot_spatial_betas( bktr_reg, plot_feature_labels, temporal_point_label, nb_cols = 1, use_dark_mode = TRUE, show_figure = TRUE, zoom = 11, google_token = NULL, stadia_token = NULL, fig_width = 8.5, fig_height = 5.5, fig_resolution = 200 )
bktr_reg |
BKTRRegressor: BKTRRegressor object. |
plot_feature_labels |
Array: Array of feature labels to plot. |
temporal_point_label |
String: Temporal point label to plot. |
nb_cols |
Integer: The number of columns to use in the facet grid. |
use_dark_mode |
Boolean: Whether to use a dark mode for the geographic map or not. Defaults to TRUE. |
show_figure |
Boolean: Whether to show the figure. Defaults to True. |
zoom |
Integer: Zoom level for the geographic map. Defaults to 11. |
google_token |
String or NULL: Google API token to use for the geographic map. Defaults to NULL. |
stadia_token |
String or NULL: Stadia API token to use for the geographic map. Defaults to NULL. |
fig_width |
Numeric: Figure width when figure is shown. Defaults to 8.5. |
fig_height |
Numeric: Figure height when figure is shown. Defaults to 5.5. |
fig_resolution |
Numeric: Figure resolution PPI. Defaults to 200. |
ggplot or NULL: ggplot object or NULL if show_figure is set to FALSE.
# Launch MCMC sampling on a light version of the BIXI dataset bixi_data <- BixiData$new(is_light = TRUE) bktr_regressor <- BKTRRegressor$new( data_df <- bixi_data$data_df, spatial_positions_df = bixi_data$spatial_positions_df, temporal_positions_df = bixi_data$temporal_positions_df, burn_in_iter = 5, sampling_iter = 10) # For example only (too few iterations) bktr_regressor$mcmc_sampling() # Plot spatial beta coefficients for the first time point and the two features # Use your own Google API token instead of 'GOOGLE_API_TOKEN' plot_spatial_betas( bktr_regressor, plot_feature_labels = c('mean_temp_c', 'area_park'), temporal_point_label = bixi_data$temporal_positions_df$time[1], google_token = 'GOOGLE_API_TOKEN') # We can also use light mode and plot the maps side by side # Use your own Stadia API token instead of 'STADIA_API_TOKEN' plot_spatial_betas( bktr_regressor, plot_feature_labels = c('mean_temp_c', 'area_park', 'total_precip_mm'), temporal_point_label = bixi_data$temporal_positions_df$time[10], use_dark_mode = FALSE, nb_cols = 3, stadia_token = 'STADIA_API_TOKEN')# Launch MCMC sampling on a light version of the BIXI dataset bixi_data <- BixiData$new(is_light = TRUE) bktr_regressor <- BKTRRegressor$new( data_df <- bixi_data$data_df, spatial_positions_df = bixi_data$spatial_positions_df, temporal_positions_df = bixi_data$temporal_positions_df, burn_in_iter = 5, sampling_iter = 10) # For example only (too few iterations) bktr_regressor$mcmc_sampling() # Plot spatial beta coefficients for the first time point and the two features # Use your own Google API token instead of 'GOOGLE_API_TOKEN' plot_spatial_betas( bktr_regressor, plot_feature_labels = c('mean_temp_c', 'area_park'), temporal_point_label = bixi_data$temporal_positions_df$time[1], google_token = 'GOOGLE_API_TOKEN') # We can also use light mode and plot the maps side by side # Use your own Stadia API token instead of 'STADIA_API_TOKEN' plot_spatial_betas( bktr_regressor, plot_feature_labels = c('mean_temp_c', 'area_park', 'total_precip_mm'), temporal_point_label = bixi_data$temporal_positions_df$time[10], use_dark_mode = FALSE, nb_cols = 3, stadia_token = 'STADIA_API_TOKEN')
Create a plot of the beta values through time for a given spatial point and a set of feature labels.
plot_temporal_betas( bktr_reg, plot_feature_labels, spatial_point_label, date_format = "%Y-%m-%d", show_figure = TRUE, fig_width = 8.5, fig_height = 5.5, fig_resolution = 200 )plot_temporal_betas( bktr_reg, plot_feature_labels, spatial_point_label, date_format = "%Y-%m-%d", show_figure = TRUE, fig_width = 8.5, fig_height = 5.5, fig_resolution = 200 )
bktr_reg |
BKTRRegressor: BKTRRegressor object. |
plot_feature_labels |
Array: Array of feature labels to plot. |
spatial_point_label |
String: Spatial point label to plot. |
date_format |
String: Format of the date to use in bktr dataframes for the time. Defaults to '%Y-%m-%d'. |
show_figure |
Boolean: Whether to show the figure. Defaults to True. |
fig_width |
Numeric: Figure width when figure is shown. Defaults to 8.5. |
fig_height |
Numeric: Figure height when figure is shown. Defaults to 5.5. |
fig_resolution |
Numeric: Figure resolution PPI. Defaults to 200. |
ggplot or NULL: ggplot object or NULL if show_figure is set to FALSE.
# Launch MCMC sampling on a light version of the BIXI dataset bixi_data <- BixiData$new(is_light = TRUE) bktr_regressor <- BKTRRegressor$new( data_df <- bixi_data$data_df, spatial_positions_df = bixi_data$spatial_positions_df, temporal_positions_df = bixi_data$temporal_positions_df, burn_in_iter = 5, sampling_iter = 10) # For example only (too few iterations) bktr_regressor$mcmc_sampling() # Plot temporal beta coefficients for the first station and the two features plot_temporal_betas( bktr_regressor, plot_feature_labels = c('mean_temp_c', 'area_park'), spatial_point_label = bixi_data$spatial_positions_df$location[1])# Launch MCMC sampling on a light version of the BIXI dataset bixi_data <- BixiData$new(is_light = TRUE) bktr_regressor <- BKTRRegressor$new( data_df <- bixi_data$data_df, spatial_positions_df = bixi_data$spatial_positions_df, temporal_positions_df = bixi_data$temporal_positions_df, burn_in_iter = 5, sampling_iter = 10) # For example only (too few iterations) bktr_regressor$mcmc_sampling() # Plot temporal beta coefficients for the first station and the two features plot_temporal_betas( bktr_regressor, plot_feature_labels = c('mean_temp_c', 'area_park'), spatial_point_label = bixi_data$spatial_positions_df$location[1])
Plot y estimates vs observed y values.
plot_y_estimates( bktr_reg, show_figure = TRUE, fig_width = 5, fig_height = 5, fig_resolution = 200, fig_title = "y estimates vs observed y values" )plot_y_estimates( bktr_reg, show_figure = TRUE, fig_width = 5, fig_height = 5, fig_resolution = 200, fig_title = "y estimates vs observed y values" )
bktr_reg |
BKTRRegressor: BKTRRegressor object. |
show_figure |
Boolean: Whether to show the figure. Defaults to True. |
fig_width |
Numeric: Figure width when figure is shown. Defaults to 5. |
fig_height |
Numeric: Figure height when figure is shown. Defaults to 5. |
fig_resolution |
Numeric: Figure resolution PPI when figure is shown. Defaults to 200. |
fig_title |
String or NULL: Figure title if provided. Defaults to 'y estimates vs observed y values' |
ggplot or NULL: ggplot object or NULL if show_figure is set to FALSE.
# Launch MCMC sampling on a light version of the BIXI dataset bixi_data <- BixiData$new(is_light = TRUE) bktr_regressor <- BKTRRegressor$new( data_df <- bixi_data$data_df, spatial_positions_df = bixi_data$spatial_positions_df, temporal_positions_df = bixi_data$temporal_positions_df, burn_in_iter = 5, sampling_iter = 10) # For example only (too few iterations) bktr_regressor$mcmc_sampling() # Plot Y estimates vs observed y values plot_y_estimates(bktr_regressor)# Launch MCMC sampling on a light version of the BIXI dataset bixi_data <- BixiData$new(is_light = TRUE) bktr_regressor <- BKTRRegressor$new( data_df <- bixi_data$data_df, spatial_positions_df = bixi_data$spatial_positions_df, temporal_positions_df = bixi_data$temporal_positions_df, burn_in_iter = 5, sampling_iter = 10) # For example only (too few iterations) bktr_regressor$mcmc_sampling() # Plot Y estimates vs observed y values plot_y_estimates(bktr_regressor)
Print the summary of a BKTRRegressor instance
## S3 method for class 'BKTRRegressor' print(x, ...)## S3 method for class 'BKTRRegressor' print(x, ...)
x |
A BKTRRegressor instance |
... |
Additional arguments to comply with generic function |
Function used to transform covariates coming from two dataframes one for spatial and one for temporal into a single dataframe with the right shape for the BKTR Regressor. This is useful when the temporal covariates do not vary trough space and the spatial covariates do not vary trough time (Like in the BIXI example). The function also adds a column for the target variable at the beginning of the dataframe.
reshape_covariate_dfs(spatial_df, temporal_df, y_df, y_column_name = "y")reshape_covariate_dfs(spatial_df, temporal_df, y_df, y_column_name = "y")
spatial_df |
data.table: Spatial covariates dataframe with an index named location and a shape of (n_locations, n_spatial_covariates) |
temporal_df |
data.table: Temporal covariates dataframe with an index named time and a shape of (n_times, n_temporal_covariates) |
y_df |
data.table: The dataframe containing the target variable. It should have a shape of (n_locations, n_times). The columns and index names of this dataframe should be correspond to the one of the spatial_df and temporal_df. |
y_column_name |
string: The name of the target variable column in y_df. Default to 'y'. |
data.table: The reshaped covariates dataframe with a shape of (n_locations * n_times, 1 + n_spatial_covariates + n_temporal_covariates). The first two columns are the indexes (location, time), the following column is the target variable and the other columns are the covariates.
# Let's reshape the BIXI dataframes without normalization new_data_df <- reshape_covariate_dfs( spatial_df = BKTR::bixi_spatial_features, temporal_df = BKTR::bixi_temporal_features, y_df = BKTR::bixi_station_departures, y_column_name = 'whole_nb_departure') # The resulting dataframe has the right shape to be a BKTRRegressor data_df head(new_data_df)# Let's reshape the BIXI dataframes without normalization new_data_df <- reshape_covariate_dfs( spatial_df = BKTR::bixi_spatial_features, temporal_df = BKTR::bixi_temporal_features, y_df = BKTR::bixi_station_departures, y_column_name = 'whole_nb_departure') # The resulting dataframe has the right shape to be a BKTRRegressor data_df head(new_data_df)
Simulate Spatiotemporal Data Using Kernel Covariances.
simulate_spatiotemporal_data( nb_locations, nb_time_points, nb_spatial_dimensions, spatial_scale, time_scale, spatial_covariates_means, temporal_covariates_means, spatial_kernel, temporal_kernel, noise_variance_scale )simulate_spatiotemporal_data( nb_locations, nb_time_points, nb_spatial_dimensions, spatial_scale, time_scale, spatial_covariates_means, temporal_covariates_means, spatial_kernel, temporal_kernel, noise_variance_scale )
nb_locations |
Integer: Number of spatial locations |
nb_time_points |
Integer: Number of time points |
nb_spatial_dimensions |
Integer: Number of spatial dimensions |
spatial_scale |
Numeric: Spatial scale |
time_scale |
Numeric: Time scale |
spatial_covariates_means |
Vector: Spatial covariates means |
temporal_covariates_means |
Vector: Temporal covariates means |
spatial_kernel |
Kernel: Spatial kernel |
temporal_kernel |
Kernel: Temporal kernel |
noise_variance_scale |
Numeric: Noise variance scale |
A list containing 4 dataframes: - 'data_df' contains the response variable and the covariates - 'spatial_positions_df' contains the spatial locations and their coordinates - 'temporal_positions_df' contains the time points and their coordinates - 'beta_df' contains the true beta coefficients
# Simulate data with 20 locations, 30 time points, in 2D spatial dimensions # with 3 spatial covariates and 2 temporal covariates simu_data <- simulate_spatiotemporal_data( nb_locations=20, nb_time_points=30, nb_spatial_dimensions=2, spatial_scale=10, time_scale=10, spatial_covariates_means=c(0, 2, 4), temporal_covariates_means=c(1, 3), spatial_kernel=KernelMatern$new(), temporal_kernel=KernelSE$new(), noise_variance_scale=1) # The dataframes are similar to bixi_data, we have: # - data_df head(simu_data$data_df) # - spatial_positions_df head(simu_data$spatial_positions_df) # - temporal_positions_df head(simu_data$temporal_positions_df) # We also obtain the true beta coefficients used to simulate the data head(simu_data$beta_df)# Simulate data with 20 locations, 30 time points, in 2D spatial dimensions # with 3 spatial covariates and 2 temporal covariates simu_data <- simulate_spatiotemporal_data( nb_locations=20, nb_time_points=30, nb_spatial_dimensions=2, spatial_scale=10, time_scale=10, spatial_covariates_means=c(0, 2, 4), temporal_covariates_means=c(1, 3), spatial_kernel=KernelMatern$new(), temporal_kernel=KernelSE$new(), noise_variance_scale=1) # The dataframes are similar to bixi_data, we have: # - data_df head(simu_data$data_df) # - spatial_positions_df head(simu_data$spatial_positions_df) # - temporal_positions_df head(simu_data$temporal_positions_df) # We also obtain the true beta coefficients used to simulate the data head(simu_data$beta_df)
Summarize a BKTRRegressor instance
## S3 method for class 'BKTRRegressor' summary(object, ...)## S3 method for class 'BKTRRegressor' summary(object, ...)
object |
A BKTRRegressor instance |
... |
Additional arguments to comply with generic function |
Tensor backend configuration and methods for all the tensor operations in BKTR
R6P::Singleton -> TensorOperator
fp_typeThe floating point type to use for the tensor operations
fp_deviceThe device to use for the tensor operations
new()
Initialize the tensor operator with the given floating point type and device
TensorOperator$new(fp_type = "float64", fp_device = "cpu")
fp_typeThe floating point type to use for the tensor operations (either "float64" or "float32")
fp_deviceThe device to use for the tensor operations (either "cpu" or "cuda")
A new tensor operator instance
set_params()
Set the tensor operator parameters
TensorOperator$set_params(fp_type = NULL, fp_device = NULL, seed = NULL)
fp_typeThe floating point type to use for the tensor operations (either "float64" or "float32")
fp_deviceThe device to use for the tensor operations (either "cpu" or "cuda")
seedThe seed to use for the random number generator
get_default_jitter()
Get the default jitter value for the floating point type used by the tensor operator
TensorOperator$get_default_jitter()
The default jitter value for the floating point type used by the tensor operator
tensor()
Create a tensor from a vector or matrix of data with the tensor operator dtype and device
TensorOperator$tensor(tensor_data)
tensor_dataThe vector or matrix of data to create the tensor from
A new tensor with the tensor operator dtype and device
is_tensor()
Check if a provided object is a tensor
TensorOperator$is_tensor(tensor)
tensorThe object to check
A boolean indicating if the object is a tensor
eye()
Create a tensor with a diagonal of ones and zeros with the tensor operator dtype and device for a given dimension
TensorOperator$eye(eye_dim)
eye_dimThe dimension of the tensor to create
A new tensor with a diagonal of ones and zeros with the tensor operator dtype and device
ones()
Create a tensor of ones with the tensor operator dtype and device for a given dimension
TensorOperator$ones(tsr_dim)
tsr_dimThe dimension of the tensor to create
A new tensor of ones with the tensor operator dtype and device
zeros()
Create a tensor of zeros with the tensor operator dtype and device for a given dimension
TensorOperator$zeros(tsr_dim)
tsr_dimThe dimension of the tensor to create
A new tensor of zeros with the tensor operator dtype and device
rand()
Create a tensor of random uniform values with the tensor operator dtype and device for a given dimension
TensorOperator$rand(tsr_dim)
tsr_dimThe dimension of the tensor to create
A new tensor of random values with the tensor operator dtype and device
randn()
Create a tensor of random normal values with the tensor operator dtype and device for a given dimension
TensorOperator$randn(tsr_dim)
tsr_dimThe dimension of the tensor to create
A new tensor of random normal values with the tensor operator dtype and device
randn_like()
Create a tensor of random uniform values with the same shape as a given tensor with the tensor operator dtype and device
TensorOperator$randn_like(input_tensor)
input_tensorThe tensor to use as a shape reference
A new tensor of random uniform values with the same shape as a given tensor
arange()
Create a tensor of a range of values with the tensor operator dtype and device for a given start and end
TensorOperator$arange(start, end)
startThe start of the range
endThe end of the range
A new tensor of a range of values with the tensor operator dtype and device
rand_choice()
Choose random values from a tensor for a given number of samples
TensorOperator$rand_choice( choices_tsr, nb_sample, use_replace = FALSE, weights_tsr = NULL )
choices_tsrThe tensor to choose values from
nb_sampleThe number of samples to choose
use_replaceA boolean indicating if the sampling should be done with replacement. Defaults to FALSE
weights_tsrThe weights to use for the sampling. If NULL, the sampling is uniform. Defaults to NULL
A new tensor of randomly chosen values from a tensor
kronecker_prod()
Efficiently compute the kronecker product of two matrices in tensor format
TensorOperator$kronecker_prod(a, b)
aThe first tensor
bThe second tensor
The kronecker product of the two matrices
khatri_rao_prod()
Efficiently compute the khatri rao product of two matrices in tensor format having the same number of columns
TensorOperator$khatri_rao_prod(a, b)
aThe first tensor
bThe second tensor
The khatri rao product of the two matrices
clone()
The objects of this class are cloneable with this method.
TensorOperator$clone(deep = FALSE)
deepWhether to make a deep clone.
# Set the seed, setup the tensor floating point type and device TSR$set_params(fp_type='float64', fp_device='cpu', seed=42) # Create a tensor from a vector TSR$tensor(c(1, 2, 3)) # Create a tensor from a matrix TSR$tensor(matrix(c(1, 2, 3, 4), nrow=2)) # Create a 3x3 tensor with a diagonal of ones and zeros elsewhere TSR$eye(3) # Create a tensor of ones (with 6 elements, 2 rows and 3 columns) TSR$ones(c(2, 3)) # Create a tensor of zeros (with 12 elements, 3 rows and 4 columns) TSR$zeros(c(3, 4)) # Create a tensor of random uniform values (with 6 elements) TSR$rand(c(2, 3)) # Create a tensor of random normal values (with 6 elements) TSR$randn(c(2, 3)) # Create a tensor of random normal values with the same shape as a given tensor tsr_a <- TSR$randn(c(2, 3)) TSR$randn_like(tsr_a) # Create a tensor of a range of values (1, 2, 3, 4) TSR$arange(1, 4) # Choose two random values from a given tensor without replacement tsr_b <- TSR$rand(6) TSR$rand_choice(tsr_b, 2) # Use the tensor operator to compute the kronecker product of two 2x2 matrices tsr_c <- TSR$tensor(matrix(c(1, 2, 3, 4), nrow=2)) tsr_d <- TSR$tensor(matrix(c(5, 6, 7, 8), nrow=2)) TSR$kronecker_prod(tsr_c, tsr_d) # Returns a 4x4 tensor # Use the tensor operator to compute the khatri rao product of two 2x2 matrices TSR$khatri_rao_prod(tsr_c, tsr_d) # Returns a 4x2 tensor # Check if a given object is a tensor TSR$is_tensor(tsr_d) # Returns TRUE TSR$is_tensor(TSR$eye(2)) # Returns TRUE TSR$is_tensor(1) # Returns FALSE# Set the seed, setup the tensor floating point type and device TSR$set_params(fp_type='float64', fp_device='cpu', seed=42) # Create a tensor from a vector TSR$tensor(c(1, 2, 3)) # Create a tensor from a matrix TSR$tensor(matrix(c(1, 2, 3, 4), nrow=2)) # Create a 3x3 tensor with a diagonal of ones and zeros elsewhere TSR$eye(3) # Create a tensor of ones (with 6 elements, 2 rows and 3 columns) TSR$ones(c(2, 3)) # Create a tensor of zeros (with 12 elements, 3 rows and 4 columns) TSR$zeros(c(3, 4)) # Create a tensor of random uniform values (with 6 elements) TSR$rand(c(2, 3)) # Create a tensor of random normal values (with 6 elements) TSR$randn(c(2, 3)) # Create a tensor of random normal values with the same shape as a given tensor tsr_a <- TSR$randn(c(2, 3)) TSR$randn_like(tsr_a) # Create a tensor of a range of values (1, 2, 3, 4) TSR$arange(1, 4) # Choose two random values from a given tensor without replacement tsr_b <- TSR$rand(6) TSR$rand_choice(tsr_b, 2) # Use the tensor operator to compute the kronecker product of two 2x2 matrices tsr_c <- TSR$tensor(matrix(c(1, 2, 3, 4), nrow=2)) tsr_d <- TSR$tensor(matrix(c(5, 6, 7, 8), nrow=2)) TSR$kronecker_prod(tsr_c, tsr_d) # Returns a 4x4 tensor # Use the tensor operator to compute the khatri rao product of two 2x2 matrices TSR$khatri_rao_prod(tsr_c, tsr_d) # Returns a 4x2 tensor # Check if a given object is a tensor TSR$is_tensor(tsr_d) # Returns TRUE TSR$is_tensor(TSR$eye(2)) # Returns TRUE TSR$is_tensor(1) # Returns FALSE
Singleton instance of the TensorOperator class that contains
all informations related the tensor API; tensor methods, used data type and used device.
TSRTSR
An object of class TensorOperator (inherits from Singleton, R6) of length 19.