Vespucci  1.0.0
VespucciDataset Class Reference

The VespucciDataset class This is the main class for dealing with hyperspectral data. This handles the import and export of spectra, and the creation of maps. Images are handled by the MapData class. This class is intended to be allocated on the heap inside of a smart pointer, there is no copy constructor. More...

#include <vespuccidataset.h>

Inheritance diagram for VespucciDataset:
AbstractDataAnalyzer MetaDataset

Public Member Functions

 VespucciDataset (const QString &h5_filename, MainWindow *main_window, QSharedPointer< VespucciWorkspace > ws)
 VespucciDataset::VespucciDataset. More...
 
 VespucciDataset (QString text_filename, MainWindow *main_window, QString *directory, QString name, QString x_axis_description, QString y_axis_description, bool swap_spatial, std::string format)
 VespucciDataset::VespucciDataset. More...
 
 VespucciDataset (map< pair< int, int >, string > text_filenames, MainWindow *main_window, QString *directory, QString name, QString x_axis_description, QString y_axis_description, int rows, int cols)
 
 VespucciDataset (QString name, MainWindow *main_window, QString *directory, QSharedPointer< VespucciDataset > original, uvec indices)
 VespucciDataset::VespucciDataset. More...
 
 VespucciDataset (QString name, MainWindow *main_window, QString *directory)
 VespucciDataset::VespucciDataset. More...
 
 ~VespucciDataset ()
 VespucciDataset::~VespucciDataset. More...
 
void Undo ()
 VespucciDataset::Undo Swap spectra_ and spectra_old_ to undo an action. Calling this function again re-does the action that was undid. More...
 
void CropSpectra (double x_min, double x_max, double y_min, double y_max, double wl_min, double wl_max)
 VespucciDataset::CropSpectra Crops spectra_ based on. More...
 
mat ZScoreNormCopy ()
 VespucciDataset::ZScoreNormCopy For when you want to Z-score normalize without changing spectra_. More...
 
void MinMaxNormalize ()
 VespucciDataset::MinMaxNormalize Normalizes each column of spectra_ so that the highest value is 1 and lowest value is 0. More...
 
void VectorNormalize (uword norm)
 VespucciDataset::VectorNormalize. More...
 
void MeanCenter ()
 VespucciDataset::MeanCenter Subtract the mean of all spectra in the dataset at each wavelength from each wavelength (i.e. in preparation for PCA) More...
 
void ZScoreNormalize ()
 VespucciDataset::ZScoreNormalize Normalizes each row using the score of the normal distribution. More...
 
void SNVNormalize (double offset, bool center)
 
void AbsoluteValue ()
 
void PeakIntensityNormalize (double peak_position)
 
void PeakIntensityNormalize (double left_bound, double right_bound)
 VespucciDataset::PeakIntensityNormalize. More...
 
void Booleanize (double min, double max, bool keep_inside, bool oneify)
 VespucciDataset::Booleanize. More...
 
void Clamp (double min, double max)
 VespucciDataset::Clamp. More...
 
void ShedZeroSpectra ()
 
void ShedZeroWavelengths ()
 
void SubtractBackground (const QStringList &data_keys)
 VespucciDataset::SubtractBackground. More...
 
void MedianFilter (unsigned int window_size)
 VespucciDataset::MedianFilter performs median filtering on the spectral data. Entries near the boundaries of spectra are not processed. See also VespucciDataset::LinearMovingAverage. More...
 
void LinearMovingAverage (unsigned int window_size)
 VespucciDataset::LinearMovingAverage Performs moving average filtering on the spectral data. Entries near the boundaries of spectra are not processed. See also VespucciDataset::MedianFilter. More...
 
void SavitzkyGolay (unsigned int derivative_order, unsigned int polynomial_order, unsigned int window_size)
 VespucciDataset::SavitzkyGolay Performs derivatization/Savitzky-Golay smoothing. More...
 
void SingularValue (unsigned int singular_values)
 VespucciDataset::SingularValue Denoises the spectra matrix using a truncated singular value decomposition. The first singular_values singular values are used to "reconstruct" the spectra matrix. The function used to find the truncated SVD is Vespucci::Math::DimensionReduction::svds. More...
 
int QUIC_SVD (double epsilon)
 VespucciDataset::QUIC_SVD. More...
 
void MFBaseline (int window_size, int iterations)
 VespucciDataset::MFBaseline Baseline-adjusts the data. This function uses a median filter with a large window to determine the baseline on the assumption that the median value is more likely to be basline than spectrum. This will complicate things if there are many peaks. Additionally, this significantly narrows the shape of peaks. Other baseline methods may be implemented later. More...
 
void RollingBallBaseline (size_t wm, size_t ws)
 VespucciDataset::RollingBallBaseline. More...
 
void CWTBaseline (int lambda, int penalty_order, double SNR_threshold, double peak_shape_threshold)
 
void IModPolyBaseline (const uword poly_order, const uword max_it, double threshold)
 
void RemoveClippedSpectra (double threshold)
 VespucciDataset::RemoveClippedSpectra. More...
 
void RemoveFlatSpectra (double threshold)
 VespucciDataset::RemoveFlatSpectra. More...
 
void ZeroClippedSpectra (double threshold)
 VespucciDataset::ZeroClippedSpectra. More...
 
void ZeroFlatSpectra (double threshold)
 VespucciDataset::ZeroFlatSpectra. More...
 
void Scale (double scaling_factor)
 
void ShedSpectrum (const uword index)
 VespucciDataset::ShedSpectrum. More...
 
void ZeroSpectrum (const uword index)
 
int HySime ()
 VespucciDataset::HySime. More...
 
void TransformAbscissa (QString input_units, double input_factor, QString output_units, double output_factor, QString description)
 VespucciDataset::TransformAbscissa. More...
 
void InterpolateToNewAbscissa (const vec &new_abscissa, unsigned int poly_order, unsigned int window_size)
 VespucciDataset::InterpolateToNewAbscissa. More...
 
void InterpolateToNewAbscissa (const vec &new_abscissa)
 VespucciDataset::InterpolateToNewAbscissa. More...
 
void FourierTransform (int n)
 VespucciDataset::FourierTransform. More...
 
void InverseFourierTransform (int n)
 
void ApplyFTWeight (QString type, double param)
 VespucciDataset::ApplyFTWeight. More...
 
void ApplyFTWeight (double start_offset, double end_offset, double power)
 VespucciDataset::ApplyFTWeight. More...
 
vec PointSpectrum (uword index) const
 VespucciDataset::PointSpectrum. More...
 
QVector< double > WavelengthQVector () const
 
int FindIndex (int x_value, int y_value) const
 
uword FindIndex (double abscissa_value) const
 
