Public

Public functions

Documentation for JWAS.jl's public interface. Below are functions available to general users.

Index

Public Interface

JWAS.build_modelFunction.
build_model(model_equations::AbstractString,R;df::Float64=4.0)
  • Build models from model equations with residual varainces R and degree of freedom for residual variance df defaulting to 4.0.

  • By default, all variabels in model_equations are fixed and factors. Set variables to be covariates or random using functions set_covariate() or set_random().

#single-trait
model_equations = "BW = intercept + age + sex"
R               = 6.72
models          = build_model(model_equations,R);

#multi-trait
model_equations = "BW = intercept + age + sex
                   CW = intercept + litter";
R               = [6.72   24.84
                   24.84  708.41]
models          = build_model(model_equations,R);
source
JWAS.set_covariateFunction.
set_covariate(mme::MME,variables::AbstractString...)
  • set variables as covariates; mme is the output of function build_model().

#After running build_model, variabels age and year can be set to be covariates as
set_covariate(models,"age","year")
#or
set_covariate(models,"age year")
source
JWAS.set_randomFunction.
set_random(mme::MME,randomStr::AbstractString,ped::Pedigree, G;df=4)
  • set variables as random polygenic effects with pedigree information ped, variances G whose degree of freedom df defaults to 4.0.

#single-trait (example 1)
model_equation  = "y = intercept + Age + Animal"
model           = build_model(model_equation,R)
ped             = get_pedigree(pedfile)
G               = 1.6
set_random(model,"Animal Animal*Age", ped,G)

#single-trait (example 2)
model_equation  = "y = intercept + Age + Animal + Animal*Age"
model           = build_model(model_equation,R)
ped             = get_pedigree(pedfile)
G               = [1.6   0.2
                   0.2  1.0]
set_random(model,"Animal Animal*Age", ped,G)

#multi-trait
model_equations = "BW = intercept + age + sex + Animal
                   CW = intercept + age + sex + Animal"
model           = build_model(model_equations,R);
ped             = get_pedigree(pedfile);
G               = [6.72   2.84
                   2.84  8.41]
set_random(model,"Animal", ped,G)
source
set_random(mme::MME,randomStr::AbstractString,G;df=4)
  • set variables as i.i.d random effects with variances G whose degree of freedom df defaults to 4.0.

#single-trait (example 1)
model_equation  = "y = intercept + litter + sex"
model           = build_model(model_equation,R)
G               = 0.6
set_random(model,"litter",G)

#multi-trait
model_equations = "BW = intercept + litter + sex
                   CW = intercept + litter + sex"
model           = build_model(model_equations,R);
G               = [3.72  1.84
                   1.84  3.41]
set_random(model,"litter",G)
source
JWAS.get_pedigreeFunction.
get_pedigree(pedfile::AbstractString;header=false,separator=' ')
  • Get pedigree informtion from a pedigree file.

  • File format:

a 0 0
b 0 0
c a b
d a c
source
JWAS.add_genotypesFunction.
add_genotypes(mme::MME,file,G;separator=' ',header=true,center=true,G_is_marker_variance=false,df=4.0)
  • Get marker informtion from a genotype file (same order as the phenotype file).

  • G defaults to the genetic variance with degree of freedom df=4.0.

  • File format:

Animal,marker1,marker2,marker3,marker4,marker5
S1,1,0,1,1,1
D1,2,0,2,2,1
O1,1,2,0,1,0
O3,0,0,2,1,1
source
JWAS.add_markersFunction.
same to add_genotypes
source
JWAS.runMCMCFunction.
runMCMC(mme,df;Pi=0.0,estimatePi=false,chain_length=1000,burnin = 0,starting_value=false,printout_frequency=100,missing_phenotypes=false,constraint=false,methods="conventional (no markers)",output_samples_frequency::Int64 = 0)

Run MCMC (marker information included or not) with sampling of variance components.

  • available methods include "conventional (no markers)", "RR-BLUP", "BayesB", "BayesC".

  • save MCMC samples every output_samples_frequency iterations to files output_file defaulting to MCMC_samples.

  • the first burnin iterations are discarded at the beginning of an MCMC run

  • Pi for single-trait analyses is a number; Pi for multi-trait analyses is a dictionary such as Pi=Dict([1.0; 1.0]=>0.7,[1.0; 0.0]=>0.1,[0.0; 1.0]=>0.1,[0.0; 0.0]=>0.1),

    • if Pi (Π) is not provided in multi-trait analysis, it will be generated assuming all markers have effects on all traits.

  • starting_value can be provided as a vector for all location parameteres except marker effects.

  • print out the monte carlo mean in REPL with printout_frequency

  • constraint=true if constrain residual covariances between traits to be zeros.

source
outputMCMCsamples(mme::MME,trmStr::AbstractString...)

Get MCMC samples for specific location parameters.

source
JWAS.showMMEFunction.
showMME(mme::MME,df::DataFrame)
  • Show left-hand side and right-hand side of mixed model equations (no markers).

source
JWAS.solveFunction.
solve(mme::MME,df::DataFrame;solver="default",printout_frequency=100,tolerance = 0.000001,maxiter = 5000)
  • Solve the mixed model equations (no marker information) without estimating variance components.

Available solvers includes default,Jacobi,GaussSeidel,Gibbs sampler.

source
JWAS.misc.GWASFunction.
GWAS(marker_effects_file;header=false)

Compute the model frequency for each marker (the probability the marker is included in the model) using samples of marker effects stored in marker_effects_file.

source
GWAS(marker_effects_file,map_file,model;header=false,window_size="1 Mb",threshold=0.001)

run genomic window-based GWAS

  • MCMC samples of marker effects are stored in marker_effects_file

  • map_file has the marker position information

source
get_additive_genetic_variances(model::MME,files...;header=true)
  • Get MCMC samples for additive genetic variances using samples for marker effects stored in files.

  • Return a vector for single-trait analysis and an array of matrices for multi-trait analysis

source
get_heritability(samples_for_genetic_variances::Array{Array{Float64,2},1},samples_for_residual_variances::Array{Array{Float64,2},1}))
  • Get MCMC samples for heritabilities using MCMC samples for genetic variances and residual variances for multi-trait analysis

source
get_heritability(samples_for_genetic_variances::Array{Float64,1},samples_for_residual_variances::Array{Float64,1}))
  • Get MCMC samples for heritabilities using MCMC samples for genetic variances and residual variances for single-trait analysis

source
get_correlations(samples_for_genetic_variances::Array{Array{Float64,2},1})
  • Get MCMC samples for correlations using MCMC samples for covariance matrces

source
get_breeding_values(model::MME,files...;header=true)
  • Get esitimated breeding values and prediction error variances using samples of marker effects stored in files.

source
JWAS.misc.reformatFunction.
reformat(G::Array{Float64,2},ntraits=false)
  • convert the format from a Matrix to a Array of Matrices, when the matrix is

read from text file storing samples for covariance matrices or heritabilities

  • the size of the Matrix is

    • ntraits x nsamples for MCMC samples for heritbilities

    • (ntraits^2) x nsamples for MCMC samples for covariance matrices

  • the Array of Matrices is an array (length=nsamples) of

    • matrices of size ntraits X ntraits for covariance matrices

    • vectors of length ntraits for heritabilities

source
reformat(G::Array{Array{Float64,2},1},1})
  • convert the format from a Array of Matrices to a Matrix, then the Matrix chance

be write to text files to save samples for covariance matrices or heritabilities

  • the Array of Matrices is an array (length=nsamples) of

    • matrices of size ntraits X ntraits for covariance matrices (Array{Array{Float64,2},1})

    • vectors of length ntraits for heritabilities (Array{Array{Float64,1},1})

  • the size of the Matrix is

    • ntraits x nsamples for MCMC samples for heritbilities

    • (ntraits^2) x nsamples for MCMC samples for covariance matrices

source
JWAS.misc.reportFunction.
report(X::Array{Array{Float64,2},1};index=false)
  • show summary statistics for MCMC samples (matrices or index [i,j] of matrices)

source
report(X::Array{Array{Float64,1},1};index=false)
  • show summary statistics for MCMC samples (vectors or index i of vectors)

source
JWAS.misc.QCFunction.
QC(infile,outfile;separator=' ',header=true,missing=false,MAF=0.1)
  • Quality control for input file infile, then write out to output file: outfile.

    • Delete loci with minor allele frequency < MAF.

    • missing genotypes are replaced by column means.

  • File format (header=true,separator=',',missing=9):

Animal,marker1,marker2,marker3,marker4,marker5
S1,1,0,1,1,1
D1,2,0,9,2,1
O1,1,2,0,1,0
O3,0,0,2,1,1
source