C++

Class Declaration

class bqreg_t

Bayesian Quantile Regression class.

Public Functions

~bqreg_t() = default

Default destructor.

bqreg_t(const bqreg_t &obj_inp)

Constructor using a bqreg_t object to copy from.

Parameters

obj_inp – an object of type bqreg_t whose elements will be used to initalize a new object of the same type.

bqreg_t(bqreg_t &&obj_inp)

Constructor using a moveable bqreg_t object.

Parameters

obj_inp – a moveable object of type bqreg_t whose elements will be used to initalize a new object of the same type.

inline explicit bqreg_t(const ColVec_t &Y_inp, const Mat_t &X_inp)

Constructor using data matrices.

Parameters
  • Y_inp – an n x 1 vector defining the target variable

  • X_inp – an n x K matrix of features

inline explicit bqreg_t(ColVec_t &&Y_inp, Mat_t &&X_inp)

Constructor using moveable data matrices.

Parameters
  • Y_inp – an n x 1 vector defining the target variable

  • X_inp – an n x K matrix of features

inline bqreg_t &operator=(const bqreg_t &obj_inp)

Assignment operator.

Parameters

obj_inp – an object of type bqreg_t whose elements will be used to initalize a new object of the same type.

Returns

an object of type bqreg_t

inline bqreg_t &operator=(bqreg_t &&obj_inp)

Assignment operator.

Parameters

obj_inp – a moveable object of type bqreg_t whose elements will be used to initalize a new object of the same type.

Returns

an object of type bqreg_t

inline int get_omp_n_threads()

OpenMP threads check.

Returns

the number of OpenMP threads to use in the Gibbs sampler stage.

inline void set_omp_n_threads(const int omp_n_threads_inp)

OpenMP threads set.

Parameters

omp_n_threads_inp – the number of OpenMP threads to use in the Gibbs sampler stage.

inline void set_seed_value(const size_t seed_val_inp)

RNG engine seeding.

Parameters

seed_val_inp – the seed value for the RNG engine.

inline void load_data(const ColVec_t &Y_inp, const Mat_t &X_inp)

Load data.

Parameters
  • Y_inp – an n x 1 vector defining the target variable

  • X_inp – an n x K matrix of features

inline void set_quantile_target(const fp_t tau_inp)

Set the target quantile value.

Parameters

tau_inp – a real value between zero and one

inline void set_prior_params(const ColVec_t &prior_beta_mean_inp, const Mat_t &prior_beta_var_inp, const fp_t prior_sigma_shape_inp, const fp_t prior_sigma_scale_inp)

Set the Prior.

Set the parameters of the distributions that define the prior

Parameters
  • prior_beta_mean_inp – mean of the prior distribution for \( \beta \)

  • prior_beta_var_inp – variance of the prior distribution for \( \beta \)

  • prior_sigma_shape_inp – shape parameter of the prior distribution for \( \sigma \)

  • prior_sigma_scale_inp – scale parameter of the prior distribution for \( \sigma \)

inline ColVec_t get_initial_beta_draw()

Initial value check for the Gibbs sampler.

Returns

a vector defining the initial values for \( \beta \) in the Gibbs sampler.

inline void set_initial_beta_draw(const ColVec_t &beta_initial_draw_inp)

Set the initial values for the Gibbs sampler.

Parameters

beta_initial_draw_inp – a vector defining the initial values for \( \beta \) in the Gibbs sampler.

inline void gibbs(const size_t n_burnin_draws, const size_t n_keep_draws, const size_t thinning_factor, Mat_t &beta_draws, Mat_t &z_draws, ColVec_t &sigma_draws)

Run the Gibbs sampler.

Parameters
  • n_burnin_draws – the number of burnin draws

  • n_keep_draws – the number of draws to keep, post burnin

  • thinning_factor – the number of draws to skip between keep draws

  • beta_draws – a writable matrix to store the draws of \( \beta \)

  • z_draws – a writable matrix to store the draws of \( z \)

  • sigma_draws – a writable vector to store the draws of \( \sigma \)

Public Members

ColVec_t Y

An n x 1 vector defining the target variable

Mat_t X

An n x K matrix of features

fp_t tau

The target quantile value

ColVec_t prior_beta_mean

Mean of the prior distribution for \( \beta \)

Mat_t prior_beta_var

Variance of the prior distribution for \( \beta \)

fp_t prior_sigma_shape

Shape parameter of the prior distribution for \( \sigma \)

fp_t prior_sigma_scale

Scale parameter of the prior distribution for \( \sigma \)


Example

Code to run a median regression example is given below.

#include "bqreg.hpp"

inline
Eigen::MatrixXd
eigen_randn(const size_t nr, const size_t nc)
{
    static std::mt19937 gen{ std::random_device{}() };
    static std::normal_distribution<> dist;

    return Eigen::MatrixXd{ nr, nc }.unaryExpr([&](double x) { (void)(x); return dist(gen); });
}

int main()
{
    // generate data
    const int n = 500;
    const int K = 3;

    Eigen::MatrixXd X = eigen_randn(n, K);
    X.col(0).setConstant(1);

    Eigen::VectorXd beta0(3);
    beta0 << 5.0, 1.3, 1.8;

    Eigen::VectorXd Y = (X * beta0).array() + eigen_randn(n, 1).array();

    // initialize object
    bqreg::bqreg_t obj = bqreg::bqreg_t(Y, X);

    // set prior pars
    Eigen::VectorXd beta_bar = Eigen::VectorXd::Zero(K);

    Eigen::MatrixXd V0 = Eigen::MatrixXd::Zero(K, K);
    V0.diagonal().array() = 1.0 / 0.001;

    int prior_shape = 3;
    int prior_scale = 3;

    obj.set_prior_params(beta_bar, V0, prior_shape, prior_scale);

    // (optional) set the number of OpenMP threads
    obj.set_omp_n_threads(4);
    std::cout << "Number of OpenMP threads: " << obj.get_omp_n_threads() << std::endl;

    // (optional) set the seed value of the Gibbs sampler RNG
    obj.set_seed_value(1111);

    // set the target quantile
    obj.set_quantile_target(0.5);

    // run Gibbs sampler
    size_t n_burnin_draws = 10000;
    size_t n_keep_draws = 10000;
    size_t thinning_factor = 0;

    Eigen::MatrixXd beta_draws;
    Eigen::MatrixXd z_draws;
    Eigen::VectorXd sigma_draws;

    obj.gibbs(n_burnin_draws, n_keep_draws, thinning_factor, beta_draws, z_draws, sigma_draws);

    //

    std::cout << "\n - beta posterior mean:\n" << beta_draws.rowwise().mean() << std::endl;

    //

    return 0;
}