9 void UnivariateData::Apply(
double left_bound,
double right_bound, uword bound_window,
const mat &spectra,
const vec &abscissa)
11 cout <<
"UnivariateData::Apply\n";
15 field<mat> inflection_baselines;
19 left_bound, right_bound,
20 bound_window, baselines,
21 inflection_baselines);
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();
34 QStringList keys({
"Peak Centers",
36 "Adjusted Peak Intensities",
39 "Area Between Inflection Points",
40 "Adjusted Area Between Inflection Points",
41 "Empirical Full-Width at Half-Maximum"});
44 AddMatrix(
"Inflection Points", inflection_points);
47 AddField(
"Inflection Baselines", inflection_baselines);
56 mat results, baselines, fits, params, residuals;
57 if (peak_shape ==
"Gaussian")
59 else if (peak_shape ==
"Lorentzian")
61 left_bound, right_bound,
66 else if (peak_shape ==
"Voigt")
75 else throw invalid_argument(
"UnivariateData::Apply: improper peak shape specified");
78 if (peak_shape ==
"Voigt"){
79 columns = QStringList({
"Intensity",
81 "Gaussian Full-Width at Half Maximum",
82 "Lorentzian Full-Width at Half Maximum",
83 "Full-Width at Half Maximum",
88 columns = QStringList({
"Intensity",
90 "Full-Width at Half Maximum",
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)
103 mat first_baselines, second_baselines;
104 field<mat> inflection_first_baselines, inflection_second_baselines;
107 first_left_bound, first_right_bound,
108 bound_window, first_baselines,
109 inflection_first_baselines);
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);
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();
125 mat ratios = first_results / second_results;
126 QStringList headings({
"Peak Centers",
128 "Adjusted Peak Intensities",
131 "Area Between Inflection Points",
132 "Adjusted Area Between Inflection Points",
133 "Empirical Full-Width at Half-Maximum"});
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);
154 AddMatrix(
"Correlation Coefficients", results);
168 if (column >=
GetMatrix(
"Band Ratios").n_cols || column >=
GetMatrix(
"Results").n_cols || column >=
GetMatrix(
"Correlation Coefficients").n_cols)
173 else if (
HasMatrix(
"Correlation Coefficients")) data =
GetMatrix(
"Correlation Coefficients").col(0);
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 ¶ms, 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)
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 ¶ms, 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 ¶ms, arma::mat &residuals)
Vespucci::Math::Quantification::FitLorentzianPeakMat.