MINT2
|
#include <HyperPointSet.h>
Public Member Functions | |
const int & | getDimension () const |
bool | compatible (const HyperPoint &other, bool printError=true) const |
const HyperPoint & | at (int i) const |
HyperPoint & | at (int i) |
unsigned int | size () const |
void | push_back (const HyperPoint &point) |
void | addHyperPointSet (const HyperPointSet &other) |
void | removeDuplicates () |
void | sort () |
double | getSumW () const |
double | getSumW2 () const |
double | getCorrelation (int i, int j) const |
TMatrixD | getCorrelationMatrix () const |
double | getCovarience (int i, int j) const |
TMatrixD | getCovarienceMatrix () const |
HyperPoint | getMin () const |
Get a HyperCuboid that surrounds the points. More... | |
HyperPoint | getMax () const |
Get a HyperPoint that gives the maximum in each dim. More... | |
HyperPoint | mean () const |
HyperPoint | geometricMean () const |
HyperPoint | harmonicMean () const |
void | save (TString path) |
void | load (TString path) |
void | print (std::ostream &os=std::cout) const |
HyperPointSet (int dimension) | |
HyperPointSet (const HyperPoint &point) | |
HyperPointSet (const HyperPoint &point1, const HyperPoint &point2) | |
HyperPointSet (const HyperPoint &point1, const HyperPoint &point2, const HyperPoint &point3) | |
HyperPointSet (int npoints, const HyperPoint &point) | |
HyperPointSet (TString path) | |
virtual | ~HyperPointSet () |
bool | linearlyIndependant () const |
HyperPointSet | gramSchmidtProcess () const |
void | normalise () |
Private Member Functions | |
void | save () |
Private Attributes | |
int | _dimension |
std::vector< HyperPoint > | _points |
HyperPlot, Author: Sam Harnew, sam.h , Date: Dec 2015 arne w@gma il.c om
A vector of HyperPoints
Definition at line 31 of file HyperPointSet.h.
HyperPointSet::HyperPointSet | ( | int | dimension | ) |
Standard constuctor where only the dimensionality is given
Definition at line 6 of file HyperPointSet.cpp.
HyperPointSet::HyperPointSet | ( | const HyperPoint & | point | ) |
Constuctor for a HyperPointSet that takes dimensionality from a given HyperPoint and adds that point to the HyperPointSet
Definition at line 35 of file HyperPointSet.cpp.
HyperPointSet::HyperPointSet | ( | const HyperPoint & | point1, |
const HyperPoint & | point2 | ||
) |
Constuctor for a HyperPointSet that takes dimensionality from a given HyperPoint and adds that point and another to the HyperPointSet
Definition at line 43 of file HyperPointSet.cpp.
HyperPointSet::HyperPointSet | ( | const HyperPoint & | point1, |
const HyperPoint & | point2, | ||
const HyperPoint & | point3 | ||
) |
Constuctor for a HyperPointSet that takes dimensionality from a given HyperPoint and adds that point and another 2 to the HyperPointSet
Definition at line 52 of file HyperPointSet.cpp.
HyperPointSet::HyperPointSet | ( | int | npoints, |
const HyperPoint & | point | ||
) |
Constuctor for a HyperPointSet that repeats the same point a total of npoints times.
Definition at line 24 of file HyperPointSet.cpp.
HyperPointSet::HyperPointSet | ( | TString | path | ) |
Constuctor for loading exisiting HyperPointSet from a file
Definition at line 14 of file HyperPointSet.cpp.
|
virtual |
void HyperPointSet::addHyperPointSet | ( | const HyperPointSet & | other | ) |
Add all the points from another HyperPointSet to this one.
Definition at line 202 of file HyperPointSet.cpp.
const HyperPoint & HyperPointSet::at | ( | int | i | ) | const |
HyperPoint & HyperPointSet::at | ( | int | i | ) |
bool HyperPointSet::compatible | ( | const HyperPoint & | other, |
bool | printError = true |
||
) | const |
Is this HyperPoint compatible with this HyperPointSet? (same dimension)
Definition at line 62 of file HyperPointSet.cpp.
HyperPoint HyperPointSet::geometricMean | ( | ) | const |
Get the geometric mean of all points. This does not use the weights.
Definition at line 460 of file HyperPointSet.cpp.
double HyperPointSet::getCorrelation | ( | int | i, |
int | j | ||
) | const |
Get the correlation between two variables.
Definition at line 217 of file HyperPointSet.cpp.
TMatrixD HyperPointSet::getCorrelationMatrix | ( | ) | const |
Get the correlation matrix associated with the HyperPointSet.
Definition at line 286 of file HyperPointSet.cpp.
double HyperPointSet::getCovarience | ( | int | i, |
int | j | ||
) | const |
Get the covarience between two variables.
Definition at line 245 of file HyperPointSet.cpp.
TMatrixD HyperPointSet::getCovarienceMatrix | ( | ) | const |
Get the covarience matrix associated with the HyperPointSet.
Definition at line 271 of file HyperPointSet.cpp.
|
inline |
Get the dimensionality of HyperPoints in the HyperPointSet
Definition at line 40 of file HyperPointSet.h.
HyperPoint HyperPointSet::getMax | ( | ) | const |
Get a HyperPoint that gives the maximum in each dim.
Definition at line 434 of file HyperPointSet.cpp.
HyperPoint HyperPointSet::getMin | ( | ) | const |
Get a HyperCuboid that surrounds the points.
Get a HyperPoint that gives the minimum in each dim
Definition at line 420 of file HyperPointSet.cpp.
double HyperPointSet::getSumW | ( | ) | const |
Find the Sum of weights for all HyperPoints in the HyperPointSet
Definition at line 99 of file HyperPointSet.cpp.
double HyperPointSet::getSumW2 | ( | ) | const |
Find the Sum of weights squared for all HyperPoints in the HyperPointSet
Definition at line 111 of file HyperPointSet.cpp.
HyperPointSet HyperPointSet::gramSchmidtProcess | ( | ) | const |
Use the gramSchmidtProcess to make the first nDim vectors orthanormal. If there are less than this, only size() vectors will be returned (but they will be orthanormal).
Definition at line 127 of file HyperPointSet.cpp.
HyperPoint HyperPointSet::harmonicMean | ( | ) | const |
Get the harmonic mean of all points. This does not use the weights.
Definition at line 483 of file HyperPointSet.cpp.
bool HyperPointSet::linearlyIndependant | ( | ) | const |
Check if the HyperPoints in the HyperPointSet are linearly independent (this isn't implemented yet!!!)
Definition at line 153 of file HyperPointSet.cpp.
void HyperPointSet::load | ( | TString | path | ) |
Load exisiting HyperPointSet from a file. HyperPoints are appended to the end of the
Definition at line 362 of file HyperPointSet.cpp.
HyperPoint HyperPointSet::mean | ( | ) | const |
Get the arithmatic mean of all points. This does not use the weights.
Definition at line 450 of file HyperPointSet.cpp.
void HyperPointSet::normalise | ( | ) |
void HyperPointSet::print | ( | std::ostream & | os = std::cout | ) | const |
Print out the entire HyperPointSet.
Definition at line 301 of file HyperPointSet.cpp.
void HyperPointSet::push_back | ( | const HyperPoint & | point | ) |
Add new HyperPoint to the HyperPointSet
Definition at line 91 of file HyperPointSet.cpp.
void HyperPointSet::removeDuplicates | ( | ) |
Remove all duplicate events. IMPORTANT must be followed by sort()
Definition at line 179 of file HyperPointSet.cpp.
|
private |
Save the HyperPointSet to the file currently in scope.
Definition at line 332 of file HyperPointSet.cpp.
void HyperPointSet::save | ( | TString | path | ) |
Save the HyperPointSet to a file (opened with RECREATE).
Definition at line 314 of file HyperPointSet.cpp.
unsigned int HyperPointSet::size | ( | ) | const |
void HyperPointSet::sort | ( | ) |
Sort all the HyperPoints in acending order. First sort by the first element, then by the second etc.
Definition at line 171 of file HyperPointSet.cpp.
|
private |
The dimensionality of HyperPoints in the HyperPointSet
Definition at line 33 of file HyperPointSet.h.
|
private |
std::vector containing the HyperPoints
Definition at line 34 of file HyperPointSet.h.