diff options
Diffstat (limited to 'src')
-rw-r--r-- | src/lm.cpp | 8 | ||||
-rw-r--r-- | src/lmm.cpp | 202 | ||||
-rw-r--r-- | src/param.h | 5 |
3 files changed, 84 insertions, 131 deletions
@@ -362,7 +362,7 @@ void LM::AnalyzeGene(const gsl_matrix *W, const gsl_vector *x) { time_opt += (clock() - time_start) / (double(CLOCKS_PER_SEC) * 60.0); // Store summary data. - SUMSTAT SNPs = {beta, se, 0.0, 0.0, p_wald, p_lrt, p_score}; + SUMSTAT SNPs = {beta, se, 0.0, 0.0, p_wald, p_lrt, p_score, -0.0 }; sumStat.push_back(SNPs); } cout << endl; @@ -587,7 +587,7 @@ void LM::Analyzebgen(const gsl_matrix *W, const gsl_vector *y) { time_opt += (clock() - time_start) / (double(CLOCKS_PER_SEC) * 60.0); // Store summary data. - SUMSTAT SNPs = {beta, se, 0.0, 0.0, p_wald, p_lrt, p_score}; + SUMSTAT SNPs = {beta, se, 0.0, 0.0, p_wald, p_lrt, p_score, -0.0}; sumStat.push_back(SNPs); } cout << endl; @@ -702,7 +702,7 @@ void LM::AnalyzeBimbam(const gsl_matrix *W, const gsl_vector *y) { time_opt += (clock() - time_start) / (double(CLOCKS_PER_SEC) * 60.0); // Store summary data. - SUMSTAT SNPs = {beta, se, 0.0, 0.0, p_wald, p_lrt, p_score}; + SUMSTAT SNPs = {beta, se, 0.0, 0.0, p_wald, p_lrt, p_score, -0.0}; sumStat.push_back(SNPs); } cout << endl; @@ -844,7 +844,7 @@ void LM::AnalyzePlink(const gsl_matrix *W, const gsl_vector *y) { p_lrt, p_score); // store summary data - SUMSTAT SNPs = {beta, se, 0.0, 0.0, p_wald, p_lrt, p_score}; + SUMSTAT SNPs = {beta, se, 0.0, 0.0, p_wald, p_lrt, p_score, -0.0}; sumStat.push_back(SNPs); time_opt += (clock() - time_start) / (double(CLOCKS_PER_SEC) * 60.0); diff --git a/src/lmm.cpp b/src/lmm.cpp index 37f2f5b..e2f23a2 100644 --- a/src/lmm.cpp +++ b/src/lmm.cpp @@ -95,6 +95,7 @@ void LMM::CopyToParam(PARAM &cPar) { } void LMM::WriteFiles() { + string file_str; file_str = path_out + "/" + file_out; file_str += ".assoc.txt"; @@ -105,150 +106,99 @@ void LMM::WriteFiles() { return; } - if (!file_gene.empty()) { - outfile << "geneID" - << "\t"; - - if (a_mode == 1) { - outfile << "beta" - << "\t" - << "se" - << "\t" - << "l_remle" - << "\t" + auto common_header = [&] () { + outfile << "beta" << "\t" + << "se" << "\t"; + + outfile << "logl_H1" << "\t"; // we may make this an option + + switch(a_mode) { + case 1: + outfile << "l_remle" << "\t" << "p_wald" << endl; - } else if (a_mode == 2) { - outfile << "l_mle" - << "\t" + break; + case 2: + outfile << "l_mle" << "\t" << "p_lrt" << endl; - } else if (a_mode == 3) { - outfile << "beta" - << "\t" - << "se" - << "\t" + break; + case 3: + outfile << "p_score" << endl; + break; + case 4: + outfile << "l_remle" << "\t" + << "l_mle" << "\t" + << "p_wald" << "\t" + << "p_lrt" << "\t" << "p_score" << endl; - } else if (a_mode == 4) { - outfile << "beta" - << "\t" - << "se" - << "\t" - << "l_remle" - << "\t" - << "l_mle" - << "\t" - << "p_wald" - << "\t" - << "p_lrt" - << "\t" - << "p_score" << endl; - } else { + break; } + }; + + auto sumstats = [&] (SUMSTAT st) { + outfile << scientific << setprecision(6) << st.beta << "\t" + << st.se << "\t"; + + outfile << st.logl_H1 << "\t"; + + switch(a_mode) { + case 1: + outfile << st.lambda_remle << "\t" + << st.p_wald << endl; + break; + case 2: + outfile << st.lambda_mle << "\t" + << st.p_lrt << endl; + break; + case 3: + outfile << st.p_score << endl; + break; + case 4: + outfile << st.lambda_remle << "\t" + << st.lambda_mle << "\t" + << st.p_wald << "\t" + << st.p_lrt << "\t" + << st.p_score << endl; + break; + } + }; + + + if (!file_gene.empty()) { + outfile << "geneID" << "\t"; + + common_header(); for (vector<SUMSTAT>::size_type t = 0; t < sumStat.size(); ++t) { outfile << snpInfo[t].rs_number << "\t"; - - if (a_mode == 1) { - outfile << scientific << setprecision(6) << sumStat[t].beta << "\t" - << sumStat[t].se << "\t" << sumStat[t].lambda_remle << "\t" - << sumStat[t].p_wald << endl; - } else if (a_mode == 2) { - outfile << scientific << setprecision(6) << sumStat[t].lambda_mle - << "\t" << sumStat[t].p_lrt << endl; - } else if (a_mode == 3) { - outfile << scientific << setprecision(6) << sumStat[t].beta << "\t" - << sumStat[t].se << "\t" << sumStat[t].p_score << endl; - } else if (a_mode == 4) { - outfile << scientific << setprecision(6) << sumStat[t].beta << "\t" - << sumStat[t].se << "\t" << sumStat[t].lambda_remle << "\t" - << sumStat[t].lambda_mle << "\t" << sumStat[t].p_wald << "\t" - << sumStat[t].p_lrt << "\t" << sumStat[t].p_score << endl; - } else { - } + sumstats(sumStat[t]); } } else { bool process_gwasnps = setGWASnps.size(); - outfile << "chr" - << "\t" - << "rs" - << "\t" - << "ps" - << "\t" - << "n_miss" - << "\t" - << "allele1" - << "\t" - << "allele0" - << "\t" - << "af" - << "\t"; - - if (a_mode == 1) { - outfile << "beta" - << "\t" - << "se" - << "\t" - << "l_remle" - << "\t" - << "p_wald" << endl; - } else if (a_mode == 2) { - outfile << "l_mle" - << "\t" - << "p_lrt" << endl; - } else if (a_mode == 3) { - outfile << "beta" - << "\t" - << "se" - << "\t" - << "p_score" << endl; - } else if (a_mode == 4) { - outfile << "beta" - << "\t" - << "se" - << "\t" - << "l_remle" - << "\t" - << "l_mle" - << "\t" - << "p_wald" - << "\t" - << "p_lrt" - << "\t" - << "p_score" << endl; - } else { - } + + outfile << "chr" << "\t" + << "rs" << "\t" + << "ps" << "\t" + << "n_miss" << "\t" + << "allele1" << "\t" + << "allele0" << "\t" + << "af" << "\t"; + + common_header(); size_t t = 0; for (size_t i = 0; i < snpInfo.size(); ++i) { - if (indicator_snp[i] == 0) continue; auto snp = snpInfo[i].rs_number; if (process_gwasnps && setGWASnps.count(snp) == 0) continue; // cout << t << endl; - outfile << snpInfo[i].chr << "\t" << snpInfo[i].rs_number << "\t" << snpInfo[i].base_position << "\t" << snpInfo[i].n_miss << "\t" << snpInfo[i].a_minor << "\t" << snpInfo[i].a_major << "\t" << fixed << setprecision(3) << snpInfo[i].maf << "\t"; - if (a_mode == 1) { - outfile << scientific << setprecision(6) << sumStat[t].beta << "\t" - << sumStat[t].se << "\t" << sumStat[t].lambda_remle << "\t" - << sumStat[t].p_wald << endl; - } else if (a_mode == 2) { - outfile << scientific << setprecision(6) << sumStat[t].lambda_mle - << "\t" << sumStat[t].p_lrt << endl; - } else if (a_mode == 3) { - outfile << scientific << setprecision(6) << sumStat[t].beta << "\t" - << sumStat[t].se << "\t" << sumStat[t].p_score << endl; - } else if (a_mode == 4) { - outfile << scientific << setprecision(6) << sumStat[t].beta << "\t" - << sumStat[t].se << "\t" << sumStat[t].lambda_remle << "\t" - << sumStat[t].lambda_mle << "\t" << sumStat[t].p_wald << "\t" - << sumStat[t].p_lrt << "\t" << sumStat[t].p_score << endl; - } else { - } + sumstats(sumStat[t]); t++; } } @@ -1299,7 +1249,7 @@ void LMM::AnalyzeGene(const gsl_matrix *U, const gsl_vector *eval, time_opt += (clock() - time_start) / (double(CLOCKS_PER_SEC) * 60.0); // Store summary data. - SUMSTAT SNPs = {beta, se, lambda_remle, lambda_mle, p_wald, p_lrt, p_score}; + SUMSTAT SNPs = {beta, se, lambda_remle, lambda_mle, p_wald, p_lrt, p_score, logl_H1}; sumStat.push_back(SNPs); } cout << endl; @@ -1400,7 +1350,7 @@ void LMM::AnalyzeBimbam(const gsl_matrix *U, const gsl_vector *eval, // Store summary data. SUMSTAT SNPs = {beta, se, lambda_remle, lambda_mle, - p_wald, p_lrt, p_score}; + p_wald, p_lrt, p_score, logl_H1}; sumStat.push_back(SNPs); } }; @@ -1653,7 +1603,7 @@ void LMM::AnalyzePlink(const gsl_matrix *U, const gsl_vector *eval, // Store summary data. SUMSTAT SNPs = {beta, se, lambda_remle, lambda_mle, - p_wald, p_lrt, p_score}; + p_wald, p_lrt, p_score, logl_H1}; sumStat.push_back(SNPs); } } @@ -1930,7 +1880,7 @@ void LMM::Analyzebgen(const gsl_matrix *U, const gsl_vector *eval, // Store summary data. SUMSTAT SNPs = {beta, se, lambda_remle, lambda_mle, - p_wald, p_lrt, p_score}; + p_wald, p_lrt, p_score, logl_H1}; sumStat.push_back(SNPs); } } @@ -2411,7 +2361,7 @@ void LMM::AnalyzeBimbamGXE(const gsl_matrix *U, const gsl_vector *eval, time_opt += (clock() - time_start) / (double(CLOCKS_PER_SEC) * 60.0); // Store summary data. - SUMSTAT SNPs = {beta, se, lambda_remle, lambda_mle, p_wald, p_lrt, p_score}; + SUMSTAT SNPs = {beta, se, lambda_remle, lambda_mle, p_wald, p_lrt, p_score, logl_H1}; sumStat.push_back(SNPs); } cout << endl; @@ -2589,7 +2539,7 @@ void LMM::AnalyzePlinkGXE(const gsl_matrix *U, const gsl_vector *eval, time_opt += (clock() - time_start) / (double(CLOCKS_PER_SEC) * 60.0); // Store summary data. - SUMSTAT SNPs = {beta, se, lambda_remle, lambda_mle, p_wald, p_lrt, p_score}; + SUMSTAT SNPs = {beta, se, lambda_remle, lambda_mle, p_wald, p_lrt, p_score, logl_H1}; sumStat.push_back(SNPs); } cout << endl; diff --git a/src/param.h b/src/param.h index 08b1e10..ff279bd 100644 --- a/src/param.h +++ b/src/param.h @@ -56,6 +56,9 @@ public: double p_wald; // p value from a Wald test. double p_lrt; // p value from a likelihood ratio test. double p_score; // p value from a score test. + double logl_H1; // log likelihood under the alternative + // hypothesis as a measure of goodness of fit, + // see https://github.com/genetics-statistics/GEMMA/issues/81 }; // Results for mvLMM. @@ -118,7 +121,7 @@ public: bool mode_debug = false; uint issue; // enable tests for issue on github tracker - int a_mode; // Analysis mode, 1/2/3/4 for Frequentist tests + uint a_mode; // Analysis mode, 1/2/3/4 for Frequentist tests int k_mode; // Kinship read mode: 1: n by n matrix, 2: id/id/k_value; vector<size_t> p_column; // Which phenotype column needs analysis. size_t d_pace; // Display pace |