Title: | Dynamic Bayesian Network Learning and Inference |
---|---|
Description: | Learning and inference over dynamic Bayesian networks of arbitrary Markovian order. Extends some of the functionality offered by the 'bnlearn' package to learn the networks from data and perform exact inference. It offers three structure learning algorithms for dynamic Bayesian networks: Trabelsi G. (2013) <doi:10.1007/978-3-642-41398-8_34>, Santos F.P. and Maciel C.D. (2014) <doi:10.1109/BRC.2014.6880957>, Quesada D., Bielza C. and LarraƱaga P. (2021) <doi:10.1007/978-3-030-86271-8_14>. It also offers the possibility to perform forecasts of arbitrary length. A tool for visualizing the structure of the net is also provided via the 'visNetwork' package. |
Authors: | David Quesada [aut, cre], Gabriel Valverde [ctb] |
Maintainer: | David Quesada <[email protected]> |
License: | GPL-3 |
Version: | 0.7.9 |
Built: | 2025-02-14 05:12:22 UTC |
Source: | https://github.com/dkesada/dbnr |
Generic parameter replacement method for "dbn.fit" S3 objects. Calls bnlearn underneath.
## S3 replacement method for class 'dbn.fit' x[[name]] <- value
## S3 replacement method for class 'dbn.fit' x[[name]] <- value
x |
the fitted network |
name |
name of the node to replace its parameters |
value |
the new parameters |
the modified network
Generic parameter replacement method for "dbn.fit" S3 objects. Calls bnlearn underneath.
## S3 replacement method for class 'dbn.fit' x$name <- value
## S3 replacement method for class 'dbn.fit' x$name <- value
x |
the fitted network |
name |
name of the node to replace its parameters |
value |
the new parameters |
the modified network
Generic method for calculating the Akaike information criterion (AIC) of a
"dbn" S3 object given some data. Calls bnlearn's AIC
underneath.
## S3 method for class 'dbn' AIC(object, ..., k)
## S3 method for class 'dbn' AIC(object, ..., k)
object |
the structure of the network |
... |
additional parameters for the network scoring |
k |
the penalty parameter |
the AIC score of the network
Generic method for calculating the Akaike information criterion (AIC) of a
"dbn.fit" S3 object given some data. Calls bnlearn's AIC
underneath.
## S3 method for class 'dbn.fit' AIC(object, ..., k)
## S3 method for class 'dbn.fit' AIC(object, ..., k)
object |
the fitted network |
... |
additional parameters for the network scoring |
k |
the penalty parameter |
the AIC score of the network
Generic method for checking the equality of two
"dbn" S3 objects. Calls bnlearn's all.equal
underneath.
## S3 method for class 'equal.dbn' all(target, current, ...)
## S3 method for class 'equal.dbn' all(target, current, ...)
target |
"dbn" object |
current |
the other "dbn" object |
... |
additional parameters |
boolean result of the comparison
Generic method for checking the equality of two
"dbn.fit" S3 objects. Calls bnlearn's all.equal
underneath.
## S3 method for class 'equal.dbn.fit' all(target, current, ...)
## S3 method for class 'equal.dbn.fit' all(target, current, ...)
target |
"dbn.fit" object |
current |
the other "dbn.fit" object |
... |
additional parameters |
boolean result of the comparison
Generic method for converting a "dbn" S3 object into a string.
Calls bnlearn's as.character
underneath.
## S3 method for class 'dbn' as.character(x, ...)
## S3 method for class 'dbn' as.character(x, ...)
x |
a "dbn" object |
... |
additional parameters |
string representing the DBN model
Generic method for calculating the Bayesian information criterion (BIC) of a
"dbn" S3 object given some data. Calls bnlearn's BIC
underneath.
## S3 method for class 'dbn' BIC(object, ...)
## S3 method for class 'dbn' BIC(object, ...)
object |
the structure of the network |
... |
additional parameters for the network scoring |
the BIC score of the network
Generic method for calculating the Bayesian information criterion (BIC) of a
"dbn.fit" S3 object given some data. Calls bnlearn's BIC
underneath.
## S3 method for class 'dbn.fit' BIC(object, ...)
## S3 method for class 'dbn.fit' BIC(object, ...)
object |
the fitted network |
... |
additional parameters for the network scoring |
the BIC score of the network
Given a "bn.fit" or a "dbn.fit" object, calculate the mu vector of the equivalent multivariate Gaussian distribution. Front end of a C++ function.
calc_mu(fit)
calc_mu(fit)
fit |
a bn.fit or dbn.fit object |
a named numeric vector of the means of each variable
dt_train <- dbnR::motor[200:2500] net <- bnlearn::mmhc(dt_train) fit <- bnlearn::bn.fit(net, dt_train, method = "mle-g") mu <- dbnR::calc_mu(fit) f_dt_train <- dbnR::fold_dt(dt_train, size = 2) net <- dbnR::learn_dbn_struc(dt_train, size = 2) fit <- dbnR::fit_dbn_params(net, f_dt_train) mu <- dbnR::calc_mu(fit)
dt_train <- dbnR::motor[200:2500] net <- bnlearn::mmhc(dt_train) fit <- bnlearn::bn.fit(net, dt_train, method = "mle-g") mu <- dbnR::calc_mu(fit) f_dt_train <- dbnR::fold_dt(dt_train, size = 2) net <- dbnR::learn_dbn_struc(dt_train, size = 2) fit <- dbnR::fit_dbn_params(net, f_dt_train) mu <- dbnR::calc_mu(fit)
Given a "bn.fit" or a "dbn.fit" object, calculate the sigma covariance matrix of the equivalent multivariate Gaussian distribution. Front end of a C++ function.
calc_sigma(fit)
calc_sigma(fit)
fit |
a bn.fit or dbn.fit object |
a named numeric covariance matrix of the nodes
dt_train <- dbnR::motor[200:2500] net <- bnlearn::mmhc(dt_train) fit <- bnlearn::bn.fit(net, dt_train, method = "mle-g") sigma <- dbnR::calc_sigma(fit) f_dt_train <- dbnR::fold_dt(dt_train, size = 2) net <- dbnR::learn_dbn_struc(dt_train, size = 2) fit <- dbnR::fit_dbn_params(net, f_dt_train) sigma <- dbnR::calc_sigma(fit)
dt_train <- dbnR::motor[200:2500] net <- bnlearn::mmhc(dt_train) fit <- bnlearn::bn.fit(net, dt_train, method = "mle-g") sigma <- dbnR::calc_sigma(fit) f_dt_train <- dbnR::fold_dt(dt_train, size = 2) net <- dbnR::learn_dbn_struc(dt_train, size = 2) fit <- dbnR::fit_dbn_params(net, f_dt_train) sigma <- dbnR::calc_sigma(fit)
Generic method for "dbn.fit" S3 objects. Calls bnlearn underneath.
## S3 method for class 'dbn.fit' coef(object, ...)
## S3 method for class 'dbn.fit' coef(object, ...)
object |
the fitted network |
... |
additional parameters |
the coefficients of the network
#' Generic method for calculating the degree of a list of nodes in a
BN or a DBN. Calls bnlearn's degree
underneath.
I have to redefine the generic and mask the original for it to work on both
bn and dbn objects without the user having to import bnlearn.
degree(object, Nodes, ...)
degree(object, Nodes, ...)
object |
a "bn", "dbn", "bn.fit" or "dbn.fit" object |
Nodes |
which nodes to check |
... |
additional parameters |
the degree of the nodes
Given an id variable that labels the different instances of a time series inside a dataset, discard the rows that have values from more than 1 id.
filter_same_cycle(f_dt, size, id_var)
filter_same_cycle(f_dt, size, id_var)
f_dt |
folded data.table |
size |
the size of the data.table |
id_var |
the variable that labels each individual instance of the time series |
the filtered data.table
dt <- dbnR::motor[201:2500] dt[, n_sec := rep(seq(46), each = 50)] # I'll create secuences of 50 instances each f_dt <- dbnR::fold_dt(dt, size = 2) f_dt[50, .SD, .SDcols = c("n_sec_t_0", "n_sec_t_1")] f_dt <- dbnR::filter_same_cycle(f_dt, size = 2, id_var = "n_sec") f_dt[50, .SD, .SDcols = c("n_sec_t_0", "n_sec_t_1")]
dt <- dbnR::motor[201:2500] dt[, n_sec := rep(seq(46), each = 50)] # I'll create secuences of 50 instances each f_dt <- dbnR::fold_dt(dt, size = 2) f_dt[50, .SD, .SDcols = c("n_sec_t_0", "n_sec_t_1")] f_dt <- dbnR::filter_same_cycle(f_dt, size = 2, id_var = "n_sec") f_dt[50, .SD, .SDcols = c("n_sec_t_0", "n_sec_t_1")]
If the dataset that is going to be folded contains several different time series instances of the same process, folding it could introduce false rows with data from different time series. Given an id variable that labels the different instances of a time series inside a dataset and a desired size, this function folds the dataset and avoids mixing data from different origins in the same instance.
filtered_fold_dt(dt, size, id_var, clear_id_var = TRUE)
filtered_fold_dt(dt, size, id_var, clear_id_var = TRUE)
dt |
data.table to be folded |
size |
the size of the data.table |
id_var |
the variable that labels each individual instance of the time series |
clear_id_var |
boolean that decides whether or not the id_var column is deleted |
the filtered data.table
dt <- dbnR::motor[201:2500] dt[, n_sec := rep(seq(46), each = 50)] # I'll create secuences of 50 instances each f_dt <- dbnR::fold_dt(dt, size = 2) dim(f_dt) f_dt <- dbnR::filtered_fold_dt(dt, size = 2, id_var = "n_sec") dim(f_dt) # The filtered folded dt has a row less for each independent secuence
dt <- dbnR::motor[201:2500] dt[, n_sec := rep(seq(46), each = 50)] # I'll create secuences of 50 instances each f_dt <- dbnR::fold_dt(dt, size = 2) dim(f_dt) f_dt <- dbnR::filtered_fold_dt(dt, size = 2, id_var = "n_sec") dim(f_dt) # The filtered folded dt has a row less for each independent secuence
Fits the parameters of the DBN via MLE. The "mu" vector of means and the "sigma" covariance matrix are set as attributes of the dbn.fit object for future exact inference.
fit_dbn_params(net, f_dt, ...)
fit_dbn_params(net, f_dt, ...)
net |
the structure of the DBN |
f_dt |
a folded data.table |
... |
additional parameters for the |
a "dbn.fit" S3 object with the fitted net
size = 3 dt_train <- dbnR::motor[200:2500] net <- learn_dbn_struc(dt_train, size) f_dt_train <- fold_dt(dt_train, size) fit <- fit_dbn_params(net, f_dt_train, method = "mle-g")
size = 3 dt_train <- dbnR::motor[200:2500] net <- learn_dbn_struc(dt_train, size) f_dt_train <- fold_dt(dt_train, size) fit <- fit_dbn_params(net, f_dt_train, method = "mle-g")
Generic method for "dbn.fit" S3 objects. Calls bnlearn underneath.
## S3 method for class 'dbn.fit' fitted(object, ...)
## S3 method for class 'dbn.fit' fitted(object, ...)
object |
the fitted network |
... |
additional parameters |
the fitted values of the network
This function will widen the dataset to put the t previous time slices in each row, so that it can be used to learn temporal arcs in the second phase of the dmmhc.
fold_dt(dt, size)
fold_dt(dt, size)
dt |
the data.table to be treated |
size |
number of time slices to unroll. Markovian 1 would be size 2 |
the extended data.table
data(motor) size <- 3 f_dt <- fold_dt(motor, size)
data(motor) size <- 3 f_dt <- fold_dt(motor, size)
Given a dbn.fit object, the size of the net and a folded dataset, performs a forecast over the initial evidence taken from the dataset.
forecast_ts( dt, fit, size = NULL, obj_vars, ini = 1, len = dim(dt)[1] - ini, rep = 1, num_p = 50, print_res = TRUE, plot_res = TRUE, mode = "exact", prov_ev = NULL )
forecast_ts( dt, fit, size = NULL, obj_vars, ini = 1, len = dim(dt)[1] - ini, rep = 1, num_p = 50, print_res = TRUE, plot_res = TRUE, mode = "exact", prov_ev = NULL )
dt |
data.table object with the TS data |
fit |
dbn.fit object |
size |
number of time slices of the net. Deprecated, will be removed in the future |
obj_vars |
variables to be predicted |
ini |
starting point in the dataset to forecast. |
len |
length of the forecast |
rep |
number of times to repeat the approximate forecasting |
num_p |
number of particles in the approximate forecasting |
print_res |
if TRUE prints the mae and sd metrics of the forecast |
plot_res |
if TRUE plots the results of the forecast |
mode |
"exact" for exact inference, "approx" for approximate |
prov_ev |
variables to be provided as evidence in each forecasting step |
a list with the original time series values and the results of the forecast
size = 3 data(motor) dt_train <- motor[200:900] dt_val <- motor[901:1000] obj <- c("pm_t_0") net <- learn_dbn_struc(dt_train, size) f_dt_train <- fold_dt(dt_train, size) f_dt_val <- fold_dt(dt_val, size) fit <- fit_dbn_params(net, f_dt_train, method = "mle-g") res <- suppressWarnings(forecast_ts(f_dt_val, fit, obj_vars = obj, len = 10, print_res = FALSE, plot_res = FALSE))
size = 3 data(motor) dt_train <- motor[200:900] dt_val <- motor[901:1000] obj <- c("pm_t_0") net <- learn_dbn_struc(dt_train, size) f_dt_train <- fold_dt(dt_train, size) f_dt_val <- fold_dt(dt_val, size) fit <- fit_dbn_params(net, f_dt_train, method = "mle-g") res <- suppressWarnings(forecast_ts(f_dt_val, fit, obj_vars = obj, len = 10, print_res = FALSE, plot_res = FALSE))
This function generates both a random DBN and a dataset that can be used to learn its structure from data. It's intended for experimental use.
generate_random_network_exp( n_vars, size, min_mu, max_mu, min_sd, max_sd, min_coef, max_coef, seed = NULL )
generate_random_network_exp( n_vars, size, min_mu, max_mu, min_sd, max_sd, min_coef, max_coef, seed = NULL )
n_vars |
number of desired variables per time-slice |
size |
desired size of the networks |
min_mu |
minimum mean allowed for the variables |
max_mu |
maximum mean allowed for the variables |
min_sd |
minimum standard deviation allowed for the variables |
max_sd |
maximum standard deviation allowed for the variables |
min_coef |
minimum coefficient allowed for the parent nodes |
max_coef |
maximum coefficient allowed for the parent nodes |
seed |
the seed of the experiment |
a list with the original network structure and the sampled dataset
Learns a gaussian dynamic Bayesian network from a dataset. It allows the creation of markovian n nets rather than only markov 1.
learn_dbn_struc(dt, size = 2, method = "dmmhc", f_dt = NULL, ...)
learn_dbn_struc(dt, size = 2, method = "dmmhc", f_dt = NULL, ...)
dt |
the data.frame or data.table to be used |
size |
number of time slices of the net. Markovian 1 would be size 2 |
method |
the structure learning method of choice to use |
f_dt |
previously folded dataset, in case some specific rows have to be removed after the folding |
... |
additional parameters for |
a "dbn" S3 object with the structure of the network
data("motor") net <- learn_dbn_struc(motor, size = 3)
data("motor") net <- learn_dbn_struc(motor, size = 3)
Generic method for calculating the log-likelihood of a
"dbn" S3 object given some data. Calls bnlearn's logLik
underneath.
## S3 method for class 'dbn' logLik(object, dt, ...)
## S3 method for class 'dbn' logLik(object, dt, ...)
object |
the structure of the network |
dt |
the dataset to calculate the score of the network |
... |
additional parameters for the network scoring |
the log-likelihood score of the network
Generic method for calculating the log-likelihood of a
"dbn.fit" S3 object given some data. Calls bnlearn's logLik
underneath.
## S3 method for class 'dbn.fit' logLik(object, dt, ...)
## S3 method for class 'dbn.fit' logLik(object, dt, ...)
object |
the fitted network |
dt |
the dataset to calculate the score of the network |
... |
additional parameters for the network scoring |
the log-likelihood score of the network
Generic method for "dbn.fit" S3 objects. Calls bnlearn underneath.
## S3 method for class 'dbn.fit' mean(x, ...)
## S3 method for class 'dbn.fit' mean(x, ...)
x |
the fitted network |
... |
additional parameters |
the averaged parameters
Data from several sensors on an electric motor that records different benchmark sessions of measurements at 2 Hz. The dataset is reduced to 3000 instances from the 60th session in order to include it in the package for testing purposes. For the complete dataset, refer to the source.
data(motor)
data(motor)
An object of class data.table
(inherits from data.frame
) with 3000 rows and 11 columns.
Kaggle, <https://www.kaggle.com/wkirgsn/electric-motor-temperature>
Given some evidence, this function performs inference over a multivariate normal
distribution. After converting a Gaussian linear network to its MVN form, this
kind of inference can be performed. It's recommended to use
predict_dt
functions instead unless you need a more flexible
inference method.
mvn_inference(mu, sigma, evidence)
mvn_inference(mu, sigma, evidence)
mu |
the mean vector |
sigma |
the covariance matrix |
evidence |
a single row data.table or a named vector with the values and names of the variables given as evidence |
a list with the posterior mean and covariance matrix
size = 3 data(motor) dt_train <- motor[200:2500] dt_val <- motor[2501:3000] obj <- c("pm_t_0") net <- learn_dbn_struc(dt_train, size) f_dt_train <- fold_dt(dt_train, size) f_dt_val <- fold_dt(dt_val, size) ev <- f_dt_val[1, .SD, .SDcols = obj] fit <- fit_dbn_params(net, f_dt_train, method = "mle-g") pred <- mvn_inference(calc_mu(fit), calc_sigma(fit), ev)
size = 3 data(motor) dt_train <- motor[200:2500] dt_val <- motor[2501:3000] obj <- c("pm_t_0") net <- learn_dbn_struc(dt_train, size) f_dt_train <- fold_dt(dt_train, size) f_dt_val <- fold_dt(dt_val, size) ev <- f_dt_val[1, .SD, .SDcols = obj] fit <- fit_dbn_params(net, f_dt_train, method = "mle-g") pred <- mvn_inference(calc_mu(fit), calc_sigma(fit), ev)
Generic method for obtaining the names of the nodes in a BN or a DBN.
Calls bnlearn's nodes
underneath.
I have to redefine the generic and mask the original for it to work on both
bn and dbn objects without the user having to import bnlearn.
nodes(object, ...)
nodes(object, ...)
object |
a "bn", "dbn", "bn.fit" or "dbn.fit" object |
... |
additional parameters |
the names of the nodes
Generic method for renaming the nodes in a BN or a DBN.
Calls bnlearn's nodes<-
underneath.
I have to redefine the generic and mask the original for it to work on both
bn and dbn objects without the user having to import bnlearn.
nodes(object) <- value
nodes(object) <- value
object |
a "bn", "dbn", "bn.fit" or "dbn.fit" object |
value |
a list with the new names |
the modified object
To plot the DBN, this method first computes a hierarchical structure for a time slice and replicates it for each slice. Then, it calculates the relative position of each node with respect to his equivalent in the first slice. The result is a net where each time slice is ordered and separated from one another, where the leftmost slice is the oldest and the rightmost represents the present time. This function is also called by the generic plot function of "dbn" and "dbn.fit" S3 objects.
plot_dynamic_network( structure, offset = 200, subset_nodes = NULL, reverse = FALSE )
plot_dynamic_network( structure, offset = 200, subset_nodes = NULL, reverse = FALSE )
structure |
the structure or fit of the network. |
offset |
the blank space between time slices |
subset_nodes |
a vector containing the names of the subset of nodes to plot |
reverse |
reverse to the classic naming convention of the nodes. The oldest time-slice will now be t_0 and the most recent one t_n. Only for visualization purposes, the network is unmodified underneath. If using subset_nodes, remember that t_0 is now the oldest time-slice. |
the visualization of the DBN
size = 3 dt_train <- dbnR::motor[200:2500] net <- learn_dbn_struc(dt_train, size) plot_dynamic_network(net)
size = 3 dt_train <- dbnR::motor[200:2500] net <- learn_dbn_struc(dt_train, size) plot_dynamic_network(net)
This function calculates the levels of each node and then plots them in a hierarchical layout in visNetwork. Can be used in place of the generic plot function offered by bnlearn for "bn" and "bn.fit" S3 objects.
plot_static_network(structure)
plot_static_network(structure)
structure |
the structure or fit of the network. |
dt_train <- dbnR::motor[200:2500] net <- bnlearn::mmhc(dt_train) plot_static_network(net) fit <- bnlearn::bn.fit(net, dt_train, method = "mle-g") plot_static_network(fit) # Works for both the structure and the fitted net
dt_train <- dbnR::motor[200:2500] net <- bnlearn::mmhc(dt_train) plot_static_network(net) fit <- bnlearn::bn.fit(net, dt_train, method = "mle-g") plot_static_network(fit) # Works for both the structure and the fitted net
Generic method for plotting the "dbn" S3 objects. Calls
plot_dynamic_network
underneath.
## S3 method for class 'dbn' plot(x, ...)
## S3 method for class 'dbn' plot(x, ...)
x |
the structure of the network. |
... |
additional parameters for the visualization of a DBN |
Generic method for plotting the "dbn.fit" S3 objects. Calls
plot_dynamic_network
underneath.
## S3 method for class 'dbn.fit' plot(x, ...)
## S3 method for class 'dbn.fit' plot(x, ...)
x |
the structure of the network. |
... |
additional parameters for the visualization of a DBN |
Performs inference over a Gaussian BN. It's thought to be used in a map for
a data.table, to use as evidence each separate row. If not specifically
needed, it's recommended to use the function predict_dt
instead.
This function is deprecated and will be removed in a future version.
predict_bn(fit, evidence)
predict_bn(fit, evidence)
fit |
the fitted bn |
evidence |
values of the variables used as evidence for the net |
a data.table with the predictions
size = 3 data(motor) dt_train <- motor[200:2500] dt_val <- motor[2501:3000] net <- learn_dbn_struc(dt_train, size) f_dt_train <- fold_dt(dt_train, size) f_dt_val <- fold_dt(dt_val, size) fit <- fit_dbn_params(net, f_dt_train, method = "mle-g") res <- f_dt_val[, predict_bn(fit, .SD), .SDcols = c("pm_t_0", "coolant_t_0"), by = 1:nrow(f_dt_val)]
size = 3 data(motor) dt_train <- motor[200:2500] dt_val <- motor[2501:3000] net <- learn_dbn_struc(dt_train, size) f_dt_train <- fold_dt(dt_train, size) f_dt_val <- fold_dt(dt_val, size) fit <- fit_dbn_params(net, f_dt_train, method = "mle-g") res <- f_dt_val[, predict_bn(fit, .SD), .SDcols = c("pm_t_0", "coolant_t_0"), by = 1:nrow(f_dt_val)]
This function performs inference over each row of a folded data.table,
plots the results and gives metrics of the accuracy of the predictions. Given
that only a single row is predicted, the horizon of the prediction is at most 1.
This function is also called by the generic predict method for "dbn.fit"
objects. For long term forecasting, please refer to the
forecast_ts
function.
predict_dt(fit, dt, obj_nodes, verbose = T, look_ahead = F)
predict_dt(fit, dt, obj_nodes, verbose = T, look_ahead = F)
fit |
the fitted bn |
dt |
the test dataset |
obj_nodes |
the nodes that are going to be predicted. They are all predicted at the same time |
verbose |
if TRUE, displays the metrics and plots the real values against the predictions |
look_ahead |
boolean that defines whether or not the values of the variables in t_0 should be used when predicting, even if they are not present in obj_nodes. This decides if look-ahead bias is introduced or not. |
a data.table with the prediction results for each row
size = 3 data(motor) dt_train <- motor[200:900] dt_val <- motor[901:1000] # With a DBN obj <- c("pm_t_0") net <- learn_dbn_struc(dt_train, size) f_dt_train <- fold_dt(dt_train, size) f_dt_val <- fold_dt(dt_val, size) fit <- fit_dbn_params(net, f_dt_train, method = "mle-g") res <- suppressWarnings(predict_dt(fit, f_dt_val, obj_nodes = obj, verbose = FALSE)) # With a Gaussian BN directly from bnlearn obj <- c("pm") net <- bnlearn::mmhc(dt_train) fit <- bnlearn::bn.fit(net, dt_train, method = "mle-g") res <- suppressWarnings(predict_dt(fit, dt_val, obj_nodes = obj, verbose = FALSE))
size = 3 data(motor) dt_train <- motor[200:900] dt_val <- motor[901:1000] # With a DBN obj <- c("pm_t_0") net <- learn_dbn_struc(dt_train, size) f_dt_train <- fold_dt(dt_train, size) f_dt_val <- fold_dt(dt_val, size) fit <- fit_dbn_params(net, f_dt_train, method = "mle-g") res <- suppressWarnings(predict_dt(fit, f_dt_val, obj_nodes = obj, verbose = FALSE)) # With a Gaussian BN directly from bnlearn obj <- c("pm") net <- bnlearn::mmhc(dt_train) fit <- bnlearn::bn.fit(net, dt_train, method = "mle-g") res <- suppressWarnings(predict_dt(fit, dt_val, obj_nodes = obj, verbose = FALSE))
Generic method for predicting a dataset with a "dbn.fit" S3 objects. Calls
predict_dt
underneath.
## S3 method for class 'dbn.fit' predict(object, ...)
## S3 method for class 'dbn.fit' predict(object, ...)
object |
a "dbn.fit" object |
... |
additional parameters for the inference process |
a data.table with the prediction results
Generic print method for "dbn" S3 objects. Calls bnlearn's print underneath
## S3 method for class 'dbn' print(x, ...)
## S3 method for class 'dbn' print(x, ...)
x |
the "dbn" object |
... |
additional parameters |
Generic print method for "dbn.fit" S3 objects. Calls bnlearn's print underneath
## S3 method for class 'dbn.fit' print(x, ...)
## S3 method for class 'dbn.fit' print(x, ...)
x |
the "dbn.fit" object |
... |
additional parameters |
Generic method for "dbn.fit" S3 objects.
Calls bnlearn's rbn
underneath.
rbn.dbn.fit(x, n, ...)
rbn.dbn.fit(x, n, ...)
x |
the fitted network |
n |
number of samples |
... |
additional parameters |
the sampled dataset
In a time series dataset, there is a time difference between one row and the next one. This function reduces the number of rows from its current frequency to the desired one by averaging batches of rows. Instead of the frequency in Hz, the number of seconds between rows is asked (Hz = 1/s).
reduce_freq(dt, obj_freq, curr_freq, id_var = NULL)
reduce_freq(dt, obj_freq, curr_freq, id_var = NULL)
dt |
the original data.table |
obj_freq |
the desired number of seconds between rows |
curr_freq |
the number of seconds between rows in the original dataset |
id_var |
optional variable that labels different time series in a dataset, to avoid averaging values from different processes |
the data.table with the desired frequency
# Let's assume that the dataset has a frequency of 4Hz, 0.25 seconds between rows dt <- dbnR::motor dim(dt) # Let's change the frequency to 2Hz, 0.5 seconds between rows dt <- reduce_freq(dt, obj_freq = 0.5, curr_freq = 0.2) dim(dt)
# Let's assume that the dataset has a frequency of 4Hz, 0.25 seconds between rows dt <- dbnR::motor dim(dt) # Let's change the frequency to 2Hz, 0.5 seconds between rows dt <- reduce_freq(dt, obj_freq = 0.5, curr_freq = 0.2) dim(dt)
Generic method for "dbn.fit" S3 objects. Calls bnlearn underneath.
## S3 method for class 'dbn.fit' residuals(object, ...)
## S3 method for class 'dbn.fit' residuals(object, ...)
object |
the fitted network |
... |
additional parameters |
the residuals of fitting the network
Generic method for computing the score of a BN or a DBN.
Calls bnlearn's nodes
underneath.
I have to redefine the generic and mask the original for it to work on both
bn and dbn objects without the user having to import bnlearn.
score(object, ...)
score(object, ...)
object |
a "bn" or "dbn" object |
... |
additional parameters |
the score of the network
This function moves the values in t_0, t_1, ..., t_n-1 in a folded dataset row to t_1, t_2, ..., t_n. All the variables in t_0 will be inputed with NAs and the obtained row can be used to forecast up to any desired point.
shift_values(f_dt, row)
shift_values(f_dt, row)
f_dt |
a folded dataset |
row |
the index of the row that is going to be processed |
a one row data.table the shifted values
dt <- dbnR::motor f_dt <- dbnR::fold_dt(dt, size = 2) s_row <- dbnR::shift_values(f_dt, row = 500)
dt <- dbnR::motor f_dt <- dbnR::fold_dt(dt, size = 2) s_row <- dbnR::shift_values(f_dt, row = 500)
Generic method for "dbn.fit" S3 objects. Calls bnlearn underneath.
## S3 method for class 'dbn.fit' sigma(object, ...)
## S3 method for class 'dbn.fit' sigma(object, ...)
object |
the fitted network |
... |
additional parameters |
the standard deviation residuals of fitting the network
Given a dbn.fit object, the size of the net and a folded dataset, performs a smoothing of a trajectory. Smoothing is the opposite of forecasting: given a starting point, predict backwards in time to obtain the time series that generated that point.
smooth_ts( dt, fit, size = NULL, obj_vars, ini = dim(dt)[1], len = ini - 1, print_res = TRUE, plot_res = TRUE, prov_ev = NULL )
smooth_ts( dt, fit, size = NULL, obj_vars, ini = dim(dt)[1], len = ini - 1, print_res = TRUE, plot_res = TRUE, prov_ev = NULL )
dt |
data.table object with the TS data |
fit |
dbn.fit object |
size |
number of time slices of the net. Deprecated, will be removed in the future |
obj_vars |
variables to be predicted. Should be in the oldest time step |
ini |
starting point in the dataset to smooth |
len |
length of the smoothing |
print_res |
if TRUE prints the mae and sd metrics of the smoothing |
plot_res |
if TRUE plots the results of the smoothing |
prov_ev |
variables to be provided as evidence in each smoothing step. Should be in the oldest time step |
a list with the original values and the results of the smoothing
size = 3 data(motor) dt_train <- motor[200:900] dt_val <- motor[901:1000] obj <- c("pm_t_2") net <- learn_dbn_struc(dt_train, size) f_dt_train <- fold_dt(dt_train, size) f_dt_val <- fold_dt(dt_val, size) fit <- fit_dbn_params(net, f_dt_train, method = "mle-g") res <- suppressWarnings(smooth_ts(f_dt_val, fit, obj_vars = obj, len = 10, print_res = FALSE, plot_res = FALSE))
size = 3 data(motor) dt_train <- motor[200:900] dt_val <- motor[901:1000] obj <- c("pm_t_2") net <- learn_dbn_struc(dt_train, size) f_dt_train <- fold_dt(dt_train, size) f_dt_val <- fold_dt(dt_val, size) fit <- fit_dbn_params(net, f_dt_train, method = "mle-g") res <- suppressWarnings(smooth_ts(f_dt_val, fit, obj_vars = obj, len = 10, print_res = FALSE, plot_res = FALSE))
This will rename the columns in a data.table so that they end in '_t_0', which will be needed when folding the data.table. If any of the columns already ends in '_t_0', a warning will be issued and no further operation will be done. There is no need to use this function to learn a DBN unless some operation with the variable names wants to be done prior to folding a dataset.
time_rename(dt)
time_rename(dt)
dt |
the data.table to be treated |
the renamed data.table
data("motor") dt <- time_rename(motor)
data("motor") dt <- time_rename(motor)