Library
LiteQTL.calculate_nrLiteQTL.calculate_nrLiteQTL.calculate_pxLiteQTL.cpurunLiteQTL.filter_mafLiteQTL.find_max_idx_valueLiteQTL.get_geno_dataLiteQTL.get_pheno_block_sizeLiteQTL.get_pheno_dataLiteQTL.get_standardized_matrixLiteQTL.get_standardized_matrix_gpuLiteQTL.gpu_square_lodLiteQTL.gpurunLiteQTL.is_corr_in_rangeLiteQTL.lod2pLiteQTL.lod_kernelLiteQTL.lod_score_multithreadLiteQTL.lod_score_multithreadLiteQTL.reduce_kernelLiteQTL.scanLiteQTL.set_blas_threads
LiteQTL.get_geno_data — Methodget_geno_data(file, datatype)
returns the genotype data. Will skip every other column because genotype probability is duplicated.
LiteQTL.get_pheno_data — Methodget_pheno_data(file, datatype; transposed)
returns the phenotype data. If transposed=true, then the data will be transposed.
LiteQTL.scan — Functionscan(Y, G)
scan(Y, G, X; maf_threshold, export_matrix, usegpu, lod_or_pval, timing_file)
This function is the main API for eQTL scans. EQTL scan process includes If no covariate: - filter maf if the threshold is greater than 0 - standardizing phenotype matrix (Y) and genotype matrix (G) - calculate correlation (R) matrix. - computes log of odds (LOD) score matrix, or p-value if lod_or_pval is set to lod - calculate maximum LOD score if export_matrix is set to false.
If covariates exists:
- calculate px
- computes Y hat and G hat by matrix multiplication.
- substract Y hat and G hat from Y and G respectively.
- filter maf if the threshold is greater than 0
- standardizing phenotype matrix (Y) and genotype matrix (G)
- calculate correlation (R) matrix.
- computes log of odds (LOD) score matrix, or p-value if `lod_or_pval` is set to `lod`
- calculate maximum LOD score if `export_matrix` is set to false.Arguments:
Y: a matrix of phenotypes.G: a matrix of genotypes.X: a matrix of covariates. Default isnothing. Ifnothing, scan is run without covariates.maf_threshold: a floating point number to indicate the maf_threshold. Default is 0.05. Set to 0 if no maf filtering should be done.export_matrix: a boolean value that determines whether the result should be the maximum value of LOD score of each phenotype and its corresponding index, or the whole LOD score matrix.usegpu: a boolean value that indicates whether to run scan function on GPU or CPU. Default is false, which runs scan on CPU.lod_or_pval: a string value of eitherlodorpvalto indicate the desired output.timing_file: a string that indicates the file location for the timing outputs. Default is nothing.
Output:
returns LOD score or pval, in vector or matrix format depending on value of export_matrix.
LiteQTL.calculate_nr — Methodcalculate_nr(a, b)
Computes correlation matrix (R) * n. n is not removed because of performance choice.
Arguments
a: standardized phenotype matrixb: standardized genotype matrix
Output:
returns correlation matrix R
LiteQTL.calculate_nr — Methodcalculate_nr(a, b)
Computes correlation matrix (R) * n on GPU with cuBLAS library. n is not removed because of performance choice.
Arguments
a: standardized phenotype matrixb: standardized genotype matrix
Output:
returns correlation matrix R on GPU device (of type CuArray).
LiteQTL.calculate_px — Methodcalculate_px(x)
Arguments
x: input matrix
Output:
returns px
LiteQTL.cpurun — Functioncpurun(pheno, geno)
cpurun(pheno, geno, X; maf_threshold, export_matrix, lod_or_pval, timing_file)
Arguments:
Y: a matrix of phenotypes.G: a matrix of genotypes.X: a matrix of covariates. Default isnothing. Ifnothing, scan is run without covariates.maf_threshold: a floating point number to indicate the maf_threshold. Default is 0.05. Set to 0 if no maf filtering should be done.export_matrix: a boolean value that determines whether the result should be the maximum value of LOD score of each phenotype and its corresponding index, or the whole LOD score matrix.lod_or_pval: a string value of eitherlodorpvalto indicate the desired output.timing_file: a string that indicates the file location for the timing outputs. Default is nothing.
Output:
returns LOD score or pval, in vector or matrix format depending on value of export_matrix.
LiteQTL.filter_maf — Methodfilter_maf(genotype; maf_threshold)
Filter genotype data with a minor allele frequency threshold.
Arguments
genotype: genotype matrixmaf_threshold: default value is 0.05.
Output:
returns filtered genotype matrix.
LiteQTL.find_max_idx_value — Methodfind_max_idx_value(lod)
Computes the index of maximum, and maximum value of each row of a matrix. Optimized with multi-threading.
!Notes: Set the thread number with env JULIANUMTHREADS to your desired number of threads. For example: JULIA_NUM_THREADS=16 julia
Arguments
lod: input matrix.
Output:
returns a matrix with two columns, first column is the index of maximum, second column is the maximum value.
LiteQTL.get_pheno_block_size — Methodget_pheno_block_size(n, m, p, datatype)
Computes the number of blocks for chuncking large phenotype matrix, to fit into GPU memory.
Arguments
n: number of individualsm: number of phenotypesp: number of genotype markers- datatype: datatype of phenotype. (Float32, or Float64)
Output:
returns
num_block: number of blocks to disect phenotypeblock_size: size of each block.
LiteQTL.get_standardized_matrix — Methodget_standardized_matrix(mat)
Computes standardized matrix on CPU utilizing multi-threads.
Arguments
m: a matrix to be standardized
Output:
returns the standardized matrix of m
LiteQTL.get_standardized_matrix_gpu — Methodget_standardized_matrix_gpu(m)
Computes standardized matrix on GPU. Theoretically this function works with CPU as well. But it is not as fast as the get_standardized_matrix function, which utilizes multi threads on CPU.
Arguments
m: a matrix to be standardized
Output:
returns the standardized matrix of m
LiteQTL.gpu_square_lod — Methodgpu_square_lod(d_nr, n, m, p, export_matrix)
Sets up GPU environment and calls GPU kernels.
Arguments
d_nr: Output matrix (correlation matrix) in CuArray formatn: number of individualsm: number of phenotypesp: number of genotype markersexport_matrix: boolean value to indicate user desire output
Output:
returns correlation matrix R, or maximum of R if export_matrix is set to false
LiteQTL.gpurun — Functiongpurun(pheno, geno)
gpurun(pheno, geno, X; maf_threshold, export_matrix, timing_file)
Arguments:
Y: a matrix of phenotypes.G: a matrix of genotypes.X: a matrix of covariates. Default isnothing. Ifnothing, scan is run without covariates.maf_threshold: a floating point number to indicate the maf_threshold. Default is 0.05. Set to 0 if no maf filtering should be done.export_matrix: a boolean value that determines whether the result should be the maximum value of LOD score of each phenotype and its corresponding index, or the whole LOD score matrix.lod_or_pval: a string value of eitherlodorpvalto indicate the desired output.timing_file: a string that indicates the file location for the timing outputs. Default is nothing.
Output:
returns LOD score, in vector or matrix format depending on value of export_matrix. P value for GPU is not currently implementated.
LiteQTL.is_corr_in_range — Methodis_corr_in_range(r, min, max)
Checks whether correlation matrix is in range. Only runs if DEBUG flag is turned on from julia commandline.
Arguments
r: a matrix to be standardized- min : minimum value of range
- max : maximum value of range
Output:
errors message will show if correlation matrix is not in range.
LiteQTL.lod2p — Methodlod2p(lod)
Computes p-value based on log of odds score.
Arguments
lod: LOD matrix
Output:
returns p-value.
LiteQTL.lod_kernel — Methodlod_kernel(nr, MAX, n)
GPU kernel that computes LOD scores.
Arguments
nr: correlation matrix. will be updated with output results. (Type CuArray)MAX: number of total GPU threads, also called ndrange.n: number of individuals
Output:
Results written back to input matrix nr. Returns the LOD scores.
LiteQTL.lod_score_multithread — Methodlod_score_multithread(m, nr)
Computes log of odds (LOD) score. Optimized for correlation matrix type is Float32 (single precision).
!Notes: Set the thread number with env JULIANUMTHREADS to your desired number of threads. For example: JULIA_NUM_THREADS=16 julia
Arguments
m: number of individuals.nr: correlation matrix R times n. N will be removed during this step.
Output:
returns LOD score.
LiteQTL.lod_score_multithread — Methodlod_score_multithread(m, nr)
Computes log of odds (LOD) score. Optimized for correlation matrix type is Float64 (double precision).
!Notes: Set the thread number with env JULIANUMTHREADS to your desired number of threads. For example: JULIA_NUM_THREADS=16 julia
Arguments
m: number of individuals.nr: correlation matrix R times n. N will be removed during this step.
Output:
returns LOD score.
LiteQTL.reduce_kernel — Methodreduce_kernel(input, rows, cols)
Computes the maximum LOD score of each row. Results written back to the first two columns.
Arguments
input: LOD scorerows: number of rows forinputcols: number of cols forinput
Output:
Output is written to the input matrix, first column of input matrix is the index of maximum, and second column is the maximum value of that row.
LiteQTL.set_blas_threads — Methodset_blas_threads(nthread)
Sets the number of threads that will be used by BLAS library. !Notes: Set the thread number with env JULIANUMTHREADS to your desired number of threads. For example: JULIA_NUM_THREADS=16 julia
Arguments
nthread: desired number of threads
Output:
Shows number of threads set.