MuGen
Multitrait genetics
Public Member Functions | Protected Attributes | List of all members
MVnormBetaFt Class Reference

Regression with multiple predictors. More...

#include <MuGen.h>

Inheritance diagram for MVnormBetaFt:
[legend]
Collaboration diagram for MVnormBetaFt:
[legend]

Public Member Functions

 MVnormBetaFt ()
 Default constructor.
 
 MVnormBetaFt (gsl_vector *b, const gsl_vector *sd, gsl_matrix *pred, const size_t &iCl, gsl_matrix *allFt, const size_t &begRw, const gsl_rng *r)
 Univariate random constructor with vector. More...
 
 MVnormBetaFt (gsl_vector *b, const gsl_vector *sd, gsl_matrix *pred, const size_t &iCl, gsl_matrix *allFt, const size_t &begRw, const gsl_rng *r, const size_t &up)
 Univariate random constructor with vector and prior index. More...
 
 MVnormBetaFt (gsl_vector *b, const gsl_matrix *Sig, gsl_matrix *pred, const size_t &iCl, gsl_matrix *allFt, const size_t &begRw, const gsl_rng *r)
 Multivariate random constructor with vector. More...
 
 MVnormBetaFt (gsl_vector *b, const gsl_matrix *Sig, gsl_matrix *pred, const size_t &iCl, gsl_matrix *allFt, const size_t &begRw, const gsl_rng *r, const size_t &up)
 Multivariate random constructor with vector and prior index. More...
 
 MVnormBetaFt (const gsl_matrix *resp, gsl_matrix *pred, const size_t &iCl, gsl_matrix *allFt, const size_t &begRw, const gsl_matrix *Sig, const gsl_rng *r, gsl_matrix *bet, const size_t &iRw)
 Multivariate random constructor with response matrix. More...
 
 MVnormBetaFt (const gsl_matrix *resp, gsl_matrix *pred, const size_t &iCl, gsl_matrix *allFt, const size_t &begRw, const gsl_matrix *Sig, const gsl_rng *r, const size_t &up, gsl_matrix *bet, const size_t &iRw)
 Multivariate random constructor with response matrix and prior index. More...
 
 MVnormBetaFt (const gsl_matrix *resp, gsl_matrix *pred, const size_t &iCl, vector< double > &eaFt, const gsl_matrix *Sig, const gsl_rng *r, gsl_matrix *bet, const size_t &iRw)
 Multivariate random constructor with response matrix. More...
 
 MVnormBetaFt (const gsl_matrix *resp, gsl_matrix *pred, const size_t &iCl, vector< double > &eaFt, const gsl_matrix *Sig, const gsl_rng *r, const size_t &up, gsl_matrix *bet, const size_t &iRw)
 Multivariate random constructor with response matrix and prior index. More...
 
 MVnormBetaFt (const Grp &resp, gsl_matrix *pred, const size_t &iCl, gsl_matrix *allFt, const size_t &begRw, const gsl_matrix *Sig, const gsl_rng *r, gsl_matrix *bet, const size_t &iRw)
 Multivariate random constructor with Grp type response. More...
 
 MVnormBetaFt (const Grp &resp, gsl_matrix *pred, const size_t &iCl, gsl_matrix *allFt, const size_t &begRw, const gsl_matrix *Sig, const gsl_rng *r, const size_t &up, gsl_matrix *bet, const size_t &iRw)
 Multivariate random constructor with Grp type response and prior index. More...
 
 MVnormBetaFt (const Grp &resp, gsl_matrix *pred, const size_t &iCl, vector< double > &eaFt, const gsl_matrix *Sig, const gsl_rng *r, gsl_matrix *bet, const size_t &iRw)
 Multivariate random constructor with Grp type response. More...
 
 MVnormBetaFt (const Grp &resp, gsl_matrix *pred, const size_t &iCl, vector< double > &eaFt, const gsl_matrix *Sig, const gsl_rng *r, const size_t &up, gsl_matrix *bet, const size_t &iRw)
 Multivariate random constructor with Grp type response and prior index. More...
 
 MVnormBetaFt (gsl_matrix *pred, const size_t &iCl, vector< double > &eaFt, gsl_matrix *bet, const size_t &iRw)
 Deterministic constructor. More...
 
 MVnormBetaFt (gsl_matrix *pred, const size_t &iCl, vector< double > &eaFt, const size_t &up, gsl_matrix *bet, const size_t &iRw)
 Deterministic constructor with a prior index. More...
 
 MVnormBetaFt (const MVnormBetaFt &)
 Deterministic copy constructor. More...
 
