MuGen
Multitrait genetics
|
C++ classes for hierarchical Bayesian Multi-trait quantitative-genetic models. More...
#include <gsl/gsl_matrix.h>
#include <gsl/gsl_vector.h>
#include <gsl/gsl_rng.h>
#include <gsl/gsl_permutation.h>
#include <vector>
#include <string>
#include <list>
Go to the source code of this file.
Classes | |
class | MVnorm |
The abstract base class for location parameter rows. More... | |
class | MVnormMu |
Individual vector of means. More... | |
class | MVnormMuPEX |
Individual vector of means with parameter expansion. More... | |
class | MVnormBetaPEX |
Individual regression with parameter expansion. More... | |
class | MVnormMuMiss |
Individual vector of means with missing data. More... | |
class | MVnormBeta |
Generic regression. More... | |
class | MVnormBetaFt |
Regression with multiple predictors. More... | |
class | MVnormMuBlk |
Individual vector of means with blocks of traits. More... | |
class | MVnormBetaBlk |
Individual vector of regression coefficients with blocks of traits. More... | |
class | MVnormBetaFtBlk |
Individual vector of conditional regression coefficients with blocks of traits. More... | |
class | RanIndex |
Generic index class. More... | |
class | RanIndexVS |
BVSR index class. More... | |
class | Apex |
Multiplicative redundant parameter. More... | |
class | Grp |
Base location parameter group class. More... | |
class | MuGrp |
Hierarchical mean. More... | |
class | MuGrpPEX |
"Random effect" with parameter expansion More... | |
class | MuGrpMiss |
Data with missing phenotype values. More... | |
class | MuGrpEE |
Data with measurement error. More... | |
class | MuGrpEEmiss |
Data with measurement error and missing phenotypes. More... | |
class | BetaGrpFt |
Multivariate multiple regression. More... | |
class | BetaGrpPEX |
Multivariate multiple regression with parameter expansion. More... | |
class | BetaGrpPC |
Relationship matrix regression. More... | |
class | BetaGrpPCpex |
Multiplicative parameter expansion for PC regression. More... | |
class | BetaGrpSnp |
Simple single-SNP regression class. More... | |
class | BetaGrpSnpCV |
Single-SNP regression with conditional variance. More... | |
class | BetaGrpPSR |
Single-SNP regression with partial effects. More... | |
class | BetaGrpSnpMiss |
Simple single-SNP regression class with missing data. More... | |
class | BetaGrpSnpMissCV |
Single-SNP regression with conditional variance and missing data. More... | |
class | BetaGrpPSRmiss |
Single-SNP regression with partial effects and missing data. More... | |
class | BetaGrpBVSR |
Bayesian variable selection regression. More... | |
class | MuBlk |
Hierarchical mean with independent blocks of traits. More... | |
class | BetaBlk |
Regression with independent blocks of traits. More... | |
class | SigmaI |
Basic inverse-covariance. More... | |
class | SigmaIblk |
Block-diagonal inverse-covariance. More... | |
class | SigmaIpex |
PEX inverse-covariance. More... | |
class | Qgrp |
Standard Student- \(t\) weights. More... | |
class | QgrpPEX |
Student- \(t\) weights for PEX. More... | |
class | MixP |
Dirichlet-multinomial mixture prior. More... | |
Functions | |
void | MVgauss (const gsl_vector *, const gsl_matrix *, const gsl_rng *, gsl_vector *) |
Multivariate Gaussian sampling function. More... | |
void | Wishart (const gsl_matrix *, const int &, const gsl_rng *, gsl_matrix *) |
void | Wishart (const gsl_matrix *, const size_t &, const gsl_rng *, gsl_matrix *) |
size_t | rtgeom (const double &, const size_t &, const gsl_rng *) |
Truncated geometric distribution. More... | |
double | mhl (const gsl_vector *beta, const gsl_matrix *SigI) |
Mahalanobis distance. More... | |
void | colCenter (gsl_matrix *inplace) |
Matrix centering in-place. More... | |
void | colCenter (const gsl_matrix *source, gsl_matrix *res) |
Matrix centering with copy. More... | |
void | colCenter (gsl_matrix *inplace, const double &absLab) |
Matrix centering in-place with missing values. More... | |
void | colCenter (const gsl_matrix *source, gsl_matrix *res, const double &absLab) |
Matrix centering with copy and missing values. More... | |
void | vecCenter (gsl_vector *inplace) |
Vector centering in-place. More... | |
void | vecCenter (const gsl_vector *source, gsl_vector *res) |
Vector centering with copy. More... | |
void | meanImpute (gsl_matrix *inplace, const double &absLab) |
Mean imputation without centering. More... | |
void | printMat (const gsl_matrix *m) |
Print matrix to screen. More... | |
unsigned long long | rdtsc () |
Accessing the processor RTDSC instruction. More... | |
MuGrp | operator+ (const Grp &m1, const Grp &m2) |
Addition operator. More... | |
MuGrp | operator+ (const MuGrp &m1, const Grp &m2) |
Addition operator. More... | |
MuGrp | operator+ (const Grp &m1, const MuGrp &m2) |
Addition operator. More... | |
MuGrp | operator- (const Grp &m1, const Grp &m2) |
Subtraction operator. More... | |
MuGrp | operator- (const MuGrp &m1, const Grp &m2) |
Subtraction operator. More... | |
MuGrp | operator- (const Grp &m1, const MuGrp &m2) |
Subtraction operator. More... | |
C++ classes for hierarchical Bayesian Multi-trait quantitative-genetic models.
This is the project header file containing class definitions and interface documentation.