fat2Lpoly {fat2Lpoly} | R Documentation |
Two-locus Family-based Association Test with Polytomous Outcome
Description
Performs family-based association tests with a polytomous outcome under 2-locus and 1-locus models as described in reference [1]. Various functions design.constraint
to create design matrices are provided in this package. When SNP pairs are specified, the tested SNP is the second one of each pair, while the first one is considered the conditioning SNP. The function may also perform one-locus tests if individual SNPs are specified instead of SNP pairs.
Usage
fat2Lpoly(pedfilenames, datfilenames, freq.data, ibdfilenames = NULL,
snp.names.mat, ibd.loci = NULL, joint.tests = NULL,
contingency.file = FALSE, design.constraint,
par.constrained, constraints, pairweights=calcule.poids.alphafixe,
lc = NULL, alpha = NULL)
Arguments
pedfilenames |
vector of 1 or 2 (the number of loci involved in the |
datfilenames |
vector of 1 or 2 (the number of loci involved in the |
freq.data |
Either
(1) a vector of 1 or 2 (the number of loci involved in the |
ibdfilenames |
vector of 1 or 2 (the number of loci involved in the |
snp.names.mat |
matrix of one or two columns giving the names of the SNPs (if one column matrix) or pairs of SNPs (if two columns matrix) to be analyzed. These SNPs represent all or part of the SNPs in the data files |
ibd.loci |
matrix of the same dimensions as |
joint.tests |
list of vectors of numbers between 1 and the total number of parameters in the |
contingency.file |
if 'TRUE' (default is 'FALSE'), then a file called descriptive_statistics'date_and_time'.txt is created and contingency tables with the numbers of subjects per level are progressively added to this file. |
design.constraint |
function building the design matrices WITHIN each category, for constraints specific to each category. It also returns the design matrices comprising only the loci main effects that are used for computing the covariances. An attribute |
par.constrained |
Optional matrix of dimensions ( |
constraints |
Optional matrix of dimensions ( |
pairweights |
function calculating the weights of the observation pair differences when conditioning on the first SNP in the test of the second SNP in a SNP pair. Default is calcule.poids.alphafixe, implementing the weighting function of equation (6) of reference [1]. An alternative is calcule.poids.Chen, implementing the weighting function of equation (7) of reference [1]. |
lc |
numerical identifier of the SNP (locus) on which to condition when testing model terms. Defaults to NULL, or no conditioning. |
alpha |
vector of length |
Details
All subjects included in the pedigree files must also be found in the IBD files.
All fields in the pedigree files must be numeric. No letters allowed, even for family and subject ID's.
Families whose genotyped subjects are all in the same category (phenotype combination), are uninformative and will be excluded.
Conditioning on the first SNP in a SNP pair is implemented by weighting the observation pair differences according to a model of the polytomous outcome as a function of the first SNP genotypes. The function converting the coefficients of this regression model into weights is specified by the argument pairweights
. The default function calcule.poids.alphafixe
provided satisfactory power in simulations described in reference [1].
File "descriptive_statistics'date_and_time'.txt" (will be created if contingency.file='TRUE'): For each tested SNP, it shows contingency tables of the subjects in the 2 or 4 different categories, first for all families together and then for each individual family.
If one or both of the arguments ibd.loci and ibdfilenames are left unspecified (or NULL, their default), then we use the kinship coefficients multiplied by two, instead of the expectation of the IBD probabilities, in the computation of the score statistics. The kinship coefficients are obtained using the function kinship
from the package kinship2
.
Value
returns a list of five objects:
scores.covs.all.SNPs |
list of length 'nrow( |
p.values.scores |
data frame of p-values for all the SNPs or SNP pairs in |
MA.table |
data frame giving the minor allele numbers of all the SNPs contained in the allele frequency files. |
y1 |
affection name extracted from first line of the data file(s) |
y2 |
affection name extracted from second line of the data file(s) |
Author(s)
Alexandre Bureau and Jordie Croteau
References
1. Bureau A., Croteau J., Chagnon, Y.C., Roy, M.-A. and Maziade, M. Extension of the Generalized Disequilibrium Test to polytomous phenotypes and two locus models. Frontiers in Genetics, 5: Article 258. 2. http://www.sph.umich.edu/csg/abecasis/Merlin/tour/input_files.html
See Also
Examples
path.data=paste(.libPaths()[which(unlist(lapply(.libPaths(),
function(x) length(grep("fat2Lpoly",dir(x)))))>0)],"/fat2Lpoly/extdata/",sep="")
if(length(path.data)>1) path.data=path.data[length(path.data)]
snps.anal=c("snp3.loc2","snp4.loc2")
microsat.names.loc2=c("2_3_mrk:","2_4_mrk:")
############ design.endo2disease with conditioning on locus 1 ################
## Not run:
joint.tests=list(c(2,5))
snp.names.mat=cbind(rep("snp4.loc1",length(snps.anal)),snps.anal)
microsat.names.mat=cbind(rep("1_4_mrk:",length(snps.anal)),microsat.names.loc2)
test=fat2Lpoly(pedfilenames=paste(path.data,c("loc1.ped","loc2.ped"),sep=""),
datfilenames=paste(path.data,c("loc1.dat","loc2.dat"),sep=""),
freq.data=paste(path.data,c("loc1.freq","loc2.freq"),sep=""),
ibdfilenames=paste(path.data,c("loc1.ibd","loc2.ibd"),sep=""),
snp.names.mat=snp.names.mat,ibd.loci=microsat.names.mat,
joint.tests=joint.tests,contingency.file=TRUE,
design.constraint=design.endo2disease,lc=1)
test$p.values.scores
## End(Not run)
###############################################################################
################### design.endo2disease without conditioning ##################
joint.tests=list(c(2,5))
snp.names.mat=cbind(rep("snp4.loc1",length(snps.anal)),snps.anal)
microsat.names.mat=cbind(rep("1_4_mrk:",length(snps.anal)),microsat.names.loc2)
test=fat2Lpoly(pedfilenames=paste(path.data,c("loc1.ped","loc2.ped"),sep=""),
datfilenames=paste(path.data,c("loc1.dat","loc2.dat"),sep=""),
freq.data=paste(path.data,c("loc1.freq","loc2.freq"),sep=""),
ibdfilenames=paste(path.data,c("loc1.ibd","loc2.ibd"),sep=""),
snp.names.mat=snp.names.mat,ibd.loci=microsat.names.mat,
joint.tests=joint.tests,contingency.file=FALSE,
design.constraint=design.endo2disease)
test$p.values.scores
###############################################################################
################# design.full with conditioning on locus 1 ##################
## Not run:
joint.tests=list(c(2,3),c(5,6),c(8,9),c(2,3,5,6,8,9))
snp.names.mat=cbind(rep("snp4.loc1",length(snps.anal)),snps.anal)
microsat.names.mat=cbind(rep("1_4_mrk:",length(snps.anal)),microsat.names.loc2)
test=fat2Lpoly(pedfilenames=paste(path.data,c("loc1.ped","loc2.ped"),sep=""),
datfilenames=paste(path.data,c("loc1.dat","loc2.dat"),sep=""),
freq.data=paste(path.data,c("loc1.freq","loc2.freq"),sep=""),
ibdfilenames=paste(path.data,c("loc1.ibd","loc2.ibd"),sep=""),
snp.names.mat=snp.names.mat,ibd.loci=microsat.names.mat,
joint.tests=joint.tests,
design.constraint=design.full,lc=1)
test$p.values.scores
## End(Not run)
##############################################################################
############################# design.1locus #################################
snp.names.mat=as.matrix(snps.anal)
microsat.names.mat=as.matrix(microsat.names.loc2)
test=fat2Lpoly(pedfilenames=paste(path.data,"loc2.ped",sep=""),
datfilenames=paste(path.data,"loc2.dat",sep=""),
freq.data=paste(path.data,"loc2.freq",sep=""),
ibdfilenames=paste(path.data,"loc2.ibd",sep=""),
snp.names.mat=snp.names.mat,ibd.loci=microsat.names.mat,
design.constraint=design.1locus)
test$p.values.scores
##############################################################################
############# design.dichotomous with conditioning on locus 1 ##############
## Not run:
joint.tests=list(c(2,3))
snp.names.mat=cbind(rep("snp4.loc1",length(snps.anal)),snps.anal)
microsat.names.mat=cbind(rep("1_4_mrk:",length(snps.anal)),microsat.names.loc2)
test=fat2Lpoly(pedfilenames=paste(path.data,c("loc1.ped","loc2.ped"),sep=""),
datfilenames=paste(path.data,c("loc1.dat","loc2.dat"),sep=""),
freq.data=paste(path.data,c("loc1.freq","loc2.freq"),sep=""),
ibdfilenames=paste(path.data,c("loc1.ibd","loc2.ibd"),sep=""),
snp.names.mat=snp.names.mat,ibd.loci=microsat.names.mat,
joint.tests=joint.tests,
design.constraint=design.dichotomous,lc=1)
test$p.values.scores
## End(Not run)
##############################################################################