MVnormBetaFtoperator= (const MVnormBetaFt &)
 Deterministic assignment operator. More...
 
virtual ~MVnormBetaFt ()
 Destructor.
 
void update (const Grp &resp, const SigmaI &SigIb, const gsl_rng *r)
 Gaussian likelihood. More...
 
void update (const Grp &resp, const Qgrp &q, const SigmaI &SigIb, const gsl_rng *r)
 Sudent- \(t\) likelihood. More...
 
void update (const Grp &resp, const SigmaI &SigIb, const SigmaI &SigIp, const gsl_rng *r)
 Gaussian likelihood, Gaussian prior. More...
 
void update (const Grp &resp, const SigmaI &SigIb, const double &qPr, const SigmaI &SigIp, const gsl_rng *r)
 Gaussian likelihood, Student- \(t\) prior. More...
 
void update (const Grp &resp, const Qgrp &q, const SigmaI &SigIb, const SigmaI &SigIp, const gsl_rng *r)
 Student- \(t\) likelihood, Gaussian prior. More...
 
void update (const Grp &resp, const Qgrp &q, const SigmaI &SigIb, const double &qPr, const SigmaI &SigIp, const gsl_rng *r)
 Student- \(t\) likelihood, Student- \(t\) prior. More...
 
void update (const Grp &resp, const SigmaI &SigIb, const Grp &muPr, const SigmaI &SigIp, const gsl_rng *r)
 Gaussian likelihood, Gaussian prior. More...
 
void update (const Grp &resp, const SigmaI &SigIb, const Grp &muPr, const double &qPr, const SigmaI &SigIp, const gsl_rng *r)
 Gaussian likelihood, Student- \(t\) prior. More...
 
void update (const Grp &resp, const Qgrp &q, const SigmaI &SigIb, const Grp &muPr, const SigmaI &SigIp, const gsl_rng *r)
 Student- \(t\) likelihood, Gaussian prior. More...
 
void update (const Grp &resp, const Qgrp &q, const SigmaI &SigIb, const Grp &muPr, const double &qPr, const SigmaI &SigIp, const gsl_rng *r)
 Student- \(t\) likelihood, Student- \(t\) prior. More...
 
- Public Member Functions inherited from MVnormBeta
 MVnormBeta ()
 Default constructor.
 
 MVnormBeta (const size_t d)
 Dimension-only constructor. More...
 
 MVnormBeta (gsl_vector *b, const gsl_vector *sd, gsl_matrix *pred, const size_t &iCl, const gsl_rng *r)
 Univariate random constructor with vector. More...
 
 MVnormBeta (gsl_vector *b, const gsl_vector *sd, gsl_matrix *pred, const size_t &iCl, const gsl_rng *r, const size_t &up)
 Univariate random constructor with vector and pointer to prior. More...
 
 MVnormBeta (gsl_vector *b, const gsl_matrix *Sig, gsl_matrix *pred, const size_t &iCl, const gsl_rng *r)
 Multivariate random constructor with vector. More...
 
 MVnormBeta (gsl_vector *b, const gsl_matrix *Sig, gsl_matrix *pred, const size_t &iCl, const gsl_rng *r, const size_t &up)
 Multivariate random constructor with vector and pointer to prior. More...
 
 MVnormBeta (const gsl_matrix *resp, gsl_matrix *pred, const size_t &iCl, const gsl_matrix *Sig, const gsl_rng *r, gsl_matrix *bet, const size_t &iRw)
 Multivariate constructor with response matrix. More...
 
 MVnormBeta (const gsl_matrix *resp, gsl_matrix *pred, const size_t &iCl, const gsl_matrix *Sig, const gsl_rng *r, const size_t &up, gsl_matrix *bet, const size_t &iRw)
 Multivariate constructor with response matrix and index for a prior. More...
 
 MVnormBeta (const Grp &resp, gsl_matrix *pred, const size_t &iCl, const gsl_matrix *Sig, const gsl_rng *r, gsl_matrix *bet, const size_t &iRw)
 Multivariate constructor with a Grp type response. More...
 
 MVnormBeta (const Grp &resp, gsl_matrix *pred, const size_t &iCl, const gsl_matrix *Sig, const gsl_rng *r, const size_t &up, gsl_matrix *bet, const size_t &iRw)
 Multivariate constructor with a Grp type response and index for a prior. More...
 
 MVnormBeta (gsl_matrix *pred, const size_t &iCl, gsl_matrix *bet, const size_t &iRw)
 Deterministic constructor. More...
 
 MVnormBeta (gsl_matrix *pred, const size_t &iCl, const size_t &up, gsl_matrix *bet, const size_t &iRw)
 Deterministic constructor with prior index. More...
 
 MVnormBeta (const MVnormBeta &)
 Deterministic copy constructor. More...
 
MVnormBetaoperator= (const MVnormBeta &)
 Deterministic assignement operator. More...
 
virtual ~MVnormBeta ()
 Destructor.
 
const size_t * up () const
 Points to the prior. More...
 
double scalePar () const
 Scale parameter. More...
 
- Public Member Functions inherited from MVnorm
 MVnorm (const MVnorm &)
 Copy constructor. More...
 
MVnormoperator= (const MVnorm &)
 Assignement operator. More...
 
virtual ~MVnorm ()
 Virtual destructor. More...
 
virtual double mhl (const MVnorm *x, const SigmaI &SigI)
 Mahalanobis distance to a vector. More...
 
virtual double mhl (const MVnorm *x, const SigmaI &SigI) const
 Mahalanobis distance to a vector. More...
 
virtual double mhl (const gsl_vector *x, const SigmaI &SigI)
 Mahalanobis distance to a vector. More...
 
virtual double mhl (const gsl_vector *x, const SigmaI &SigI) const
 Mahalanobis distance to a vector. More...
 
virtual double mhl (const SigmaI &SigI)
 Mahalanobis distance to zero. More...
 
virtual double mhl (const SigmaI &SigI) const
 Mahalanobis distance to zero. More...
 
double density (const gsl_vector *theta, const SigmaI &SigI)
 Multivariate Gaussian density. More...
 
double density (const gsl_vector *theta, const SigmaI &SigI) const
 Multivariate Gaussian density. More...
 
double density (const MVnorm *theta, const SigmaI &SigI)
 Multivariate Gaussian density. More...
 
double density (const MVnorm *theta, const SigmaI &SigI) const
 Multivariate Gaussian density. More...
 
void save (const string &fileNam, const char *how="a")
 Save function. More...
 
void save (FILE *fileStr)
 Save function. More...
 
double operator[] (const size_t i) const
 Subscript operator. More...
 
void valSet (const size_t i, const double x)
 Setting an element to a value. More...
 
const gsl_vector * getVec () const
 Access the location vector. More...
 
size_t len () const
 Length of the location vector. More...
 
virtual size_t nMissP () const
 Number of missing values. More...
 
virtual const vector< size_t > getMisPhen () const
 Indexes of missing values. More...
 
virtual const vector< size_t > * down () const
 Points to the corresponding data. More...
 

Protected Attributes

gsl_matrix_view _fitted
 Matrix of already fitted values. More...
 
- Protected Attributes inherited from MVnormBeta
gsl_vector_view _X
 Predictor. More...
 
double _scale
 Scale parameter. More...
 
size_t _N
 Length of the predictor.
 
const size_t * _upLevel
 Row index of the prior. More...
 
- Protected Attributes inherited from MVnorm
gsl_vector_view _vec
 Data vector. More...
 
size_t _d
 Length of the data vector.
 

Additional Inherited Members

- Protected Member Functions inherited from MVnorm
 MVnorm ()
 Default constructor. More...
 
 MVnorm (const size_t &d)
 Dimension-only constructor. More...
 
 MVnorm (gsl_vector *mn)
 Dimension and vector value constructor. More...
 
 MVnorm (gsl_vector *mn, const gsl_vector *sd, const gsl_rng *r)
 Univariate Gaussian constructor. More...
 
 MVnorm (gsl_vector *mn, const gsl_matrix *Sig, const gsl_rng *r)
 Multivariate Gaussian constructor. More...
 
 MVnorm (gsl_matrix *mn, const size_t &iRw)
 Dimension and vector value constructor. More...
 
 MVnorm (gsl_matrix *mn, const size_t &iRw, const gsl_vector *sd, const gsl_rng *r)
 Univariate Gaussian constructor with a matrix. More...
 
 MVnorm (gsl_matrix *mn, const size_t &iRw, const gsl_matrix *Sig, const gsl_rng *r)
 Multivariate Gaussian constructor with a matrix. More...
 

Detailed Description

Regression with multiple predictors.

Full implementation of a multiple regression using one-at-a-time updating. Updates regression coefficients for the current predictor, while accounting for the effects of other predictors

Constructor & Destructor Documentation

◆ MVnormBetaFt() [1/15]

MVnormBetaFt::MVnormBetaFt ( gsl_vector *  b,
const gsl_vector *  sd,
gsl_matrix *  pred,
const size_t &  iCl,
gsl_matrix *  allFt,
const size_t &  begRw,
const gsl_rng *  r 
)
inline

Univariate random constructor with vector.

Sets up _vec to point to the provided vector of values, and the predictor to a column in the provided matrix of predictors.
The _fitted member points to the indicated submatrix of a matrix that contains all the partial fitted matrices.

Parameters
[in]gsl_vector*vector of values
[in]gsl_vector*vector of standard deviations
[in]gsl_matrix*matrix of predictors
[in]size_t&column index for the predictor matrix
[in]gsl_matrix*matrix with all the partial fitted matrices stacked
[in]size_t&index of the row where the relevant fitted matrix begins
[in]gsl_rng*pointer to a PNG

◆ MVnormBetaFt() [2/15]

MVnormBetaFt::MVnormBetaFt ( gsl_vector *  b,
const gsl_vector *  sd,
gsl_matrix *  pred,
const size_t &  iCl,
gsl_matrix *  allFt,
const size_t &  begRw,
const gsl_rng *  r,
const size_t &  up 
)
inline

Univariate random constructor with vector and prior index.

Sets up _vec to point to the provided vector of values, and the predictor to a column in the provided matrix of predictors. The _fitted member points to the indicated submatrix of a matrix that contains all the partial fitted matrices. The _upLevel member is set to point to a row in the prior matrix.

Parameters
[in]gsl_vector*vector of values
[in]gsl_vector*vector of standard deviations
[in]gsl_matrix*matrix of predictors
[in]size_t&column index for the predictor matrix
[in]gsl_matrix*matrix with all the partial fitted matrices stacked
[in]size_t&index of the row where the relevant fitted matrix begins
[in]gsl_rng*pointer to a PNG
[in]size_t&row index of the prior matrix

◆ MVnormBetaFt() [3/15]

MVnormBetaFt::MVnormBetaFt ( gsl_vector *  b,
const gsl_matrix *  Sig,
gsl_matrix *  pred,
const size_t &  iCl,
gsl_matrix *  allFt,
const size_t &  begRw,
const gsl_rng *  r 
)
inline

Multivariate random constructor with vector.

Sets up _vec to point to the provided vector of values, and the predictor to a column in the provided matrix of predictors. The _fitted member points to the indicated submatrix of a matrix that contains all the partial fitted matrices.

Parameters
[in]gsl_vector*vector of values
[in]gsl_matrix*covariance matrix
[in]gsl_matrix*matrix of predictors
[in]size_t&column index for the predictor matrix
[in]gsl_matrix*matrix with all the partial fitted matrices stacked
[in]size_t&index of the row where the relevant fitted matrix begins
[in]gsl_rng*pointer to a PNG

◆ MVnormBetaFt() [4/15]

MVnormBetaFt::MVnormBetaFt ( gsl_vector *  b,
const gsl_matrix *  Sig,
gsl_matrix *  pred,
const size_t &  iCl,
gsl_matrix *  allFt,
const size_t &  begRw,
const gsl_rng *  r,
const size_t &  up 
)
inline

Multivariate random constructor with vector and prior index.

Sets up _vec to point to the provided vector of values, and the predictor to a column in the provided matrix of predictors. The _fitted member points to the indicated submatrix of a matrix that contains all the partial fitted matrices. The _upLevel member is set to point to a row in the prior matrix.

Parameters
[in]gsl_vector*vector of values
[in]gsl_matrix*covariance matrix
[in]gsl_matrix*matrix of predictors
[in]size_t&column index for the predictor matrix
[in]gsl_matrix*matrix with all the partial fitted matrices stacked
[in]size_t&index of the row where the relevant fitted matrix begins
[in]gsl_rng*pointer to a PNG
[in]size_t&row index of the prior matrix

◆ MVnormBetaFt() [5/15]

MVnormBetaFt::MVnormBetaFt ( const gsl_matrix *  resp,
gsl_matrix *  pred,
const size_t &  iCl,
gsl_matrix *  allFt,
const size_t &  begRw,
const gsl_matrix *  Sig,
const gsl_rng *  r,
gsl_matrix *  bet,
const size_t &  iRw 
)
inline

Multivariate random constructor with response matrix.

Sets up _vec to point to the provided vector of values, and the predictor to a column in the provided matrix of predictors. Initializes regression coefficients using the provided response matrix. The _fitted member points to the indicated submatrix of a matrix that contains all the partial fitted matrices.

Parameters
[in]gsl_matrix*response matrix
[in]gsl_matrix*predictor matrix
[in]size_t&column index for the predictor matrix
[in]gsl_matrix*matrix with all the partial fitted matrices stacked
[in]size_t&index of the row where the relevant fitted matrix begins
[in]gsl_matrix*covariance matrix
[in]gsl_rng*pointer to a PNG
[in]gsl_matrix*matrix of regression coefficients
[in]size_t&row index of the regression coefficient matrix

◆ MVnormBetaFt() [6/15]

MVnormBetaFt::MVnormBetaFt ( const gsl_matrix *  resp,
gsl_matrix *  pred,
const size_t &  iCl,
gsl_matrix *  allFt,
const size_t &  begRw,
const gsl_matrix *  Sig,
const gsl_rng *  r,
const size_t &  up,
gsl_matrix *  bet,
const size_t &  iRw 
)
inline

Multivariate random constructor with response matrix and prior index.

Sets up _vec to point to the provided vector of values, and the predictor to a column in the provided matrix of predictors. Initializes regression coefficients using the provided response matrix. The _fitted member points to the indicated submatrix of a matrix that contains all the partial fitted matrices. The _upLevel member is set to point to a row in the prior matrix.

Parameters
[in]gsl_matrix*response matrix
[in]gsl_matrix*predictor matrix
[in]size_t&column index for the predictor matrix
[in]gsl_matrix*matrix with all the partial fitted matrices stacked
[in]size_t&index of the row where the relevant fitted matrix begins
[in]gsl_matrix*covariance matrix
[in]gsl_rng*pointer to a PNG
[in]size_t&row index of the prior matrix
[in]gsl_matrix*matrix of regression coefficients
[in]size_t&row index of the regression coefficient matrix

◆ MVnormBetaFt() [7/15]

MVnormBetaFt::MVnormBetaFt ( const gsl_matrix *  resp,
gsl_matrix *  pred,
const size_t &  iCl,
vector< double > &  eaFt,
const gsl_matrix *  Sig,
const gsl_rng *  r,
gsl_matrix *  bet,
const size_t &  iRw 
)
inline

Multivariate random constructor with response matrix.

Sets up _vec to point to the provided vector of values, and the predictor to a column in the provided matrix of predictors. Initializes regression coefficients using the provided response matrix. The _fitted member points to the indicated vector that has the appropriate partial fitted values.

Parameters
[in]gsl_matrix*response matrix
[in]gsl_matrix*predictor matrix
[in]size_t&column index for the predictor matrix
[in]vector<double>&vectorized partial fitted matrix
[in]gsl_matrix*covariance matrix
[in]gsl_rng*pointer to a PNG
[in]gsl_matrix*matrix of regression coefficients
[in]size_t&row index of the regression coefficient matrix

◆ MVnormBetaFt() [8/15]

MVnormBetaFt::MVnormBetaFt ( const gsl_matrix *  resp,
gsl_matrix *  pred,
const size_t &  iCl,
vector< double > &  eaFt,
const gsl_matrix *  Sig,
const gsl_rng *  r,
const size_t &  up,
gsl_matrix *  bet,
const size_t &  iRw 
)
inline

Multivariate random constructor with response matrix and prior index.

Sets up _vec to point to the provided vector of values, and the predictor to a column in the provided matrix of predictors. Initializes regression coefficients using the provided response matrix. The _fitted member points to the indicated vector that has the appropriate partial fitted values. The _upLevel member is set to point to a row in the prior matrix.

Parameters
[in]gsl_matrix*response matrix
[in]gsl_matrix*predictor matrix
[in]size_t&column index for the predictor matrix
[in]vector<double>&vectorized partial fitted matrix
[in]gsl_matrix*covariance matrix
[in]gsl_rng*pointer to a PNG
[in]size_t&row index of the prior matrix
[in]gsl_matrix*matrix of regression coefficients
[in]size_t&row index of the regression coefficient matrix

◆ MVnormBetaFt() [9/15]

MVnormBetaFt::MVnormBetaFt ( const Grp resp,
gsl_matrix *  pred,
const size_t &  iCl,
gsl_matrix *  allFt,
const size_t &  begRw,
const gsl_matrix *  Sig,
const gsl_rng *  r,
gsl_matrix *  bet,
const size_t &  iRw 
)
inline

Multivariate random constructor with Grp type response.

Sets up _vec to point to the provided vector of values, and the predictor to a column in the provided matrix of predictors. Initializes regression coefficients using the provided response matrix. The _fitted member points to the indicated submatrix of a matrix that contains all the partial fitted matrices.

Parameters
[in]Grp&response
[in]gsl_matrix*predictor matrix
[in]size_t&column index for the predictor matrix
[in]gsl_matrix*matrix with all the partial fitted matrices stacked
[in]size_t&index of the row where the relevant fitted matrix begins
[in]gsl_matrix*covariance matrix
[in]gsl_rng*pointer to a PNG
[in]gsl_matrix*matrix of regression coefficients
[in]size_t&row index of the regression coefficient matrix

◆ MVnormBetaFt() [10/15]

MVnormBetaFt::MVnormBetaFt ( const Grp resp,
gsl_matrix *  pred,
const size_t &  iCl,
gsl_matrix *  allFt,
const size_t &  begRw,
const gsl_matrix *  Sig,
const gsl_rng *  r,
const size_t &  up,
gsl_matrix *  bet,
const size_t &  iRw 
)
inline

Multivariate random constructor with Grp type response and prior index.

Sets up _vec to point to the provided vector of values, and the predictor to a column in the provided matrix of predictors. Initializes regression coefficients using the provided response matrix. The _fitted member points to the indicated submatrix of a matrix that contains all the partial fitted matrices. The _upLevel member is set to point to a row in the prior matrix.

Parameters
[in]Grp&response
[in]gsl_matrix*predictor matrix
[in]size_t&column index for the predictor matrix
[in]gsl_matrix*matrix with all the partial fitted matrices stacked
[in]size_t&index of the row where the relevant fitted matrix begins
[in]gsl_matrix*covariance matrix
[in]gsl_rng*pointer to a PNG
[in]size_t&row index of the prior matrix
[in]gsl_matrix*matrix of regression coefficients
[in]size_t&row index of the regression coefficient matrix

◆ MVnormBetaFt() [11/15]

MVnormBetaFt::MVnormBetaFt ( const Grp resp,
gsl_matrix *  pred,
const size_t &  iCl,
vector< double > &  eaFt,
const gsl_matrix *  Sig,
const gsl_rng *  r,
gsl_matrix *  bet,
const size_t &  iRw 
)
inline

Multivariate random constructor with Grp type response.

Sets up _vec to point to the provided vector of values, and the predictor to a column in the provided matrix of predictors. Initializes regression coefficients using the provided response matrix. The _fitted member points to the indicated vector that has the appropriate partial fitted values.

Parameters
[in]Grp&response
[in]gsl_matrix*predictor matrix
[in]size_t&column index for the predictor matrix
[in]vector<double>&vectorized partial fitted matrix
[in]gsl_matrix*covariance matrix
[in]gsl_rng*pointer to a PNG
[in]gsl_matrix*matrix of regression coefficients
[in]size_t&row index of the regression coefficient matrix

◆ MVnormBetaFt() [12/15]

MVnormBetaFt::MVnormBetaFt ( const Grp resp,
gsl_matrix *  pred,
const size_t &  iCl,
vector< double > &  eaFt,
const gsl_matrix *  Sig,
const gsl_rng *  r,
const size_t &  up,
gsl_matrix *  bet,
const size_t &  iRw 
)
inline

Multivariate random constructor with Grp type response and prior index.

Sets up _vec to point to the provided vector of values, and the predictor to a column in the provided matrix of predictors. Initializes regression coefficients using the provided response matrix. The _fitted member points to the indicated vector that has the appropriate partial fitted values. The _upLevel member is set to point to a row in the prior matrix.

