about summary refs log tree commit diff
path: root/src
diff options
context:
space:
mode:
authorPjotr Prins2017-10-06 11:51:18 +0000
committerPjotr Prins2017-10-13 15:27:24 +0000
commit63b6796a16ed2c82bb8b08eb4685bf8d56a9d360 (patch)
treebe1659d5f3770364d6a0773d9913e25f13d74485 /src
parentec8a139cb07b46bd0b9e5de2ea8db1f7f335a56c (diff)
downloadpangemma-63b6796a16ed2c82bb8b08eb4685bf8d56a9d360.tar.gz
Consolidate into ProgressBar into one function and related updates
Diffstat (limited to 'src')
-rw-r--r--src/io.cpp54
-rw-r--r--src/io.h4
-rw-r--r--src/lm.cpp6
-rw-r--r--src/mvlmm.cpp8
-rw-r--r--src/param.cpp4
-rw-r--r--src/param.h3
6 files changed, 27 insertions, 52 deletions
diff --git a/src/io.cpp b/src/io.cpp
index 8d77b0f..71d55c3 100644
--- a/src/io.cpp
+++ b/src/io.cpp
@@ -49,43 +49,17 @@
 using namespace std;
 
 // Print progress bar.
-void ProgressBar(string str, double p, double total) {
-  double progress = (100.0 * p / total);
-  int barsize = (int)(progress / 2.0);
-  char bar[51];
-
-  cout << str;
-  for (int i = 0; i < 50; i++) {
-    if (i < barsize) {
-      bar[i] = '=';
-    } else {
-      bar[i] = ' ';
-    }
-    cout << bar[i];
-  }
-  cout << setprecision(2) << fixed << progress << "%\r" << flush;
-
-  return;
-}
-
-// Print progress bar with acceptance ratio.
 void ProgressBar(string str, double p, double total, double ratio) {
-  double progress = (100.0 * p / total);
-  int barsize = (int)(progress / 2.0);
-  char bar[51];
-
-  cout << str;
-  for (int i = 0; i < 50; i++) {
-    if (i < barsize) {
-      bar[i] = '=';
-    } else {
-      bar[i] = ' ';
-    }
-    cout << bar[i];
-  }
-  cout << setprecision(2) << fixed << progress << "%    " << ratio << "\r"
-       << flush;
-  return;
+  assert(p<=total);
+  const double progress = (100.0 * p / total);
+  const uint barsize = (int)(progress / 2.0); // characters
+  cout << str << " ";
+  cout << std::string(barsize,'=');
+  cout << std::string(50-barsize,' ');
+  cout << setprecision(0) << fixed << progress << "%";
+  if (ratio != -1.0)
+    cout << setprecision(2) << "    " << ratio;
+  cout << "\r" << flush;
 }
 
 bool isBlankLine(char const *line) {
@@ -1407,7 +1381,7 @@ bool BimbamKin(const string file_geno, const set<string> ksnps,
     string line;
     !safeGetline(infile, line).eof();
     if (t % display_pace == 0 || t == (indicator_snp.size() - 1)) {
-      ProgressBar("Reading SNPs  ", t, indicator_snp.size() - 1);
+      ProgressBar("Reading SNPs", t, indicator_snp.size() - 1);
     }
     if (indicator_snp[t] == 0)
       continue;
@@ -1556,7 +1530,7 @@ bool PlinkKin(const string &file_bed, vector<int> &indicator_snp,
 
   for (size_t t = 0; t < indicator_snp.size(); ++t) {
     if (t % display_pace == 0 || t == (indicator_snp.size() - 1)) {
-      ProgressBar("Reading SNPs  ", t, indicator_snp.size() - 1);
+      ProgressBar("Reading SNPs", t, indicator_snp.size() - 1);
     }
     if (indicator_snp[t] == 0) {
       continue;
@@ -2713,7 +2687,7 @@ bool BimbamKinUncentered(const string &file_geno, const set<string> ksnps,
   for (size_t t = 0; t < indicator_snp.size(); ++t) {
     !safeGetline(infile, line).eof();
     if (t % display_pace == 0 || t == (indicator_snp.size() - 1)) {
-      ProgressBar("Reading SNPs  ", t, indicator_snp.size() - 1);
+      ProgressBar("Reading SNPs", t, indicator_snp.size() - 1);
     }
     if (indicator_snp[t] == 0)
       continue;
@@ -2918,7 +2892,7 @@ bool PlinkKin(const string &file_bed, const int display_pace,
 
   for (size_t t = 0; t < indicator_snp.size(); ++t) {
     if (t % display_pace == 0 || t == (indicator_snp.size() - 1)) {
-      ProgressBar("Reading SNPs  ", t, indicator_snp.size() - 1);
+      ProgressBar("Reading SNPs", t, indicator_snp.size() - 1);
     }
     if (indicator_snp[t] == 0) {
       continue;
diff --git a/src/io.h b/src/io.h
index 1c187b8..c65a2b0 100644
--- a/src/io.h
+++ b/src/io.h
@@ -32,8 +32,8 @@
 
 using namespace std;
 
-void ProgressBar(string str, double p, double total);
-void ProgressBar(string str, double p, double total, double ratio);
+void ProgressBar(string str, double p, double total, double ratio = -1.0);
+
 std::istream &safeGetline(std::istream &is, std::string &t);
 
 bool ReadFile_snps(const string file_snps, set<string> &setSnps);
diff --git a/src/lm.cpp b/src/lm.cpp
index a44bceb..ea628cf 100644
--- a/src/lm.cpp
+++ b/src/lm.cpp
@@ -331,7 +331,7 @@ void LM::AnalyzeGene(const gsl_matrix *W, const gsl_vector *x) {
   for (size_t t = 0; t < ng_total; t++) {
     getline(infile, line);
     if (t % d_pace == 0 || t == ng_total - 1) {
-      ProgressBar("Performing Analysis ", t, ng_total - 1);
+      ProgressBar("Performing Analysis", t, ng_total - 1);
     }
     ch_ptr = strtok((char *)line.c_str(), " , \t");
     rs = ch_ptr;
@@ -421,7 +421,7 @@ void LM::AnalyzeBimbam(const gsl_matrix *W, const gsl_vector *y) {
   for (size_t t = 0; t < indicator_snp.size(); ++t) {
     getline(infile, line);
     if (t % d_pace == 0 || t == (ns_total - 1)) {
-      ProgressBar("Reading SNPs  ", t, ns_total - 1);
+      ProgressBar("Reading SNPs", t, ns_total - 1);
     }
     if (indicator_snp[t] == 0) {
       continue;
@@ -547,7 +547,7 @@ void LM::AnalyzePlink(const gsl_matrix *W, const gsl_vector *y) {
 
   for (vector<SNPINFO>::size_type t = 0; t < snpInfo.size(); ++t) {
     if (t % d_pace == 0 || t == snpInfo.size() - 1) {
-      ProgressBar("Reading SNPs  ", t, snpInfo.size() - 1);
+      ProgressBar("Reading SNPs", t, snpInfo.size() - 1);
     }
     if (indicator_snp[t] == 0) {
       continue;
diff --git a/src/mvlmm.cpp b/src/mvlmm.cpp
index 88df111..c0dc143 100644
--- a/src/mvlmm.cpp
+++ b/src/mvlmm.cpp
@@ -3190,7 +3190,7 @@ void MVLMM::AnalyzeBimbam(const gsl_matrix *U, const gsl_vector *eval,
   for (size_t t = 0; t < indicator_snp.size(); ++t) {
     !safeGetline(infile, line).eof();
     if (t % d_pace == 0 || t == (ns_total - 1)) {
-      ProgressBar("Reading SNPs  ", t, ns_total - 1);
+      ProgressBar("Reading SNPs", t, ns_total - 1);
     }
     if (indicator_snp[t] == 0) {
       continue;
@@ -3639,7 +3639,7 @@ void MVLMM::AnalyzePlink(const gsl_matrix *U, const gsl_vector *eval,
   }
   for (vector<SNPINFO>::size_type t = 0; t < snpInfo.size(); ++t) {
     if (t % d_pace == 0 || t == snpInfo.size() - 1) {
-      ProgressBar("Reading SNPs  ", t, snpInfo.size() - 1);
+      ProgressBar("Reading SNPs", t, snpInfo.size() - 1);
     }
     if (indicator_snp[t] == 0) {
       continue;
@@ -4167,7 +4167,7 @@ void MVLMM::AnalyzeBimbamGXE(const gsl_matrix *U, const gsl_vector *eval,
   for (size_t t = 0; t < indicator_snp.size(); ++t) {
     !safeGetline(infile, line).eof();
     if (t % d_pace == 0 || t == (ns_total - 1)) {
-      ProgressBar("Reading SNPs  ", t, ns_total - 1);
+      ProgressBar("Reading SNPs", t, ns_total - 1);
     }
     if (indicator_snp[t] == 0) {
       continue;
@@ -4624,7 +4624,7 @@ void MVLMM::AnalyzePlinkGXE(const gsl_matrix *U, const gsl_vector *eval,
 
   for (vector<SNPINFO>::size_type t = 0; t < snpInfo.size(); ++t) {
     if (t % d_pace == 0 || t == snpInfo.size() - 1) {
-      ProgressBar("Reading SNPs  ", t, snpInfo.size() - 1);
+      ProgressBar("Reading SNPs", t, snpInfo.size() - 1);
     }
     if (indicator_snp[t] == 0) {
       continue;
diff --git a/src/param.cpp b/src/param.cpp
index 8452fb8..d637481 100644
--- a/src/param.cpp
+++ b/src/param.cpp
@@ -87,7 +87,7 @@ void trim_individuals(vector<int> &idvs, size_t ni_max, bool debug) {
 // ---- PARAM class implementation
 
 PARAM::PARAM(void)
-    : mode_silence(false), a_mode(0), k_mode(1), d_pace(100000),
+    : mode_silence(false), a_mode(0), k_mode(1), d_pace(DEFAULT_PACE),
       file_out("result"), path_out("./output/"), miss_level(0.05),
       maf_level(0.01), hwe_level(0), r2_level(0.9999), l_min(1e-5), l_max(1e5),
       n_region(10), p_nr(0.001), em_prec(0.0001), nr_prec(0.0001),
@@ -1150,7 +1150,7 @@ void PARAM::CheckData(void) {
   }
 
   // Set d_pace to 1000 for gene expression.
-  if (!file_gene.empty() && d_pace == 100000) {
+  if (!file_gene.empty() && d_pace == DEFAULT_PACE) {
     d_pace = 1000;
   }
 
diff --git a/src/param.h b/src/param.h
index 3976440..4b473c0 100644
--- a/src/param.h
+++ b/src/param.h
@@ -27,6 +27,7 @@
 #include <vector>
 
 #define K_BATCH_SIZE 10000 // #snps used for batched K
+#define DEFAULT_PACE 1000
 
 using namespace std;
 
@@ -124,7 +125,7 @@ public:
   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
+  size_t d_pace = DEFAULT_PACE;   // Display pace (-pace switch)
 
   string file_bfile, file_mbfile;
   string file_geno, file_mgeno;