Vespucci  1.0.0
univariatedata.cpp
Go to the documentation of this file.
2 
4  AnalysisResults(name, "Univariate Analysis Results")
5 {
6 
7 }
8 
9 void UnivariateData::Apply(double left_bound, double right_bound, uword bound_window, const mat &spectra, const vec &abscissa)
10 {
11  cout << "UnivariateData::Apply\n";
12 
13 
14  mat baselines;
15  field<mat> inflection_baselines;
16 
18  abscissa,
19  left_bound, right_bound,
20  bound_window, baselines,
21  inflection_baselines);
22 
23 
24  AddMetadata("Type", "One-region");
25  AddMetadata("Range", QString::number(left_bound) + "–" + QString::number(right_bound));
26  AddMetadata("Search window size", QString::number(bound_window));
27  mat inflection_points(spectra.n_cols, 2);
28  for (arma::uword i = 0; i < spectra.n_cols; ++i){
29  inflection_points(i, 0) = inflection_baselines(i).col(0).min();
30  inflection_points(i, 1) = inflection_baselines(i).col(0).max();
31  }
32 
33 
34  QStringList keys({ "Peak Centers",
35  "Peak Intensities",
36  "Adjusted Peak Intensities",
37  "Total Area",
38  "Adjusted Area",
39  "Area Between Inflection Points",
40  "Adjusted Area Between Inflection Points",
41  "Empirical Full-Width at Half-Maximum"});
42 
43  AddColumns(keys, results);
44  AddMatrix("Inflection Points", inflection_points);
45  AddMatrix("Baselines", baselines);
46  AddMatrix("Inflection Baselines", Vespucci::Math::Quantification::ConvertInflectionBaselines(inflection_baselines));
47  AddField("Inflection Baselines", inflection_baselines);
48 }
49 
50 void UnivariateData::Apply(QString peak_shape,
51  double left_bound,
52  double right_bound,
53  const mat &spectra,
54  const vec &abscissa)
55 {
56  mat results, baselines, fits, params, residuals;
57  if (peak_shape == "Gaussian")
58  results = Vespucci::Math::Quantification::FitGaussianPeakMat(spectra, abscissa, left_bound, right_bound, baselines, fits, params, residuals);
59  else if (peak_shape == "Lorentzian")
60  results = Vespucci::Math::Quantification::FitLorentzianPeakMat(spectra, abscissa,
61  left_bound, right_bound,
62  baselines,
63  fits,
64  params,
65  residuals);
66  else if (peak_shape == "Voigt")
68  abscissa,
69  left_bound,
70  right_bound,
71  baselines,
72  fits,
73  params,
74  residuals);
75  else throw invalid_argument("UnivariateData::Apply: improper peak shape specified");
76 
77  QStringList columns;
78  if (peak_shape == "Voigt"){
79  columns = QStringList({"Intensity",
80  "Area",
81  "Gaussian Full-Width at Half Maximum",
82  "Lorentzian Full-Width at Half Maximum",
83  "Full-Width at Half Maximum",
84  "Peak Centers",
85  "R²"});
86  }
87  else{
88  columns = QStringList({"Intensity",
89  "Area",
90  "Full-Width at Half Maximum",
91  "Peak Centers",
92  "R²"});
93  }
94  AddColumns(columns, results);
95  AddMatrix("Baselines", baselines);
96  AddMatrix("Fits", fits);
97  AddMatrix("Parameters", params);
98  AddMatrix("Residuals", residuals);
99 }
100 
101 void UnivariateData::Apply(double first_left_bound, double first_right_bound, double second_left_bound, double second_right_bound, uword bound_window, const mat &spectra, const vec &abscissa)
102 {
103  mat first_baselines, second_baselines;
104  field<mat> inflection_first_baselines, inflection_second_baselines;
105  mat first_results = Vespucci::Math::Quantification::QuantifyPeakMat(spectra,
106  abscissa,
107  first_left_bound, first_right_bound,
108  bound_window, first_baselines,
109  inflection_first_baselines);
110  mat second_results = Vespucci::Math::Quantification::QuantifyPeakMat(spectra,
111  abscissa,
112  second_left_bound, second_right_bound,
113  bound_window, second_baselines,
114  inflection_second_baselines);
115  mat first_inflection_points(spectra.n_cols, 2);
116  mat second_inflection_points(spectra.n_cols, 2);
117 
118  for (arma::uword i = 0; i < spectra.n_cols; ++i){
119  first_inflection_points(i, 0) = inflection_first_baselines(i).col(0).min();
120  first_inflection_points(i, 1) = inflection_first_baselines(i).col(0).max();
121  second_inflection_points(i,0) = inflection_second_baselines(i).col(0).min();
122  second_inflection_points(i,1) = inflection_second_baselines(i).col(0).max();
123  }
124 
125  mat ratios = first_results / second_results;
126  QStringList headings({ "Peak Centers",
127  "Peak Intensities",
128  "Adjusted Peak Intensities",
129  "Total Area",
130  "Adjusted Area",
131  "Area Between Inflection Points",
132  "Adjusted Area Between Inflection Points",
133  "Empirical Full-Width at Half-Maximum"});
134 
135  AddMetadata("Type", "Band ratio analysis");
136  AddMetadata("First region range", QString::number(first_left_bound) + "–" + QString::number(first_right_bound));
137  AddMetadata("Second region range", QString::number(second_left_bound) + "–" + QString::number(second_right_bound));
138  AddMatrix("First Region Results", first_results, headings);
139  AddMatrix("Second Region Results", second_results, headings);
140  AddMatrix("Band Ratios", ratios, headings);
141  AddMatrix("First Region Baselines", first_baselines);
142  AddMatrix("Second Region Baselines", second_baselines);
143  AddMatrix("First Region Inflection Points", first_inflection_points);
144  AddMatrix("Second Region Inflection Points", second_inflection_points);
145  AddMatrix("First Region Inflection Baselines", Vespucci::Math::Quantification::ConvertInflectionBaselines(inflection_first_baselines));
146  AddMatrix("Second Region Inflection Baselines", Vespucci::Math::Quantification::ConvertInflectionBaselines(inflection_second_baselines));
147 
148 }
149 
150 void UnivariateData::ApplyCorrelation(const mat &spectra, const vec &control)
151 {
152  AddMetadata("Type", "Correlation analysis");
153  mat results = Vespucci::Math::Quantification::CorrelationMat(spectra, control);
154  AddMatrix("Correlation Coefficients", results);
155 }
156 
164 void UnivariateData::Calibrate(const vec &x, const vec &y, uword column)
165 {
166  if (!HasMatrix("Band Ratios") && !HasMatrix("Results") && !HasMatrix("Correlation Coefficients"))
167  return;
168  if (column >= GetMatrix("Band Ratios").n_cols || column >= GetMatrix("Results").n_cols || column >= GetMatrix("Correlation Coefficients").n_cols)
169  return;
170 
171  vec data;
172  if (HasMatrix("Band Ratios")) data = GetMatrix("Band Ratios").col(column);
173  else if (HasMatrix("Correlation Coefficients")) data = GetMatrix("Correlation Coefficients").col(0);
174 
175 }
void AddMetadata(QString key, QString value)
const mat & GetMatrix(const QString &key)
VESPUCCI_EXPORT arma::mat FitGaussianPeakMat(const arma::mat &spectra, const arma::vec &abscissa, double &min, double &max, arma::mat &baselines, arma::mat &fits, arma::mat &params, arma::mat &residuals)
void Calibrate(const vec &x, const vec &y, uword column)
UnivariateData::Calibrate.
VESPUCCI_EXPORT arma::mat CorrelationMat(const arma::mat &X, const arma::mat &control)
Definition: correlation.cpp:22
void ApplyCorrelation(const mat &spectra, const vec &control)
void AddField(const QString &key, const field< mat > &value)
void Apply(double left_bound, double right_bound, uword bound_window, const mat &spectra, const vec &abscissa)
UnivariateData(QString name)
void AddMatrix(const QString &key, const mat &value, QStringList column_headings=QStringList())
void AddColumns(const QStringList &keys, const mat &value)
VESPUCCI_EXPORT arma::mat ConvertInflectionBaselines(const arma::field< arma::mat > &inflection_baselines)
Vespucci::Math::Quantification::ConvertInflectionBaselines.
VESPUCCI_EXPORT arma::mat FitVoigtPeakMat(const arma::mat &spectra, const arma::vec &abscissa, double &min, double &max, arma::mat &baselines, arma::mat &fits, arma::mat &params, arma::mat &residuals)
Vespucci::Math::Quantification::FitVoigtPeakMat.
The AnalysisResults class A container for a mat object that allows a mat to be copied to a heap-alloc...
bool HasMatrix(const QString &key) const
VESPUCCI_EXPORT arma::mat QuantifyPeakMat(const arma::mat &spectra, const arma::vec &abscissa, double &min, double &max, arma::uword bound_window, arma::mat &total_baselines, arma::field< arma::mat > &inflection_baselines)
VESPUCCI_EXPORT arma::mat FitLorentzianPeakMat(const arma::mat &spectra, const arma::vec &abscissa, double &min, double &max, arma::mat &baselines, arma::mat &fits, arma::mat &params, arma::mat &residuals)
Vespucci::Math::Quantification::FitLorentzianPeakMat.