Parameters
[in]Grp&response
[in]gsl_matrix*predictor matrix
[in]size_t&column index for the predictor matrix
[in]vector<double>&vectorized partial fitted matrix
[in]gsl_matrix*covariance matrix
[in]gsl_rng*pointer to a PNG
[in]size_t&row index of the prior matrix
[in]gsl_matrix*matrix of regression coefficients
[in]size_t&row index of the regression coefficient matrix

◆ MVnormBetaFt() [13/15]

MVnormBetaFt::MVnormBetaFt ( gsl_matrix *  pred,
const size_t &  iCl,
vector< double > &  eaFt,
gsl_matrix *  bet,
const size_t &  iRw 
)
inline

Deterministic constructor.

Does not initialize the regression coefficients, but simply points to the already-initialized matrix of values. The _fitted member points to the indicated vector that has the appropriate partial fitted values.

Parameters
[in]gsl_matrix*predictor matrix
[in]size_t&column index of the predictor matrix
[in]vector<double>&vectorized partial fitted matrix
[in]gsl_matrix*regression coefficient matrix
[in]size_t&row index of the regression coefficient matrix

◆ MVnormBetaFt() [14/15]

MVnormBetaFt::MVnormBetaFt ( gsl_matrix *  pred,
const size_t &  iCl,
vector< double > &  eaFt,
const size_t &  up,
gsl_matrix *  bet,
const size_t &  iRw 
)
inline

Deterministic constructor with a prior index.

Does not initialize the regression coefficients, but simply points to the already-initialized matrix of values. The _fitted member points to the indicated vector that has the appropriate partial fitted values. The _upLevel member is set to point to a row in the prior matrix.

Parameters
[in]gsl_matrix*predictor matrix
[in]size_t&column index of the predictor matrix
[in]vector<double>&vectorized partial fitted matrix
[in]size_t&row index of the prior matrix
[in]gsl_matrix*regression coefficient matrix
[in]size_t&row index of the regression coefficient matrix

◆ MVnormBetaFt() [15/15]

MVnormBetaFt::MVnormBetaFt ( const MVnormBetaFt b)

Deterministic copy constructor.

Parameters
[in]MVnormBetaFt&object to be copied

Member Function Documentation

◆ operator=()

MVnormBetaFt & MVnormBetaFt::operator= ( const MVnormBetaFt b)

Deterministic assignment operator.

Parameters
[in]MVnormBetaFt&object to be copied
Returns
Refernce to an object of class MVnormBetaFt

◆ update() [1/10]

void MVnormBetaFt::update ( const Grp ,
const Qgrp ,
const SigmaI ,
const double &  ,
const SigmaI ,
const gsl_rng *   
)
virtual

Student- \(t\) likelihood, Student- \(t\) prior.

Parameters
[in]Grp&data for the likelihood
[in]Qgrp&vector Student- \(t\) covariance scale parameter for the likelihood covariance
[in]SigmaI&inverse-covariance matrix for the likelihood
[in]double&Student- \(t\) scale parameter for the prior covariance
[in]SigmaI&prior inverse-covariance matrix
[in]gsl_rng*pointer to a PNG

Reimplemented from MVnormBeta.

◆ update() [2/10]

void MVnormBetaFt::update ( const Grp ,
const Qgrp ,
const SigmaI ,
const Grp ,
const double &  ,
const SigmaI ,
const gsl_rng *   
)
virtual

Student- \(t\) likelihood, Student- \(t\) prior.

Parameters
[in]Grp&data for the likelihood
[in]Qgrp&vector Student- \(t\) covariance scale parameter for the likelihood covariance
[in]SigmaI&inverse-covariance matrix for the likelihood
[in]Grp&prior mean
[in]double&Student- \(t\) scale parameter for the prior covariance
[in]SigmaI&prior inverse-covariance matrix
[in]gsl_rng*pointer to a PNG

Reimplemented from MVnormBeta.

◆ update() [3/10]

void MVnormBetaFt::update ( const Grp ,
const Qgrp ,
const SigmaI ,
const Grp ,
const SigmaI ,
const gsl_rng *   
)
virtual

Student- \(t\) likelihood, Gaussian prior.

Parameters
[in]Grp&data for the likelihood
[in]Qgrp&vector Student- \(t\) covariance scale parameter for the likelihood covariance
[in]SigmaI&inverse-covariance matrix for the likelihood
[in]Grp&prior mean
[in]SigmaI&prior inverse-covariance matrix
[in]gsl_rng*pointer to a PNG

Reimplemented from MVnormBeta.

◆ update() [4/10]

void MVnormBetaFt::update ( const Grp ,
const Qgrp ,
const SigmaI ,
const gsl_rng *   
)
virtual

Sudent- \(t\) likelihood.

Parameters
[in]Grp&data
[in]Qgrp&vector of Student- \(t\) covariance scale parameters for the likelihood covariance
[in]SigmaI&inverse-covariance matrix for the likelihood
[in]gsl_rng*pointer to a PNG

Reimplemented from MVnormBeta.

◆ update() [5/10]

void MVnormBetaFt::update ( const Grp ,
const Qgrp ,
const SigmaI ,
const SigmaI ,
const gsl_rng *   
)
virtual

Student- \(t\) likelihood, Gaussian prior.

Parameters
[in]Grp&data for the likelihood
[in]Qgrp&vector of Student- \(t\) covariance scale parameters for the likelihood covariance
[in]SigmaI&inverse-covariance matrix for the likelihood
[in]SigmaI&prior inverse-covariance matrix
[in]gsl_rng*pointer to a PNG

Reimplemented from MVnormBeta.

◆ update() [6/10]

void MVnormBetaFt::update ( const Grp ,
const SigmaI ,
const double &  ,
const SigmaI ,
const gsl_rng *   
)
virtual

Gaussian likelihood, Student- \(t\) prior.

Parameters
[in]Grp&data for the likelihood
[in]SigmaI&inverse-covariance matrix for the likelihood
[in]double&Student- \(t\) scale parameter for the prior covariance
[in]SigmaI&prior inverse-covariance matrix
[in]gsl_rng*pointer to a PNG

Reimplemented from MVnormBeta.

◆ update() [7/10]

void MVnormBetaFt::update ( const Grp ,
const SigmaI ,
const Grp ,
const double &  ,
const SigmaI ,
const gsl_rng *   
)
virtual

Gaussian likelihood, Student- \(t\) prior.

Parameters
[in]Grp&data for the likelihood
[in]SigmaI&inverse-covariance matrix for the likelihood
[in]Grp&prior mean
[in]double&Student- \(t\) scale parameter for the prior covariance
[in]SigmaI&prior inverse-covariance matrix
[in]gsl_rng*pointer to a PNG

Reimplemented from MVnormBeta.

◆ update() [8/10]

void MVnormBetaFt::update ( const Grp ,
const SigmaI ,
const Grp ,
const SigmaI ,
const gsl_rng *   
)
virtual

Gaussian likelihood, Gaussian prior.

Parameters
[in]Grp&data for the likelihood
[in]SigmaI&inverse-covariance matrix for the likelihood
[in]Grp&prior mean
[in]SigmaI&prior inverse-covariance matrix
[in]gsl_rng*pointer to a PNG

Reimplemented from MVnormBeta.

◆ update() [9/10]

void MVnormBetaFt::update ( const Grp ,
const SigmaI ,
const gsl_rng *   
)
virtual

Gaussian likelihood.

Parameters
[in]Grp&data
[in]SigmaI&inverse-covariance matrix for the likelihood
[in]gsl_rng*pointer to a PNG

Reimplemented from MVnormBeta.

◆ update() [10/10]

void MVnormBetaFt::update ( const Grp ,
const SigmaI ,
const SigmaI ,
const gsl_rng *   
)
virtual

Gaussian likelihood, Gaussian prior.

Parameters
[in]Grp&data for the likelihood
[in]SigmaI&inverse-covariance matrix for the likelihood
[in]SigmaI&prior inverse-covariance matrix
[in]gsl_rng*pointer to a PNG

Reimplemented from MVnormBeta.

Member Data Documentation

◆ _fitted

gsl_matrix_view MVnormBetaFt::_fitted
protected

Matrix of already fitted values.

\( \boldsymbol{XB} \) matrix for all predictors but the current one (i.e., \( \boldsymbol{X}_{\cdot -j}\boldsymbol{B}_{-j\cdot} \)). It is subtracted from the response and the regression is performed on the residual.


The documentation for this class was generated from the following files: