Vespucci  1.0.0
multianalyzer.cpp
Go to the documentation of this file.
1 /*******************************************************************************
2  Copyright (C) 2014-2016 Wright State University - All Rights Reserved
3  Daniel P. Foose - Maintainer/Lead Developer
4 
5  This file is part of Vespucci.
6 
7  Vespucci is free software: you can redistribute it and/or modify
8  it under the terms of the GNU General Public License as published by
9  the Free Software Foundation, either version 3 of the License, or
10  (at your option) any later version.
11 
12  Vespucci is distributed in the hope that it will be useful,
13  but WITHOUT ANY WARRANTY; without even the implied warranty of
14  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  GNU General Public License for more details.
16 
17  You should have received a copy of the GNU General Public License
18  along with Vespucci. If not, see <http://www.gnu.org/licenses/>.
19 *******************************************************************************/
21 #include <Math/Smoothing/denoise.h>
22 MultiAnalyzer::MultiAnalyzer(QSharedPointer<VespucciWorkspace> workspace, QStringList dataset_keys)
23  : MetaAnalyzer(workspace), dataset_keys_(dataset_keys)
24 {
25  GetData();
26 }
27 
29 
30 bool MultiAnalyzer::CheckMergability(const QStringList dataset_keys)
31 {
32  for (auto key: dataset_keys)
33  if (workspace_->GetDataset(key)->spectra_ptr()->n_rows
34  != workspace_->GetDataset(dataset_keys.first())->spectra_ptr()->n_rows)
35  return false;
36  return true;
37 }
38 
39 bool MultiAnalyzer::ConcatenateSpectra(QStringList dataset_keys)
40 {
41  if (!CheckMergability(dataset_keys)) return false;
42 
43  try{
44  start_indices_.set_size(datasets_.size());
45  end_indices_.set_size(datasets_.size());
46  data_ = datasets_[0]->spectra_ref();
47  abscissa_ = datasets_[0]->abscissa_ref();
48  start_indices_(0) = 0;
49  end_indices_(0) = datasets_[0]->spectra_ref().n_cols - 1;
50 
51  for (uword i = 1; i < start_indices_.n_rows; ++ i){
52  start_indices_(i) = data_.n_cols;
53  data_ = join_horiz(data_, datasets_[i]->spectra_ref());
54  end_indices_(i) = data_.n_cols - 1;
55  }
56  }catch(...){
57  return false;
58  }
59 
60  return true;
61 }
62 
63 void MultiAnalyzer::SNVNormalize(double offset)
64 {
66  for (uword i = 0; i < start_indices_.n_rows; ++i){
67  mat spectra = data_.cols(start_indices_(i), end_indices_(i));
68  vec abscissa = datasets_[i]->abscissa();
69  vec x = datasets_[i]->x();
70  vec y = datasets_[i]->y();
71  datasets_[i]->SetData(spectra, abscissa, x, y);
72  }
73 }
74 void MultiAnalyzer::QUIC_SVD(double epsilon)
75 {
76  uword rank;
77  mat U, V;
78  vec s;
79  data_ = Vespucci::Math::Smoothing::QUICSVDDenoise(data_, epsilon, U, s, V, rank);
80  vec SVD_vec({double(rank)});
81 
82  QSharedPointer<AnalysisResults> results(new AnalysisResults("QUIC-SVD", "QUIC-SVD"));
83  results->AddMatrix("Left Singualr Vectors", U);
84  results->AddMatrix("Right Singualr Vectors", V);
85  results->AddMatrix("Singualr Values", s);
86  results->AddMatrix("Rank", SVD_vec);
87 
88  for (uword i = 0; i < start_indices_.n_rows; ++i){
89  mat spectra = data_.cols(start_indices_(i), end_indices_(i));
90  vec abscissa = datasets_[i]->abscissa();
91  vec x = datasets_[i]->x();
92  vec y = datasets_[i]->y();
93  datasets_[i]->SetData(spectra, abscissa, x, y);
94  }
95  AddAnalysisResults(results, QStringList());
96 }
98 {
99  mat U, V;
100  vec s;
102  vec SVD_vec({double(k)});
103  QSharedPointer<AnalysisResults> results(new AnalysisResults("QUIC-SVD", "QUIC-SVD"));
104  results->AddMatrix("Left Singualr Vectors", U);
105  results->AddMatrix("Right Singualr Vectors", V);
106  results->AddMatrix("Singualr Values", s);
107  results->AddMatrix("Rank", SVD_vec);
108  for (uword i = 0; i < start_indices_.n_rows; ++i){
109  mat spectra = data_.cols(start_indices_(i), end_indices_(i));
110  vec abscissa = datasets_[i]->abscissa();
111  vec x = datasets_[i]->x();
112  vec y = datasets_[i]->y();
113  datasets_[i]->SetData(spectra, abscissa, x, y);
114  }
115 
116  AddAnalysisResults(results, QStringList());
117 }
118 
119 
126 void MultiAnalyzer::GetDatasets(QStringList keys)
127 {
128  for (auto key: keys){
129  if (workspace_->dataset_names().contains(key)){
130  datasets_.append(workspace_->GetDataset(key));
131  }
132  }
133 }
134 
138 void MultiAnalyzer::GetData()
139 {
140  GetDatasets(dataset_keys_);
141  bool ok = ConcatenateSpectra(dataset_keys_);
142  if (!ok) throw runtime_error("Could not concatenate datasets");
143 }
144 
145 void MultiAnalyzer::AddAnalysisResults(QSharedPointer<AnalysisResults> results, QStringList matrices)
146 {
147  for (uword i = 0; i < uword(datasets_.size()); ++i){
148  datasets_[i]->AddAnalysisResult(results, start_indices_(i), end_indices_(i));
149  if (!matrices.isEmpty())
150  datasets_[i]->AddAnalysisResult(results->Subset(matrices, start_indices_(i), end_indices_(i)));
151  }
152 }
153 
154 QString MultiAnalyzer::FindUniqueName(QString name)
155 {
156  QString new_name = name;
157  QStringList results_names;
158  for (auto dataset: datasets_){
159  results_names = results_names + dataset->AnalysisResultsKeys();
160  }
161  int i = 0;
162  QString basename = name;
163  while (results_names.contains(new_name)){
164  new_name = basename + " (" + QString::number(i++) + ")";
165  }
166  return new_name;
167 }
QSharedPointer< VespucciWorkspace > workspace_
workspace_ The global workspace, used to obtain access to matrices
Definition: metaanalyzer.h:100
void SVDDenoise(uword k)
void QUIC_SVD(double epsilon)
VESPUCCI_EXPORT arma::mat SNVNorm(const arma::mat &X, const double offset, bool center)
Vespucci::Math::Normalization::SNVNorm.
VESPUCCI_EXPORT arma::mat QUICSVDDenoise(const arma::mat &X, double epsilon, arma::mat &U, arma::vec &s, arma::mat &V, arma::uword &rank)
Definition: denoise.cpp:28
VESPUCCI_EXPORT arma::mat SVDDenoise(const arma::mat &X, arma::uword k, arma::mat &U, arma::vec &s, arma::mat &V)
Definition: denoise.cpp:22
The MetaAnalyzer class This class is used to perform analysis on "pseudo-datasets". This allows us perform analysis on a concatentated spectra matrix from multiple datasets (the MultiAnalyzer derived class) or on a single matrix in a single dataset (the MatrixAnalyzer derived class). This class is not concerned with how the data_ and abscissa_ members are initialized (this is handled by GetData in derived classes).
Definition: metaanalyzer.h:33
arma::vec abscissa() const
MultiAnalyzer(QSharedPointer< VespucciWorkspace > workspace, QStringList dataset_keys)
The AnalysisResults class A container for a mat object that allows a mat to be copied to a heap-alloc...
void SNVNormalize(double offset)