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

Individual vector of conditional regression coefficients with blocks of traits. More...

#include <MuGen.h>

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

Public Member Functions

 MVnormBetaFtBlk ()
 Default constructor. More...
 
 MVnormBetaFtBlk (const Grp &resp, gsl_matrix *pred, vector< double > &eaFt, const vector< size_t > &iCl, const gsl_matrix *Sig, const vector< size_t > &blkStart, const gsl_rng *r, gsl_matrix *bet, const size_t &iRw)
 Constructor with no prior index. More...
 
 MVnormBetaFtBlk (const Grp &resp, gsl_matrix *pred, vector< double > &eaFt, const vector< size_t > &iCl, const gsl_matrix *Sig, const vector< size_t > &blkStart, const gsl_rng *r, const size_t &up, gsl_matrix *bet, const size_t &iRw)
 Constructor with a prior index. More...
 
 MVnormBetaFtBlk (const gsl_matrix *resp, gsl_matrix *pred, vector< double > &eaFt, const vector< size_t > &iCl, const gsl_matrix *Sig, const vector< size_t > &blkStart, const gsl_rng *r, gsl_matrix *bet, const size_t &iRw)
 Constructor with no prior index. More...
 
 MVnormBetaFtBlk (const gsl_matrix *resp, gsl_matrix *pred, vector< double > &eaFt, const vector< size_t > &iCl, const gsl_matrix *Sig, const vector< size_t > &blkStart, const gsl_rng *r, const size_t &up, gsl_matrix *bet, const size_t &iRw)
 Constructor with a prior index. More...
 
 ~MVnormBetaFtBlk ()
 Destructor.
 
void update (const Grp &dat, const SigmaI &SigIb, const gsl_rng *r)
 Gaussian likelihood. More...
 
void update (const Grp &dat, const Qgrp &q, const SigmaI &SigIb, const gsl_rng *r)
 Sudent- \(t\) likelihood. More...
 
void update (const Grp &dat, const SigmaI &SigIb, const SigmaI &SigIp, const gsl_rng *r)
 Gaussian likelihood, Gaussian prior. More...
 
void update (const Grp &dat, const SigmaI &SigIb, const double &qPr, const SigmaI &SigIp, const gsl_rng *r)
 Gaussian likelihood, Student- \(t\) prior. More...
 
void update (const Grp &dat, const Qgrp &q, const SigmaI &SigIb, const SigmaI &SigIp, const gsl_rng *r)
 Student- \(t\) likelihood, Gaussian prior. More...
 
void update (const Grp &dat, 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 &dat, const SigmaI &SigIb, const Grp &muPr, const SigmaI &SigIp, const gsl_rng *r)
 Gaussian likelihood, Gaussian prior. More...
 