uvec FindRange (double start, double end) const
 VespucciDataset::FindRange. Finds the index of the wavelength value closest to the specified wavelength range. More...
 
uword FindOrigin () const
 VespucciDataset::FindOrigin. More...
 
int KeySize () const
 VespucciDataset::KeySize Finds the number of data points in x variable to properly set axes of QCustomPlot objects. More...
 
int ValueSize () const
 VespucciDataset::ValueSize Finds number of unique y values for properly setting QCPAxis. More...
 
QCPRange KeyRange () const
 VespucciDataset::KeyRange Finds the minima and maxima of x variable to properly set axes of QCustomPlot objects. More...
 
QCPRange ValueRange () const
 VespucciDataset::ValueRange Finds the minima and maxima of y variable to properly set axes of QCustomPlot objects. More...
 
bool Save (QString filename)
 
bool Load (QString filename)
 
bool SaveSpectrum (QString filename, uword column, file_type type)
 
void Univariate (const QString &name, double &left_bound, double &right_bound, uword bound_window)
 VespucciDataset::Univariate. More...
 
void FitPeak (const QString &name, const QString &peak_shape, double &left_bound, double &right_bound)
 
void BandRatio (const QString &name, double &first_left_bound, double &first_right_bound, double &second_left_bound, double &second_right_bound, uword bound_window)
 VespucciDataset::BandRatio. More...
 
void ClassicalLeastSquares (const QString &name, const QStringList &reference_keys)
 
void PartialLeastSquares (const QString &name, uword components)
 VespucciDataset::PartialLeastSquares. More...
 
void PLSCalibration (const QString &name, const QStringList &control_keys)
 
void TrainPLSDA (const QString &name, const QStringList &label_keys)
 
void CorrelationAnalysis (const QString &control_key, QString name)
 VespucciDataset::CorrelationAnalysis. More...
 
void VertexComponents (const QString &name, uword endmembers)
 VespucciDataset::VertexComponents. More...
 
void KMeans (const QString &name, const QString &metric_text, const QString &partition_policy, bool allow_empty, size_t clusters)
 VespucciDataset::KMeans. More...
 
void PrincipalComponents (const QString &name)
 VespucciDataset::PrincipalComponents Perform PCA without creating an image. More...
 
void PrincipalComponents (const QString &name, bool scale_data)
 VespucciDataset::PrincipalComponents. More...
 
void FindPeaks (const QString &name, double sel, double threshold, uword poly_order, uword window_size)
 VespucciDataset::FindPeaks. More...
 
void AgglomerativeClustering (const QString &name, const QString &linkage, const QString &metric)
 VespucciDataset::AgglomerativeClustering. More...
 
void CalculateRepresentativeSpectrum (const QString &name, QString statistic, QString metric)
 VespucciDataset::CalculateRepresentativeSpectrum. More...
 
vec wavelength () const
 VespucciDataset::wavelength. More...
 
vec abscissa () const
 
vec x () const
 VespucciDataset::x. More...
 
vec y () const
 VespucciDataset::y. More...
 
vec x (uvec indices) const
 VespucciDataset::x. More...
 
vec y (uvec indices) const
 VespucciDataset::y. More...
 
double x (uword index) const
 
double y (uword index) const
 
vec wavelength (uvec indices) const
 
mat spectra (uvec indices) const
 VespucciDataset::spectra. More...
 
mat spectra () const
 VespucciDataset::spectra. More...
 
cx_mat cx_spectra () const
 
cx_mat cx_spectra (uvec indices) const
 
const QString name () const
 VespucciDataset::name. More...
 
vec indices () const
 VespucciDataset::indices. More...
 
mat * indices_ptr ()
 VespucciDataset::indices_ptr. More...
 
QStringList operations ()
 
bool principal_components_calculated () const
 
bool mlpack_pca_calculated () const
 
bool vertex_components_calculated () const
 
bool partial_least_squares_calculated () const
 
QList< QSharedPointer< UnivariateData > > univariate_datas ()
 
const QString x_axis_description () const
 VespucciDataset::x_axis_description The x_axis_description is printed on the spectrum viewer. More...
 
const QString y_axis_description () const
 VespucciDataset::y_axis_description. More...
 
void SetName (QString new_name)
 VespucciDataset::SetName. More...
 
void SetData (const mat &spectra, const vec &wavelength, const vec &x, const vec &y)
 VespucciDataset::SetData. More...
 
void SetIndices (vec indices)
 VespucciDataset::SetIndices. More...
 
void AddMap (QSharedPointer< MapData > map)
 VespucciDataset::AddMap Adds a map to the list of map pointers and adds its name to relevant lists. More...
 
void RemoveMap (const QString &name)
 
int map_loading_count () const
 VespucciDataset::map_names. More...
 
void SetXDescription (QString description)
 VespucciDataset::SetXDescription Sets the value of the spectral abscissa description. More...
 
void SetYDescription (QString description)
 VespucciDataset::SetYDescription. More...
 
QCPRange WavelengthRange () const
 VespucciDataset::WavelengthRange. More...
 
QCPRange PointSpectrumRange (int i) const
 VespucciDataset::PointSpectrumRange. More...
 
QCPColorGradient GetGradient (int gradient_number) const
 VespucciDataset::GetGradient Selects the color gradient from list of presets. More...
 
QCPColorGradient GetClusterGradient (int clusters) const
 VespucciDataset::GetClusterGradient Cluster gradients are slightly different from the continuous gradients. This selects the right gradient based on the number of clusters. More...
 
bool ConstructorCancelled () const
 VespucciDataset::ConstructorCancelled Specifies whether or not the constructor has been canceled. The constructor asks this and cleans everything up in case it is canceled. More...
 
mat AverageSpectrum (bool stats) const
 VespucciDataset::AverageSpectrum Finds the average of the spectrum. This can be saved by the user. Probably not all that useful, except for determining a spectrum to use as a background spectrum for other maps. More...
 
mat * x_ptr ()
 VespucciDataset::x_ptr. More...
 
mat * y_ptr ()
 VespucciDataset::y_ptr. More...
 
mat * wavelength_ptr ()
 
mat * abscissa_ptr ()
 VespucciDataset::abscissa_ptr. More...
 
mat * spectra_ptr ()
 VespucciDataset::spectra_ptr. More...
 
mat * spectra_old_ptr (bool *ok)
 
const mat & spectra_ref ()
 VespucciDataset::spectra_ref. More...
 
const vec & abscissa_ref ()
 VespucciDataset::abscissa_ref. More...
 
const vec & x_ref ()
 VespucciDataset::x_ref. More...
 
const vec & y_ref ()
 VespucciDataset::y_ref. More...
 
size_t columns () const
 
double AbscissaMin () const
 
double AbscissaMax () const
 
bool non_spatial () const
 VespucciDataset::non_spatial. More...
 
bool meta () const
 VespucciDataset::meta. More...
 
void SetParentDatasetIndices (mat parent_dataset_indices)
 VespucciDataset::SetParentDatasetIndices. More...
 
mat * parent_dataset_indices ()
 
bool Undoable () const
 VespucciDataset::Undoable. More...
 
uword UniqueX () const
 
uword UniqueY () const
 
int UnivariateCount () const
 VespucciDataset::UnivariateCount. More...
 
const QString last_operation () const
 VespucciDataset::last_operation. More...
 
MapListModelmap_list_model ()
 
void AddAnalysisResult (QSharedPointer< AnalysisResults > analysis_result)
 
void AddAnalysisResult (QSharedPointer< AnalysisResults > analysis_result, uword start_row, uword end_row)
 
QStringList AnalysisResultsKeys () const
 VespucciDataset::AnalysisResultsKeys. More...
 
QMap< QString, QStringList > AnalysisResultsTreeStructure () const
 VespucciDataset::AnalysisResultsTreeStructure. More...
 
void ImportAuxiliaryMatrix (const QString &name, const QString &filename)
 VespucciDataset::ImportAuxiliaryMatrix. More...
 
void AddAuxiliaryMatrix (const QString &name, mat &matrix)
 VespucciDataset::AddAuxiliaryMatrix. More...
 
void AddMatrix (const QString &name, mat &matrix)
 
QStringList AuxiliaryMatrixKeys () const
 VespucciDataset::AuxiliaryMatrixKeys. More...
 
QStringList CoreMatrixKeys () const
 VespucciDataset::CoreMatrixKeys. More...
 
const mat & GetAnalysisResultMatrix (const QString &results_key, const QString &matrix_key) const
 VespucciDataset::GetAnalysisResultMatrix. More...
 
QSharedPointer< AnalysisResultsGetAnalysisResult (const QString &key)
 
const mat & GetAuxiliaryMatrix (const QString &key) const
 
const mat & GetCoreMatrix (const QString &key) const
 
bool IsCoreMatrix (const QString &key) const
 VespucciDataset::IsCoreMatrix. More...
 
QSharedPointer< MapDataGetMapData (const QString &key)
 
void CreateMap (const QString &map_name, const QString &results_key, const QString &matrix_key, uword column, QCPColorGradient gradient)
 VespucciDataset::CreateMap. More...
 
void CreateMap (const QString &map_name, const QString &matrix_key, uword column, QCPColorGradient gradient)
 
bool ShowMapViewer (const QString &map_key, bool show)
 
QStringList MapKeys () const
 
void SetOldCopies ()
 
bool Contains (const QString &key) const
 VespucciDataset::Contains. More...
 
bool IsValid () const
 
bool state_changed () const
 
QString filename () const
 
QString last_save_filename () const
 
bool saved () const
 
- Public Member Functions inherited from AbstractDataAnalyzer
virtual ~AbstractDataAnalyzer ()
 
virtual void Univariate (const QString &name, double &left_bound, double &right_bound, arma::uword bound_window)=0
 
virtual void BandRatio (const QString &name, double &first_left_bound, double &first_right_bound, double &second_left_bound, double &second_right_bound, arma::uword bound_window)=0
 
virtual void VertexComponents (const QString &name, arma::uword endmembers)=0
 
virtual void PartialLeastSquares (const QString &name, arma::uword components)=0
 
virtual arma::vec PointSpectrum (arma::uword index) const =0
 

Public Attributes

QTextStream log_text_stream_
 log_text_stream_ A log of all UI-induced functions called on this dataset More...
 

Detailed Description

The VespucciDataset class This is the main class for dealing with hyperspectral data. This handles the import and export of spectra, and the creation of maps. Images are handled by the MapData class. This class is intended to be allocated on the heap inside of a smart pointer, there is no copy constructor.

Definition at line 75 of file vespuccidataset.h.

Constructor & Destructor Documentation

VespucciDataset::VespucciDataset ( const QString &  h5_filename,
MainWindow main_window,
QSharedPointer< VespucciWorkspace ws 
)

VespucciDataset::VespucciDataset.

Parameters
name
h5_filename
main_window
directoryConstructor to load dataset from an HDF5 file. Restrictions: Top level DataSets with names "Spectra" and "Spectral Abscissa" are required. Only three levels are recognized, Dataset, AnalysisResult, Matrix All matrices other than "x", "y", "abscissa" or "spectra" will be loaded into the auxiliary_matrices_ object.

All DataSets must have 2 dimensions only (column vectors stored as Nx1 matrices)

Attributes will be used to set metadata

All analysis results groups are added as base AnalysisResults

Definition at line 214 of file vespuccidataset.cpp.

VespucciDataset::VespucciDataset ( QString  text_filename,
MainWindow main_window,
QString *  directory,
QString  name,
QString  x_axis_description,
QString  y_axis_description,
bool  swap_spatial,
std::string  format 
)

VespucciDataset::VespucciDataset.

Parameters
text_filenameThe filename of the text file from which data is imported
main_windowThe main window of the program
directoryThe current global working directory
nameThe name of the dataset displayed to the user
x_axis_descriptionA description of the spectral abscissa
y_axis_descriptionA description of the spectral ordinate
swap_spatialWhether or not y is located in the first column instead of the second (some Horiba spectrometers do this). Main function for processing data from text files to create VespucciDataset objects. Currently written to accept files in "wide" format, will be expanded to deal with different ASCII formats later with conditionals.

Definition at line 360 of file vespuccidataset.cpp.

VespucciDataset::VespucciDataset ( map< pair< int, int >, string >  text_filenames,
MainWindow main_window,
QString *  directory,
QString  name,
QString  x_axis_description,
QString  y_axis_description,
int  rows,
int  cols 
)

Definition at line 444 of file vespuccidataset.cpp.

VespucciDataset::VespucciDataset ( QString  name,
MainWindow main_window,
QString *  directory,
QSharedPointer< VespucciDataset original,
uvec  indices 
)

VespucciDataset::VespucciDataset.

Parameters
nameThe name of the dataset displayed to the user
main_windowThe main window of the program
directoryThe current global working directory
originalThe dataset from which this dataset is "extracted"
indicesThe indices of the spectra in the previous dataset that will form this one This is a constructor to create a new dataset by "extracting" spectra from a another dataset based on values on an image.

Definition at line 494 of file vespuccidataset.cpp.

VespucciDataset::VespucciDataset ( QString  name,
MainWindow main_window,
QString *  directory 
)

VespucciDataset::VespucciDataset.

Parameters
nameDataset name
main_windowThe
directoryConstructor to create a dataset without spatial and spectral data set (i.e. when using MetaDataset or a factory method). Use SetData() to instantiate dataset members

Definition at line 544 of file vespuccidataset.cpp.

VespucciDataset::~VespucciDataset ( )

Member Function Documentation

vec VespucciDataset::abscissa ( ) const
virtual

Implements AbstractDataAnalyzer.

Definition at line 2266 of file vespuccidataset.cpp.

mat * VespucciDataset::abscissa_ptr ( )

VespucciDataset::abscissa_ptr.

Returns

Definition at line 2749 of file vespuccidataset.cpp.

const vec & VespucciDataset::abscissa_ref ( )

VespucciDataset::abscissa_ref.

Returns
abscissa_ as reference

Definition at line 2675 of file vespuccidataset.cpp.

double VespucciDataset::AbscissaMax ( ) const
virtual

Implements AbstractDataAnalyzer.

Definition at line 2708 of file vespuccidataset.cpp.

double VespucciDataset::AbscissaMin ( ) const
virtual

Implements AbstractDataAnalyzer.

Definition at line 2703 of file vespuccidataset.cpp.

void VespucciDataset::AbsoluteValue ( )

Definition at line 919 of file vespuccidataset.cpp.

void VespucciDataset::AddAnalysisResult ( QSharedPointer< AnalysisResults analysis_result)

Definition at line 2819 of file vespuccidataset.cpp.

void VespucciDataset::AddAnalysisResult ( QSharedPointer< AnalysisResults analysis_result,
uword  start_row,
uword  end_row 
)

Definition at line 2833 of file vespuccidataset.cpp.

void VespucciDataset::AddAuxiliaryMatrix ( const QString &  name,
mat &  matrix 
)

VespucciDataset::AddAuxiliaryMatrix.

Parameters
nameDisplay name and internal key of the matrix to add
matrixAuxiliary matrix to insert into this dataset

Definition at line 2892 of file vespuccidataset.cpp.

void VespucciDataset::AddMap ( QSharedPointer< MapData map)

VespucciDataset::AddMap Adds a map to the list of map pointers and adds its name to relevant lists.

Parameters
map

Definition at line 2457 of file vespuccidataset.cpp.

void VespucciDataset::AddMatrix ( const QString &  name,
mat &  matrix 
)

Definition at line 2904 of file vespuccidataset.cpp.

void VespucciDataset::AgglomerativeClustering ( const QString &  name,
const QString &  linkage,
const QString &  metric 
)
virtual

VespucciDataset::AgglomerativeClustering.

Parameters
name
linkage
metriclinkage and metric must be valid arguments for Vespucci::Math::Clustering::AHCA

Implements AbstractDataAnalyzer.

Definition at line 1842 of file vespuccidataset.cpp.

QStringList VespucciDataset::AnalysisResultsKeys ( ) const

VespucciDataset::AnalysisResultsKeys.

Returns
Get a list of keys for analysis results. Used to update display model.

Definition at line 2852 of file vespuccidataset.cpp.

QMap< QString, QStringList > VespucciDataset::AnalysisResultsTreeStructure ( ) const

VespucciDataset::AnalysisResultsTreeStructure.

Returns
Container parsed by TreeModel

Definition at line 2864 of file vespuccidataset.cpp.

void VespucciDataset::ApplyFTWeight ( QString  type,
double  param 
)

VespucciDataset::ApplyFTWeight.

Parameters
type
param

Definition at line 1644 of file vespuccidataset.cpp.

void VespucciDataset::ApplyFTWeight ( double  start_offset,
double  end_offset,
double  power 
)

VespucciDataset::ApplyFTWeight.

Parameters
type
start_offset
end_offset
power

Definition at line 1666 of file vespuccidataset.cpp.

QStringList VespucciDataset::AuxiliaryMatrixKeys ( ) const

VespucciDataset::AuxiliaryMatrixKeys.

Returns
List of names of imported auxiliary matrices Used in updating display model

Definition at line 2924 of file vespuccidataset.cpp.

mat VespucciDataset::AverageSpectrum ( bool  stats) const

VespucciDataset::AverageSpectrum Finds the average of the spectrum. This can be saved by the user. Probably not all that useful, except for determining a spectrum to use as a background spectrum for other maps.

Parameters
statsWhether or not to include standard deviations on the second row.
Returns
The average spectrum

Definition at line 2600 of file vespuccidataset.cpp.

void VespucciDataset::BandRatio ( const QString &  name,
double &  first_left_bound,
double &  first_right_bound,
double &  second_left_bound,
double &  second_right_bound,
uword  bound_window 
)

VespucciDataset::BandRatio.

Parameters
first_left_bound
first_right_bound
second_left_bound
second_right_bound
bound_window

Definition at line 1735 of file vespuccidataset.cpp.

void VespucciDataset::Booleanize ( double  min,
double  max,
bool  keep_inside,
bool  oneify 
)

VespucciDataset::Booleanize.

Parameters
min
max
keep_inside
oneify

Definition at line 753 of file vespuccidataset.cpp.

void VespucciDataset::CalculateRepresentativeSpectrum ( const QString &  name,
QString  statistic,
QString  metric 
)

VespucciDataset::CalculateRepresentativeSpectrum.

Parameters
name
statisticEither "centroid" or "medoid"
metricValid for DistanceMetricWrapper statist

Definition at line 1872 of file vespuccidataset.cpp.

void VespucciDataset::Clamp ( double  min,
double  max 
)

VespucciDataset::Clamp.

Parameters
min
max

Definition at line 787 of file vespuccidataset.cpp.

void VespucciDataset::ClassicalLeastSquares ( const QString &  name,
const QStringList &  reference_keys 
)
virtual

Implements AbstractDataAnalyzer.

Definition at line 2066 of file vespuccidataset.cpp.

size_t VespucciDataset::columns ( ) const
virtual

Implements AbstractDataAnalyzer.

Definition at line 2698 of file vespuccidataset.cpp.

bool VespucciDataset::ConstructorCancelled ( ) const

VespucciDataset::ConstructorCancelled Specifies whether or not the constructor has been canceled. The constructor asks this and cleans everything up in case it is canceled.

Returns

Definition at line 2587 of file vespuccidataset.cpp.

bool VespucciDataset::Contains ( const QString &  key) const

VespucciDataset::Contains.

Parameters
key
Returns
Every AnalysisResult and auxiliary matrix must have a unique name. This is called before an object is added to verify that it doesn't already exist The caller will rename the object if it exists by appending a number that is iteratively increased until the name is unique

Definition at line 156 of file vespuccidataset.cpp.

QStringList VespucciDataset::CoreMatrixKeys ( ) const

VespucciDataset::CoreMatrixKeys.

Returns

Definition at line 2933 of file vespuccidataset.cpp.

void VespucciDataset::CorrelationAnalysis ( const QString &  control_key,
QString  name 
)

VespucciDataset::CorrelationAnalysis.

Parameters
control_key
nameObject with key control_key must be column vector or else only first column used

Definition at line 2004 of file vespuccidataset.cpp.

void VespucciDataset::CreateMap ( const QString &  map_name,
const QString &  results_key,
const QString &  matrix_key,
uword  column,
QCPColorGradient  gradient 
)

VespucciDataset::CreateMap.

Parameters
map_name
results_key
matrix_key
column
gradient
tick_countCaller should check the number of unique elements in the object before plotting it A cluster color gradient should be used for 8 or less unique elements.

Definition at line 3013 of file vespuccidataset.cpp.

void VespucciDataset::CreateMap ( const QString &  map_name,
const QString &  matrix_key,
uword  column,
QCPColorGradient  gradient 
)

Definition at line 3035 of file vespuccidataset.cpp.

void VespucciDataset::CropSpectra ( double  x_min,
double  x_max,
double  y_min,
double  y_max,
double  wl_min,
double  wl_max 
)

VespucciDataset::CropSpectra Crops spectra_ based on.

Parameters
x_minvalue of x below which spectra are deleted
x_maxvalue of x above which spectra are deleted
y_minvalue of y below which spectra are deleted
y_maxvalue of y above which spectra are deleted Removes all data points outside of the range. Cannot be undone. It is preferable to create a new dataset rather than crop an old one.

Definition at line 613 of file vespuccidataset.cpp.

void VespucciDataset::CWTBaseline ( int  lambda,
int  penalty_order,
double  SNR_threshold,
double  peak_shape_threshold 
)

Definition at line 1024 of file vespuccidataset.cpp.

cx_mat VespucciDataset::cx_spectra ( ) const

Definition at line 2373 of file vespuccidataset.cpp.

cx_mat VespucciDataset::cx_spectra ( uvec  indices) const

Definition at line 2381 of file vespuccidataset.cpp.

QString VespucciDataset::filename ( ) const

Definition at line 179 of file vespuccidataset.cpp.

int VespucciDataset::FindIndex ( int  x_value,
int  y_value 
) const
uword VespucciDataset::FindIndex ( double  abscissa_value) const

Definition at line 2145 of file vespuccidataset.cpp.

uword VespucciDataset::FindOrigin ( ) const

VespucciDataset::FindOrigin.

Returns
Find the point closest to (0,0) in the map. If the function fails, find the index in the middle of the spatial data.

Definition at line 2113 of file vespuccidataset.cpp.

void VespucciDataset::FindPeaks ( const QString &  name,
double  sel,
double  threshold,
uword  poly_order,
uword  window_size 
)

VespucciDataset::FindPeaks.

Parameters
selThe amount above the surrounding data for a peak to be identified
thresholdA value which peaks must be larger than to be maxima
poly_orderPolynomial order for Savitzky-Golay derivatization
window_sizeWindow size for Savitzky-Golay derivatization

Definition at line 1812 of file vespuccidataset.cpp.

uvec VespucciDataset::FindRange ( double  start,
double  end 
) const

VespucciDataset::FindRange. Finds the index of the wavelength value closest to the specified wavelength range.

Parameters
startthe first wavelength in the spectral region of interest
endthe second wavelength in the spectral region of interest
Returns

Definition at line 2100 of file vespuccidataset.cpp.

void VespucciDataset::FitPeak ( const QString &  name,
const QString &  peak_shape,
double &  left_bound,
double &  right_bound 
)
virtual

Implements AbstractDataAnalyzer.

Definition at line 1714 of file vespuccidataset.cpp.

void VespucciDataset::FourierTransform ( int  n)

VespucciDataset::FourierTransform.

Parameters
n

Definition at line 1597 of file vespuccidataset.cpp.

QSharedPointer< AnalysisResults > VespucciDataset::GetAnalysisResult ( const QString &  key)

Definition at line 2956 of file vespuccidataset.cpp.

const mat & VespucciDataset::GetAnalysisResultMatrix ( const QString &  results_key,
const QString &  matrix_key 
) const

VespucciDataset::GetAnalysisResultMatrix.

Parameters
results_key
matrix_key
Returns
Currently iterates over all results (O(n)) could change structure to use binary tree (O(log n)) but I don't think it should ever be necessary.

Definition at line 2945 of file vespuccidataset.cpp.

const mat & VespucciDataset::GetAuxiliaryMatrix ( const QString &  key) const

Definition at line 2968 of file vespuccidataset.cpp.

QCPColorGradient VespucciDataset::GetClusterGradient ( int  clusters) const

VespucciDataset::GetClusterGradient Cluster gradients are slightly different from the continuous gradients. This selects the right gradient based on the number of clusters.

Parameters
clustersNumber of clusters
Returns
Proper color gradient for number of clusters

Definition at line 2565 of file vespuccidataset.cpp.

const mat & VespucciDataset::GetCoreMatrix ( const QString &  key) const

Definition at line 2975 of file vespuccidataset.cpp.

QCPColorGradient VespucciDataset::GetGradient ( int  gradient_number) const

VespucciDataset::GetGradient Selects the color gradient from list of presets.

Parameters
gradient_number
Returns

Definition at line 2511 of file vespuccidataset.cpp.

QSharedPointer< MapData > VespucciDataset::GetMapData ( const QString &  key)

Definition at line 2995 of file vespuccidataset.cpp.

int VespucciDataset::HySime ( )

VespucciDataset::HySime.

Returns
Dimensionality predicted by HySime algorithm

Definition at line 1407 of file vespuccidataset.cpp.

void VespucciDataset::IModPolyBaseline ( const uword  poly_order,
const uword  max_it,
double  threshold 
)

Definition at line 1029 of file vespuccidataset.cpp.

void VespucciDataset::ImportAuxiliaryMatrix ( const QString &  name,
const QString &  filename 
)

VespucciDataset::ImportAuxiliaryMatrix.

Parameters
name
filenameLoad a file into the auxiliary matrix map

Definition at line 2877 of file vespuccidataset.cpp.

vec VespucciDataset::indices ( ) const

VespucciDataset::indices.

Returns
The indices_ vector

Definition at line 2289 of file vespuccidataset.cpp.

mat * VespucciDataset::indices_ptr ( )

VespucciDataset::indices_ptr.

Returns
A pointer to the indices_ vector

Definition at line 2298 of file vespuccidataset.cpp.

void VespucciDataset::InterpolateToNewAbscissa ( const vec &  new_abscissa,
unsigned int  poly_order,
unsigned int  window_size 
)

VespucciDataset::InterpolateToNewAbscissa.

Parameters
new_abscissa
polynomial_order
window_sizeSpline interpolation to new abscissa

Definition at line 1555 of file vespuccidataset.cpp.

void VespucciDataset::InterpolateToNewAbscissa ( const vec &  new_abscissa)

VespucciDataset::InterpolateToNewAbscissa.

Parameters
new_abscissaLinear interpolation to new abscissa

Definition at line 1576 of file vespuccidataset.cpp.

void VespucciDataset::InverseFourierTransform ( int  n)

Definition at line 1617 of file vespuccidataset.cpp.

bool VespucciDataset::IsCoreMatrix ( const QString &  key) const

VespucciDataset::IsCoreMatrix.

Parameters
key
Returns
Check whether the DatasetTreeItem child is core matrix or not DataViewer needs to know whether to use GetCoreMatrix or GetAuxiliaryMatrix after action is performed on a MatrixTreeItem that is child of DatasetTreeItem

Definition at line 2990 of file vespuccidataset.cpp.

bool VespucciDataset::IsValid ( ) const

Definition at line 162 of file vespuccidataset.cpp.

QCPRange VespucciDataset::KeyRange ( ) const

VespucciDataset::KeyRange Finds the minima and maxima of x variable to properly set axes of QCustomPlot objects.

Returns

Definition at line 2171 of file vespuccidataset.cpp.

int VespucciDataset::KeySize ( ) const

VespucciDataset::KeySize Finds the number of data points in x variable to properly set axes of QCustomPlot objects.

Returns
number of unique x values

Definition at line 2185 of file vespuccidataset.cpp.

void VespucciDataset::KMeans ( const QString &  name,
const QString &  metric_text,
const QString &  partition_policy,
bool  allow_empty,
size_t  clusters 
)
virtual

VespucciDataset::KMeans.

Implements K-means clustering using MLPACK

Parameters
clustersNumber of clusters to find
metricDistance metric
nameName of map in workspace.

Implements AbstractDataAnalyzer.

Definition at line 2041 of file vespuccidataset.cpp.

const QString VespucciDataset::last_operation ( ) const

VespucciDataset::last_operation.

Returns
Description of last pre-processing operation performed

Definition at line 2813 of file vespuccidataset.cpp.

QString VespucciDataset::last_save_filename ( ) const

Definition at line 184 of file vespuccidataset.cpp.

void VespucciDataset::LinearMovingAverage ( unsigned int  window_size)

VespucciDataset::LinearMovingAverage Performs moving average filtering on the spectral data. Entries near the boundaries of spectra are not processed. See also VespucciDataset::MedianFilter.

Parameters
window_size- an odd number representing the width of the window.

Definition at line 1225 of file vespuccidataset.cpp.

bool VespucciDataset::Load ( QString  filename)
MapListModel* VespucciDataset::map_list_model ( )
int VespucciDataset::map_loading_count ( ) const

VespucciDataset::map_names.

Returns
list of names of the maps created from this dataset VespucciDataset::map_loading_count
number of maps created for this dataset

Definition at line 2446 of file vespuccidataset.cpp.

QStringList VespucciDataset::MapKeys ( ) const

Definition at line 3067 of file vespuccidataset.cpp.

void VespucciDataset::MeanCenter ( )

VespucciDataset::MeanCenter Subtract the mean of all spectra in the dataset at each wavelength from each wavelength (i.e. in preparation for PCA)

Definition at line 712 of file vespuccidataset.cpp.

void VespucciDataset::MedianFilter ( unsigned int  window_size)

VespucciDataset::MedianFilter performs median filtering on the spectral data. Entries near the boundaries of spectra are not processed. See also VespucciDataset::LinearMovingAverage.

Parameters
window_size- an odd number representing the width of the window.

Definition at line 1201 of file vespuccidataset.cpp.

bool VespucciDataset::meta ( ) const

VespucciDataset::meta.

Returns
Whether or not this is an instance of MetaDataset

Definition at line 2790 of file vespuccidataset.cpp.

void VespucciDataset::MFBaseline ( int  window_size,
int  iterations 
)

VespucciDataset::MFBaseline Baseline-adjusts the data. This function uses a median filter with a large window to determine the baseline on the assumption that the median value is more likely to be basline than spectrum. This will complicate things if there are many peaks. Additionally, this significantly narrows the shape of peaks. Other baseline methods may be implemented later.

Parameters
method
window_size

Definition at line 977 of file vespuccidataset.cpp.

void VespucciDataset::MinMaxNormalize ( )

VespucciDataset::MinMaxNormalize Normalizes each column of spectra_ so that the highest value is 1 and lowest value is 0.

Definition at line 663 of file vespuccidataset.cpp.

bool VespucciDataset::mlpack_pca_calculated ( ) const
const QString VespucciDataset::name ( ) const

VespucciDataset::name.

Returns
member name_, the name of the dataset as seen by the user

Definition at line 2403 of file vespuccidataset.cpp.

bool VespucciDataset::non_spatial ( ) const

VespucciDataset::non_spatial.

Returns
True if map has empty x_ and y_

Definition at line 2781 of file vespuccidataset.cpp.

QStringList VespucciDataset::operations ( )

Definition at line 2303 of file vespuccidataset.cpp.

mat * VespucciDataset::parent_dataset_indices ( )

Definition at line 2804 of file vespuccidataset.cpp.

bool VespucciDataset::partial_least_squares_calculated ( ) const
void VespucciDataset::PartialLeastSquares ( const QString &  name,
uword  components 
)

VespucciDataset::PartialLeastSquares.

Parameters
componentsPerform partial least squares without creating an image

Definition at line 1940 of file vespuccidataset.cpp.

void VespucciDataset::PeakIntensityNormalize ( double  peak_position)
void VespucciDataset::PeakIntensityNormalize ( double  left_bound,
double  right_bound 
)

VespucciDataset::PeakIntensityNormalize.

Parameters
left_bound
right_bound

Definition at line 732 of file vespuccidataset.cpp.

void VespucciDataset::PLSCalibration ( const QString &  name,
const QStringList &  control_keys 
)
virtual

Implements AbstractDataAnalyzer.

Definition at line 1965 of file vespuccidataset.cpp.

vec VespucciDataset::PointSpectrum ( uword  index) const

VespucciDataset::PointSpectrum.

Parameters
index
Returns

Definition at line 2129 of file vespuccidataset.cpp.

QCPRange VespucciDataset::PointSpectrumRange ( int  i) const

VespucciDataset::PointSpectrumRange.

Parameters
icol of spectra_ containing desired spectrum
Returns
the range of y values for the point spectra at i

Definition at line 2495 of file vespuccidataset.cpp.

bool VespucciDataset::principal_components_calculated ( ) const
void VespucciDataset::PrincipalComponents ( const QString &  name)
virtual

VespucciDataset::PrincipalComponents Perform PCA without creating an image.

Implements AbstractDataAnalyzer.

Definition at line 1787 of file vespuccidataset.cpp.

void VespucciDataset::PrincipalComponents ( const QString &  name,
bool  scale_data 
)
virtual

VespucciDataset::PrincipalComponents.

Parameters
scaleData
recalculatePerform PCA using MLPACK's implementation. This is designed for dimensionality reduction, not classification, so there is no corresponding imaging function

Implements AbstractDataAnalyzer.

Definition at line 1768 of file vespuccidataset.cpp.

int VespucciDataset::QUIC_SVD ( double  epsilon)

VespucciDataset::QUIC_SVD.

Parameters
epsilonError tolerance fraction for calculated subspace Use the QUIC-SVD algorithm to denoise the spectra by finding a lower-rank approximation of the matrix.

Definition at line 1285 of file vespuccidataset.cpp.

void VespucciDataset::RemoveClippedSpectra ( double  threshold)

VespucciDataset::RemoveClippedSpectra.

Parameters
thresholdRemoves spectra with a maximum value at or above the threshold

Definition at line 1078 of file vespuccidataset.cpp.

void VespucciDataset::RemoveFlatSpectra ( double  threshold)

VespucciDataset::RemoveFlatSpectra.

Parameters
thresholdRemoves "flat" spectra whose maximum is less than the threshold

Definition at line 1109 of file vespuccidataset.cpp.

void VespucciDataset::RemoveMap ( const QString &  name)

Definition at line 2467 of file vespuccidataset.cpp.

void VespucciDataset::RollingBallBaseline ( size_t  wm,
size_t  ws 
)

VespucciDataset::RollingBallBaseline.

Parameters
wm
ws

Definition at line 1004 of file vespuccidataset.cpp.

bool VespucciDataset::Save ( QString  filename)

Definition at line 42 of file vespuccidataset.cpp.

bool VespucciDataset::saved ( ) const

Definition at line 189 of file vespuccidataset.cpp.

bool VespucciDataset::SaveSpectrum ( QString  filename,
uword  column,
file_type  type 
)

Definition at line 123 of file vespuccidataset.cpp.

void VespucciDataset::SavitzkyGolay ( unsigned int  derivative_order,
unsigned int  polynomial_order,
unsigned int  window_size 
)

VespucciDataset::SavitzkyGolay Performs derivatization/Savitzky-Golay smoothing.

Parameters
derivative_orderThe order of the derivative.
polynomial_orderThe order of the polynomial
window_sizeThe size of the filter window.

Definition at line 1320 of file vespuccidataset.cpp.

void VespucciDataset::Scale ( double  scaling_factor)

Definition at line 1346 of file vespuccidataset.cpp.

void VespucciDataset::SetData ( const mat &  spectra,
const vec &  wavelength,
const vec &  x,
const vec &  y 
)

VespucciDataset::SetData.

Parameters
spectraSpectra
wavelengthSpectral abscissa
xvec of horizontal axis spatial data
yvec of vertical axis spatial data Set the data of the dataset. Used by the MetaDataset constructor

Definition at line 2424 of file vespuccidataset.cpp.

void VespucciDataset::SetIndices ( vec  indices)

VespucciDataset::SetIndices.

Parameters
indicesA new indices_ vector

Definition at line 2312 of file vespuccidataset.cpp.

void VespucciDataset::SetName ( QString  new_name)

VespucciDataset::SetName.

Parameters
new_namenew name of dataset

Definition at line 2412 of file vespuccidataset.cpp.

void VespucciDataset::SetOldCopies ( )

Definition at line 140 of file vespuccidataset.cpp.

void VespucciDataset::SetParentDatasetIndices ( mat  parent_dataset_indices)

VespucciDataset::SetParentDatasetIndices.

Parameters
parent_dataset_indicesSet the parent_dataset_indices_ variable. Used by the MetaDataset constructor.

Definition at line 2799 of file vespuccidataset.cpp.

void VespucciDataset::SetXDescription ( QString  description)

VespucciDataset::SetXDescription Sets the value of the spectral abscissa description.

Parameters
description

Definition at line 2629 of file vespuccidataset.cpp.

void VespucciDataset::SetYDescription ( QString  description)

VespucciDataset::SetYDescription.

Parameters
descriptionSets the value of the spectral ordinate axis description

Definition at line 2638 of file vespuccidataset.cpp.

void VespucciDataset::ShedSpectrum ( const uword  index)

VespucciDataset::ShedSpectrum.

Parameters
indexIndex of spectrum to shed. Removes a spectrum from the dataset

Definition at line 1366 of file vespuccidataset.cpp.

void VespucciDataset::ShedZeroSpectra ( )

Definition at line 807 of file vespuccidataset.cpp.

void VespucciDataset::ShedZeroWavelengths ( )

Definition at line 832 of file vespuccidataset.cpp.

bool VespucciDataset::ShowMapViewer ( const QString &  map_key,
bool  show 
)

Definition at line 3053 of file vespuccidataset.cpp.

void VespucciDataset::SingularValue ( unsigned int  singular_values)

VespucciDataset::SingularValue Denoises the spectra matrix using a truncated singular value decomposition. The first singular_values singular values are used to "reconstruct" the spectra matrix. The function used to find the truncated SVD is Vespucci::Math::DimensionReduction::svds.

Parameters
singular_valuesNumber of singular values to use.

Definition at line 1252 of file vespuccidataset.cpp.

void VespucciDataset::SNVNormalize ( double  offset,
bool  center 
)

Definition at line 902 of file vespuccidataset.cpp.

mat VespucciDataset::spectra ( uvec  indices) const

VespucciDataset::spectra.

Parameters
indicesVector of indices
Returns
Submat of spectra at indices

Definition at line 2394 of file vespuccidataset.cpp.

mat VespucciDataset::spectra ( ) const

VespucciDataset::spectra.

Returns
member spectra_

Definition at line 2368 of file vespuccidataset.cpp.

mat* VespucciDataset::spectra_old_ptr ( bool *  ok)
mat * VespucciDataset::spectra_ptr ( )

VespucciDataset::spectra_ptr.

Returns

Definition at line 2657 of file vespuccidataset.cpp.

const mat & VespucciDataset::spectra_ref ( )

VespucciDataset::spectra_ref.

Returns
spectra_ as reference

Definition at line 2666 of file vespuccidataset.cpp.

bool VespucciDataset::state_changed ( ) const

Definition at line 174 of file vespuccidataset.cpp.

void VespucciDataset::SubtractBackground ( const QStringList &  data_keys)

VespucciDataset::SubtractBackground.

Parameters
data_keyslocation of matrix in workspace

Definition at line 942 of file vespuccidataset.cpp.

void VespucciDataset::TrainPLSDA ( const QString &  name,
const QStringList &  label_keys 
)
virtual

Implements AbstractDataAnalyzer.

Definition at line 1982 of file vespuccidataset.cpp.

void VespucciDataset::TransformAbscissa ( QString  input_units,
double  input_factor,
QString  output_units,
double  output_factor,
QString  description 
)

VespucciDataset::TransformAbscissa.

Parameters
input_units
input_factor
output_units
output_factor

Definition at line 1435 of file vespuccidataset.cpp.

void VespucciDataset::Undo ( )

VespucciDataset::Undo Swap spectra_ and spectra_old_ to undo an action. Calling this function again re-does the action that was undid.

Definition at line 568 of file vespuccidataset.cpp.

bool VespucciDataset::Undoable ( ) const

VespucciDataset::Undoable.

Returns
Whether or not the last operation can be undone

Definition at line 2718 of file vespuccidataset.cpp.

uword VespucciDataset::UniqueX ( ) const

Definition at line 2723 of file vespuccidataset.cpp.

uword VespucciDataset::UniqueY ( ) const

Definition at line 2729 of file vespuccidataset.cpp.

void VespucciDataset::Univariate ( const QString &  name,
double &  left_bound,
double &  right_bound,
uword  bound_window 
)

VespucciDataset::Univariate.

Parameters
left_bound
right_bound
bound_window

Definition at line 1691 of file vespuccidataset.cpp.

QList<QSharedPointer<UnivariateData> > VespucciDataset::univariate_datas ( )
int VespucciDataset::UnivariateCount ( ) const

VespucciDataset::UnivariateCount.

Returns
Number of univariate/band ratio data objects have been created

Definition at line 2740 of file vespuccidataset.cpp.

QCPRange VespucciDataset::ValueRange ( ) const

VespucciDataset::ValueRange Finds the minima and maxima of y variable to properly set axes of QCustomPlot objects.

Returns

Definition at line 2158 of file vespuccidataset.cpp.

int VespucciDataset::ValueSize ( ) const

VespucciDataset::ValueSize Finds number of unique y values for properly setting QCPAxis.

Returns
number of unique y values

Definition at line 2220 of file vespuccidataset.cpp.

void VespucciDataset::VectorNormalize ( uword  norm)

VespucciDataset::VectorNormalize.

Parameters
normIf 1, unit area, if 2 unit length

Definition at line 691 of file vespuccidataset.cpp.

bool VespucciDataset::vertex_components_calculated ( ) const
void VespucciDataset::VertexComponents ( const QString &  name,
uword  endmembers 
)

VespucciDataset::VertexComponents.

Parameters
endmembersPerform VCA without creating an image

Definition at line 1895 of file vespuccidataset.cpp.

vec VespucciDataset::wavelength ( ) const

VespucciDataset::wavelength.

Returns
member abscissa_ (spectrum key values)

Definition at line 2261 of file vespuccidataset.cpp.

vec VespucciDataset::wavelength ( uvec  indices) const

Definition at line 2271 of file vespuccidataset.cpp.

mat * VespucciDataset::wavelength_ptr ( )

Definition at line 2754 of file vespuccidataset.cpp.

QVector< double > VespucciDataset::WavelengthQVector ( ) const

Definition at line 2134 of file vespuccidataset.cpp.

QCPRange VespucciDataset::WavelengthRange ( ) const

VespucciDataset::WavelengthRange.

Returns
the range of the wavlength vector (for plotting point spectra)

Definition at line 2482 of file vespuccidataset.cpp.

vec VespucciDataset::x ( ) const

VespucciDataset::x.

Returns
member x_

Definition at line 2280 of file vespuccidataset.cpp.

vec VespucciDataset::x ( uvec  indices) const

VespucciDataset::x.

Parameters
indicesVector of indices
Returns
Subvec of x corresponding to valeus in indices

Definition at line 2322 of file vespuccidataset.cpp.

double VespucciDataset::x ( uword  index) const

Definition at line 2327 of file vespuccidataset.cpp.

const QString VespucciDataset::x_axis_description ( ) const

VespucciDataset::x_axis_description The x_axis_description is printed on the spectrum viewer.

Returns
Spectral abscissa description.

Definition at line 2619 of file vespuccidataset.cpp.

mat * VespucciDataset::x_ptr ( )

VespucciDataset::x_ptr.

Returns

Definition at line 2763 of file vespuccidataset.cpp.

const vec & VespucciDataset::x_ref ( )

VespucciDataset::x_ref.

Returns
x_ as reference

Definition at line 2684 of file vespuccidataset.cpp.

vec VespucciDataset::y ( ) const

VespucciDataset::y.

Returns
member y_

Definition at line 2341 of file vespuccidataset.cpp.

vec VespucciDataset::y ( uvec  indices) const

VespucciDataset::y.

Parameters
indicesVector of indices
Returns
Subvec of y at indices

Definition at line 2351 of file vespuccidataset.cpp.

double VespucciDataset::y ( uword  index) const

Definition at line 2356 of file vespuccidataset.cpp.

const QString VespucciDataset::y_axis_description ( ) const

VespucciDataset::y_axis_description.

Returns
The spectral ordinate axis description.

Definition at line 2647 of file vespuccidataset.cpp.

mat * VespucciDataset::y_ptr ( )

VespucciDataset::y_ptr.

Returns

Definition at line 2772 of file vespuccidataset.cpp.

const vec & VespucciDataset::y_ref ( )

VespucciDataset::y_ref.

Returns
y_ as reference

Definition at line 2693 of file vespuccidataset.cpp.

void VespucciDataset::ZeroClippedSpectra ( double  threshold)

VespucciDataset::ZeroClippedSpectra.

Parameters
threshold

Definition at line 1140 of file vespuccidataset.cpp.

void VespucciDataset::ZeroFlatSpectra ( double  threshold)

VespucciDataset::ZeroFlatSpectra.

Parameters
threshold

Definition at line 1169 of file vespuccidataset.cpp.

void VespucciDataset::ZeroSpectrum ( const uword  index)

Definition at line 1387 of file vespuccidataset.cpp.

void VespucciDataset::ZScoreNormalize ( )

VespucciDataset::ZScoreNormalize Normalizes each row using the score of the normal distribution.

Definition at line 881 of file vespuccidataset.cpp.

mat VespucciDataset::ZScoreNormCopy ( )

VespucciDataset::ZScoreNormCopy For when you want to Z-score normalize without changing spectra_.

Returns
A normalized copy of the matrix.

Definition at line 861 of file vespuccidataset.cpp.

Member Data Documentation

QTextStream VespucciDataset::log_text_stream_

log_text_stream_ A log of all UI-induced functions called on this dataset

Definition at line 301 of file vespuccidataset.h.


The documentation for this class was generated from the following files: