libatomprobe
Library for Atom Probe Tomography (APT) computation
|
Data storage and retrieval class for "ranging" a spectra, where overlapping ranges are permitted. More...
#include <multiRange.h>
Public Member Functions | |
MultiRange () | |
MultiRange (const AtomProbe::RangeFile &rng, const AtomProbe::AbundanceData &natData) | |
Do our best to initialise from a rangefile. More... | |
bool | operator== (const MultiRange &oth) const |
Check for equality with other multirange. More... | |
unsigned int | addIon (const std::set< SIMPLE_SPECIES > &molecule, const std::string &name, const RGBf &ionCol) |
Add the ion to the database returns ion ID if successful, -1 otherwise. More... | |
unsigned int | addIon (const SIMPLE_SPECIES &molecule, const std::string &name, const RGBf &ionCol) |
Add a simple ion to the database - returns ion ID if successful, -1 otherwise. More... | |
unsigned int | addRange (float start, float end, unsigned int ionID) |
Add a range to the rangefile. Returns ID number of added range if adding successful, (unsigned int)-1 otherwise. More... | |
unsigned int | addRange (const std::pair< float, float > &rng, unsigned int ionID) |
Convenience wrapper for addRange(float,float, unsigned int) More... | |
void | setRangeGroups (const std::vector< unsigned int > &groups) |
Set the groupings for the ranges : this is used when splitting overlapping ranges. More... | |
bool | write (const char *fileName, unsigned int format=MULTIRANGE_FORMAT_XML) const |
Save the structure to a file. format specifies output file format type. More... | |
bool | open (const char *fileName, unsigned int format=-1) |
read the contents of the file into class More... | |
unsigned int | getIonID (unsigned int rangeId) const |
Get the ion's ID from a specified mass. More... | |
void | setIonID (unsigned int range, unsigned int newIonId) |
Set the ion ID for a given range. More... | |
RGBf | getColour (unsigned int ionID) const |
Obtain the colour for a given ion. More... | |
void | setColour (unsigned int ionID, const RGBf &r) |
Set the colour using the ion ID. More... | |
std::string | getErrString () const |
Retrieve the human-readable error associated with the current range file state. More... | |
std::string | getIonName (unsigned int ionID) const |
Get the name of a specified ionID. More... | |
unsigned int | getNumRanges () const |
Get the number of unique ranges. More... | |
unsigned int | getNumRanges (unsigned int ionID) const |
Get the number of ranges for a given ion ID. More... | |
unsigned int | getNumIons () const |
Get the number of unique ions. More... | |
std::set< SIMPLE_SPECIES > | getMolecule (unsigned int ionID) const |
Return the molecule that is associated with this ion. More... | |
std::pair< float, float > | getRange (unsigned int rangeID) const |
Retrieve the start and end of a given range as a pair(start,end) More... | |
std::vector< unsigned int > | getRanges (float mass) const |
Retrieve all of the ranges that specify the given mass. More... | |
bool | isRanged (float mass) const |
Returns true if a specified mass is ranged. More... | |
bool | isRanged (const IonHit &) const |
Returns true if an ion is ranged. More... | |
void | range (std::vector< IonHit > &ionHits) const |
Clips out ions that are not inside the rangefile's ranges. More... | |
bool | isSelfConsistent () const |
Check to see if data structure is internally consistent. More... | |
void | clear () |
Erase the contents of the rangefile. More... | |
void | flattenToMassAxis (std::vector< FLATTENED_RANGE > &ionMapping, float tolerance=0) const |
Obtain a projection onto the mass axis of ranges that do not touch one another. More... | |
void | splitOverlapping (std::vector< MultiRange > &decomposedRanges, float massTolerance=0) const |
void | copyDataFromRange (const MultiRange &src, unsigned int srcRngId) |
Copy range from a different multirange into this one, including all dependant data. More... | |
Data storage and retrieval class for "ranging" a spectra, where overlapping ranges are permitted.
Definition at line 63 of file multiRange.h.
AtomProbe::MultiRange::MultiRange | ( | ) |
Definition at line 136 of file multiRange.cpp.
AtomProbe::MultiRange::MultiRange | ( | const AtomProbe::RangeFile & | rng, |
const AtomProbe::AbundanceData & | natData | ||
) |
Do our best to initialise from a rangefile.
ions and their child ranges may be dropped if they cannot be identified
Definition at line 141 of file multiRange.cpp.
References ASSERT, AtomProbe::SIMPLE_SPECIES::atomicNumber, AtomProbe::SIMPLE_SPECIES::count, AtomProbe::RangeFile::decomposeIonNames(), AtomProbe::AbundanceData::getAtomicNumber(), AtomProbe::RangeFile::getColour(), AtomProbe::RangeFile::getIonID(), AtomProbe::RangeFile::getName(), AtomProbe::RangeFile::getNumIons(), AtomProbe::RangeFile::getNumRanges(), AtomProbe::RangeFile::getRange(), AtomProbe::AbundanceData::getSymbolIndices(), and isSelfConsistent().
unsigned int AtomProbe::MultiRange::addIon | ( | const std::set< SIMPLE_SPECIES > & | molecule, |
const std::string & | name, | ||
const RGBf & | ionCol | ||
) |
Add the ion to the database returns ion ID if successful, -1 otherwise.
Definition at line 271 of file multiRange.cpp.
References ASSERT, and isSelfConsistent().
Referenced by addIon(), AtomProbe::computeRangeAdjacency(), copyDataFromRange(), and AtomProbe::leastSquaresOverlapSolve().
unsigned int AtomProbe::MultiRange::addIon | ( | const SIMPLE_SPECIES & | molecule, |
const std::string & | name, | ||
const RGBf & | ionCol | ||
) |
Add a simple ion to the database - returns ion ID if successful, -1 otherwise.
Definition at line 293 of file multiRange.cpp.
References addIon(), ASSERT, and isSelfConsistent().
unsigned int AtomProbe::MultiRange::addRange | ( | float | start, |
float | end, | ||
unsigned int | ionID | ||
) |
Add a range to the rangefile. Returns ID number of added range if adding successful, (unsigned int)-1 otherwise.
Note that ranges are, unlike RangeFile, allowed to overlap, and still be self consistent.
Definition at line 303 of file multiRange.cpp.
References ASSERT, and isSelfConsistent().
Referenced by AtomProbe::computeRangeAdjacency(), copyDataFromRange(), and AtomProbe::leastSquaresOverlapSolve().
unsigned int AtomProbe::MultiRange::addRange | ( | const std::pair< float, float > & | rng, |
unsigned int | ionID | ||
) |
Convenience wrapper for addRange(float,float, unsigned int)
void AtomProbe::MultiRange::clear | ( | ) |
Erase the contents of the rangefile.
Definition at line 831 of file multiRange.cpp.
Referenced by open().
void AtomProbe::MultiRange::copyDataFromRange | ( | const MultiRange & | src, |
unsigned int | srcRngId | ||
) |
Copy range from a different multirange into this one, including all dependant data.
Definition at line 983 of file multiRange.cpp.
References addIon(), addRange(), ASSERT, AtomProbe::SIMPLE_SPECIES::atomicNumber, AtomProbe::ISOTOPE_ENTRY::atomicNumber, AtomProbe::RGBf::blue, AtomProbe::SIMPLE_SPECIES::count, AtomProbe::RGBf::fromHex(), AtomProbe::AbundanceData::generateSingleAtomDist(), getIonName(), getNumIons(), getNumRanges(), AtomProbe::AbundanceData::getSymbolIndices(), AtomProbe::RGBf::green, AtomProbe::AbundanceData::isotope(), isSelfConsistent(), MASS_TOL, open(), AtomProbe::RGBf::red, setRangeGroups(), splitOverlapping(), AtomProbe::AbundanceData::symbolIndex(), TEST, and write().
Referenced by splitOverlapping().
void AtomProbe::MultiRange::flattenToMassAxis | ( | std::vector< FLATTENED_RANGE > & | ionMapping, |
float | tolerance = 0 |
||
) | const |
Obtain a projection onto the mass axis of ranges that do not touch one another.
Definition at line 915 of file multiRange.cpp.
References ASSERT, and AtomProbe::FLATTENED_RANGE::startMass.
Referenced by splitOverlapping().
RGBf AtomProbe::MultiRange::getColour | ( | unsigned int | ionID | ) | const |
Obtain the colour for a given ion.
Definition at line 785 of file multiRange.cpp.
References ASSERT.
std::string AtomProbe::MultiRange::getErrString | ( | ) | const |
Retrieve the human-readable error associated with the current range file state.
Definition at line 693 of file multiRange.cpp.
References ASSERT, and AtomProbe::MULTIRANGE_ERR_ENUM_END.
unsigned int AtomProbe::MultiRange::getIonID | ( | unsigned int | rangeId | ) | const |
Get the ion's ID from a specified mass.
Returns the ions ID if there exists a range that contains this mass. Otherwise (unsigned int)-1 is returned
Definition at line 752 of file multiRange.cpp.
References ASSERT, and range().
std::string AtomProbe::MultiRange::getIonName | ( | unsigned int | ionID | ) | const |
Get the name of a specified ionID.
Definition at line 733 of file multiRange.cpp.
References ASSERT.
Referenced by copyDataFromRange(), and AtomProbe::leastSquaresOverlapSolve().
set< SIMPLE_SPECIES > AtomProbe::MultiRange::getMolecule | ( | unsigned int | ionID | ) | const |
Return the molecule that is associated with this ion.
Definition at line 727 of file multiRange.cpp.
References ASSERT.
Referenced by AtomProbe::computeIonDistAdjacency().
unsigned int AtomProbe::MultiRange::getNumIons | ( | ) | const |
Get the number of unique ions.
Definition at line 747 of file multiRange.cpp.
Referenced by AtomProbe::computeIonDistAdjacency(), copyDataFromRange(), and AtomProbe::leastSquaresOverlapSolve().
unsigned int AtomProbe::MultiRange::getNumRanges | ( | ) | const |
Get the number of unique ranges.
Definition at line 720 of file multiRange.cpp.
Referenced by AtomProbe::computeRangeAdjacency(), copyDataFromRange(), and AtomProbe::leastSquaresOverlapSolve().
unsigned int AtomProbe::MultiRange::getNumRanges | ( | unsigned int | ionID | ) | const |
Get the number of ranges for a given ion ID.
Definition at line 739 of file multiRange.cpp.
References ASSERT.
std::pair< float, float > AtomProbe::MultiRange::getRange | ( | unsigned int | rangeID | ) | const |
Retrieve the start and end of a given range as a pair(start,end)
Definition at line 759 of file multiRange.cpp.
References ASSERT.
Referenced by AtomProbe::computeRangeAdjacency().
std::vector< unsigned int > AtomProbe::MultiRange::getRanges | ( | float | mass | ) | const |
Retrieve all of the ranges that specify the given mass.
Definition at line 765 of file multiRange.cpp.
bool AtomProbe::MultiRange::isRanged | ( | float | mass | ) | const |
Returns true if a specified mass is ranged.
Definition at line 797 of file multiRange.cpp.
References range().
Referenced by isRanged(), and range().
bool AtomProbe::MultiRange::isRanged | ( | const IonHit & | ion | ) | const |
Returns true if an ion is ranged.
Definition at line 809 of file multiRange.cpp.
References AtomProbe::IonHit::getMassToCharge(), and isRanged().
bool AtomProbe::MultiRange::isSelfConsistent | ( | ) | const |
Check to see if data structure is internally consistent.
Definition at line 230 of file multiRange.cpp.
Referenced by addIon(), addRange(), copyDataFromRange(), MultiRange(), open(), and write().
bool AtomProbe::MultiRange::open | ( | const char * | fileName, |
unsigned int | format = -1 |
||
) |
read the contents of the file into class
returns true if successful, otherwise false
Definition at line 398 of file multiRange.cpp.
References AtomProbe::SIMPLE_SPECIES::atomicNumber, AtomProbe::RGBf::blue, clear(), AtomProbe::SIMPLE_SPECIES::count, AtomProbe::RGBf::green, isSelfConsistent(), AtomProbe::MULTIRANGE_ERR_BAD_GROUP, AtomProbe::MULTIRANGE_ERR_BAD_RANGE, AtomProbe::MULTIRANGE_ERR_BAD_RANGES, AtomProbe::MULTIRANGE_ERR_MOL_BAD_COLOUR, AtomProbe::MULTIRANGE_ERR_MOL_MISSING_COLOUR, AtomProbe::MULTIRANGE_ERR_MOL_MISSING_COMPONENTS, AtomProbe::MULTIRANGE_ERR_MOL_MISSING_ENTRY, AtomProbe::MULTIRANGE_ERR_MOL_MISSING_NAME, AtomProbe::MULTIRANGE_ERR_NO_CONTENT, AtomProbe::MULTIRANGE_ERR_NO_MOLECULE_ENTRIES, AtomProbe::MULTIRANGE_ERR_NO_MOLECULES, AtomProbe::MULTIRANGE_ERR_NO_NODEPTR, AtomProbe::MULTIRANGE_ERR_NO_PARSER, AtomProbe::MULTIRANGE_ERR_NO_RANGES, AtomProbe::MULTIRANGE_ERR_NO_TIMESTAMP, AtomProbe::MULTIRANGE_ERR_NO_XML_CONTEXT, AtomProbe::MULTIRANGE_ERR_NO_XML_DOC, AtomProbe::MULTIRANGE_ERR_RANGE_NOPARENT, AtomProbe::MULTIRANGE_ERR_SELF_INCONSISTENT, AtomProbe::MULTIRANGE_NOT_VALID_ROOTNODE, AtomProbe::parseColString(), AtomProbe::RGBf::red, AtomProbe::stream_cast(), AtomProbe::XMLFreeDoc(), AtomProbe::XMLHelpFwdToElem(), and AtomProbe::XMLHelpGetProp().
Referenced by copyDataFromRange().
bool AtomProbe::MultiRange::operator== | ( | const MultiRange & | oth | ) | const |
Check for equality with other multirange.
Definition at line 208 of file multiRange.cpp.
void AtomProbe::MultiRange::range | ( | std::vector< IonHit > & | ionHits | ) | const |
Clips out ions that are not inside the rangefile's ranges.
Definition at line 814 of file multiRange.cpp.
References isRanged().
Referenced by getIonID(), isRanged(), and setIonID().
void AtomProbe::MultiRange::setColour | ( | unsigned int | ionID, |
const RGBf & | r | ||
) |
void AtomProbe::MultiRange::setIonID | ( | unsigned int | range, |
unsigned int | newIonId | ||
) |
Set the ion ID for a given range.
Definition at line 778 of file multiRange.cpp.
References ASSERT, and range().
void AtomProbe::MultiRange::setRangeGroups | ( | const std::vector< unsigned int > & | groups | ) |
Set the groupings for the ranges : this is used when splitting overlapping ranges.
When splitting overlapping ranges, "grouped" ranges are considered to be part of the same split grouping. This is typically used for example, to group isotopes that belong to an ion, but of differing charge state. As the intensities for charge state in an experiment are not linked There is no need to link these in the range group. Each range with a shared integer belongs to the same group
Definition at line 321 of file multiRange.cpp.
References ASSERT.
Referenced by copyDataFromRange(), and AtomProbe::leastSquaresOverlapSolve().
void AtomProbe::MultiRange::splitOverlapping | ( | std::vector< MultiRange > & | decomposedRanges, |
float | massTolerance = 0 |
||
) | const |
Separate out the non-interacting parts of the multi-ranges into their own multiRange entries The vector of maps converts the ionids of the child MultiRange to the parent
Definition at line 844 of file multiRange.cpp.
References ASSERT, copyDataFromRange(), flattenToMassAxis(), AtomProbe::linkIdentifiers(), and AtomProbe::pairOverlaps().
Referenced by copyDataFromRange(), and AtomProbe::leastSquaresOverlapSolve().
bool AtomProbe::MultiRange::write | ( | const char * | fileName, |
unsigned int | format = MULTIRANGE_FORMAT_XML |
||
) | const |
Save the structure to a file. format specifies output file format type.
Returns true if successful, otherwise false
Definition at line 327 of file multiRange.cpp.
References ASSERT, AtomProbe::LibVersion::getVersionStr(), isSelfConsistent(), AtomProbe::MULTIRANGE_FORMAT_XML, and AtomProbe::tabs().
Referenced by copyDataFromRange().