void update (const Grp &dat, 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 &dat, 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 &dat, 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 MVnormBetaBlk
 MVnormBetaBlk ()
 Default constructor. More...
 
 MVnormBetaBlk (const Grp &resp, gsl_matrix *pred, const vector< size_t > &iCl, const gsl_matrix *Sig, const vector< size_t > &blkStart, const gsl_rng *r, gsl_matrix *bet, const size_t &iRw)
 Constructor with no prior index. More...
 
 MVnormBetaBlk (const Grp &resp, gsl_matrix *pred, const vector< size_t > &iCl, const gsl_matrix *Sig, const vector< size_t > &blkStart, const gsl_rng *r, const size_t &up, gsl_matrix *bet, const size_t &iRw)
 Constructor with a prior index. More...
 
 MVnormBetaBlk (const gsl_matrix *resp, gsl_matrix *pred, const vector< size_t > &iCl, const gsl_matrix *Sig, const vector< size_t > &blkStart, const gsl_rng *r, gsl_matrix *bet, const size_t &iRw)
 Constructor with no prior index. More...
 
 MVnormBetaBlk (const gsl_matrix *resp, gsl_matrix *pred, const vector< size_t > &iCl, const gsl_matrix *Sig, const vector< size_t > &blkStart, const gsl_rng *r, const size_t &up, gsl_matrix *bet, const size_t &iRw)
 Constructor with a prior index. More...
 
virtual ~MVnormBetaBlk ()
 Destructor.
 
const size_t * up () const
 Points to the prior. 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...
 
virtual double scalePar () const
 Scale parameter. More...
 

Protected Attributes

gsl_matrix_view _fitted
 Matrix of fitted values. More...
 
- Protected Attributes inherited from MVnormBetaBlk
const vector< size_t > * _blkStart
 Start positions of blocks. More...
 
vector< gsl_vector_view > _eachVec
 Trait blocks. More...
 
vector< gsl_vector_view > _X
 Separate predictors. More...
 
vector< double > _scale
 Vector of scales. More...
 
const size_t * _upLevel
 Pointer to a 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

Individual vector of conditional regression coefficients with blocks of traits.

Implements separate regression models for blocks of traits. The likelihood covariance matrix is block-diagonal. Traits of the same block have to be contiguous within the vector. The regression is conditional on other regression coefficient estimates within a group (i.e., one-at-a-time updating).

Constructor & Destructor Documentation

◆ MVnormBetaFtBlk() [1/5]

MVnormBetaFtBlk::MVnormBetaFtBlk ( )
inline

Default constructor.

Simply calls the default constructor of the MVnormBetaBlk class. Pointer to the fitted matrix is not set to anything.

◆ MVnormBetaFtBlk() [2/5]

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

Constructor with no prior index.

Stochastic constructor for a regression with an improper prior.

Parameters
[in]Grp&response variable
[in]gsl_matrix*predictor matrix
[in]vector<double>&vectorized fitted matrix
[in]vector<size_t>&column indexes of the predictor matrix
[in]gsl_matrix*covariance marix for stochastic initiation
[in]vector<size_t>&vector of start indixes of each block
[in]gsl_rng*pointer to a PNG
[in]gsl_matrix*matrix of regression coefficients
[in]size_t&row index of the regression coefficient matrix

◆ MVnormBetaFtBlk() [3/5]

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

Constructor with a prior index.

Stochastic constructor for a regression with a proper prior.

Parameters
[in]Grp&response variable
[in]gsl_matrix*predictor matrix
[in]vector<double>&vectorized fitted matrix
[in]vector<size_t>&column indexes of the predictor matrix
[in]gsl_matrix*covariance marix for stochastic initiation
[in]vector<size_t>&vector of start indixes of each block
[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

◆ MVnormBetaFtBlk() [4/5]

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

Constructor with no prior index.

Stochastic constructor for a regression with an improper prior.

Parameters
[in]gsl_matrix*response matrix
[in]gsl_matrix*predictor matrix
[in]vector<double>&vectorized fitted matrix
[in]vector<size_t>&column indexes of the predictor matrix
[in]gsl_matrix*covariance marix for stochastic initiation
[in]vector<size_t>&vector of start indixes of each block
[in]gsl_rng*pointer to a PNG
[in]gsl_matrix*matrix of regression coefficients
[in]size_t&row index of the regression coefficient matrix

◆ MVnormBetaFtBlk() [5/5]

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

Constructor with a prior index.

Stochastic constructor for a regression with a proper prior.

Parameters
[in]gsl_matrix*response matrix
[in]gsl_matrix*predictor matrix
[in]vector<double>&vectorized fitted matrix
[in]vector<size_t>&column indexes of the predictor matrix
[in]gsl_matrix*covariance marix for stochastic initiation
[in]vector<size_t>&vector of start indixes of each block
[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

Member Function Documentation

◆ update() [1/10]

void MVnormBetaFtBlk::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 MVnormBetaBlk.

◆ update() [2/10]

void MVnormBetaFtBlk::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 MVnormBetaBlk.

◆ update() [3/10]

void MVnormBetaFtBlk::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 MVnormBetaBlk.

◆ update() [4/10]

void MVnormBetaFtBlk::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 MVnormBetaBlk.

◆ update() [5/10]

void MVnormBetaFtBlk::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 MVnormBetaBlk.

◆ update() [6/10]

void MVnormBetaFtBlk::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 MVnormBetaBlk.

◆ update() [7/10]

void MVnormBetaFtBlk::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 MVnormBetaBlk.

◆ update() [8/10]

void MVnormBetaFtBlk::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 MVnormBetaBlk.

◆ update() [9/10]

void MVnormBetaFtBlk::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 MVnormBetaBlk.

◆ update() [10/10]

void MVnormBetaFtBlk::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 MVnormBetaBlk.

Member Data Documentation

◆ _fitted

gsl_matrix_view MVnormBetaFtBlk::_fitted
protected

Matrix of 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: