libatomprobe
Library for Atom Probe Tomography (APT) computation
multiRange.h
Go to the documentation of this file.
1 /*
2  * multiRange.h - Atom probe rangefile class
3  * Copyright (C) 2014 D Haley
4  *
5  * This program is free software: you can redistribute it and/or modify
6  * it under the terms of the GNU General Public License as published by
7  * the Free Software Foundation, either version 3 of the License, or
8  * (at your option) any later version.
9  *
10  * This program is distributed in the hope that it will be useful,
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13  * GNU General Public License for more details.
14  *
15  * You should have received a copy of the GNU General Public License
16  * along with this program. If not, see <http://www.gnu.org/licenses/>.
17  */
18 
19 #ifndef ATOMPROBE_MULTIRANGE_H
20 #define ATOMPROBE_MULTIRANGE_H
21 
22 #include <vector>
23 #include <set>
24 #include <string>
25 #include <map>
26 
27 #include "atomprobe/io/ranges.h"
28 
29 #include "atomprobe/helper/misc.h"
31 
32 namespace AtomProbe
33 {
35 {
37  unsigned int atomicNumber;
38 
40  unsigned int count;
41 
42  bool operator==(const SIMPLE_SPECIES &oth) const;
43 };
44 
46 // a 1D mass axis
48 {
49  //Start and end of the "flattened" range
50  float startMass,endMass;
51  //The ionIDs from parent multirange that are contained between start and end mass
52  std::vector<unsigned int> containedIonIDs, containedRangeIDs;
53 };
54 
55 bool operator<(const SIMPLE_SPECIES &a, const SIMPLE_SPECIES &b);
56 
57 enum
58 {
60 };
61 
64 {
65  private:
67  std::vector<std::string> ionNames;
68 
70 
73  std::vector<std::set<SIMPLE_SPECIES> > moleculeData;
74 
76  std::vector<RGBf> colours;
77 
78 
80  std::vector<std::pair<float,float> > ranges;
81 
83  // Typically ions that are linked by isotope fingerprint share a grouping identifier.
84  // This might be used e.g. for decomposition/deconvolution
85  std::vector<unsigned int> rangeGrouping;
86 
88  std::vector<unsigned int> ionIDs;
89 
90 
92  unsigned int errState;
93 
94 
95  public:
96  MultiRange();
97 
99 
101  MultiRange(const AtomProbe::RangeFile &rng, const AtomProbe::AbundanceData &natData);
102 
104  bool operator==(const MultiRange &oth) const;
105 
107  unsigned int addIon(const std::set<SIMPLE_SPECIES> &molecule,
108  const std::string &name, const RGBf &ionCol);
109 
110 
112  unsigned int addIon(const SIMPLE_SPECIES &molecule,
113  const std::string &name, const RGBf &ionCol);
114 
116 
119  unsigned int addRange(float start, float end, unsigned int ionID);
120 
122  unsigned int addRange(const std::pair<float, float> &rng, unsigned int ionID);
123 
125 
133  void setRangeGroups(const std::vector<unsigned int> &groups);
134 
136 
138  bool write(const char *fileName, unsigned int format=MULTIRANGE_FORMAT_XML) const;
139 
141 
145  bool open(const char *fileName, unsigned int format=-1);
146 
148 
151  unsigned int getIonID(unsigned int rangeId) const;
152 
154  void setIonID(unsigned int range, unsigned int newIonId);
155 
157  RGBf getColour(unsigned int ionID) const;
159  void setColour(unsigned int ionID, const RGBf &r);
160 
162  std::string getErrString() const;
163 
165  std::string getIonName(unsigned int ionID) const;
166 
168  unsigned int getNumRanges() const;
169 
171  unsigned int getNumRanges(unsigned int ionID) const;
172 
174  unsigned int getNumIons() const;
175 
177  std::set<SIMPLE_SPECIES> getMolecule(unsigned int ionID) const;
178 
179 
181  std::pair<float,float> getRange(unsigned int rangeID) const;
182 
184  std::vector<unsigned int> getRanges(float mass) const;
185 
186 
188  bool isRanged(float mass) const;
190  bool isRanged(const IonHit &) const;
191 
193  void range(std::vector<IonHit> &ionHits) const;
194 
196  bool isSelfConsistent() const;
197 
199  void clear();
200 
201 
203  // - each ion that is in the range is listed
204  // - Tolerance is the tolerance in both the + and - directions to extend ranges when considering
205  // the existence of overlaps
206  void flattenToMassAxis(std::vector<FLATTENED_RANGE> &ionMapping,float tolerance=0) const;
207 
210  // This will only consider the ranges that exist in the multirange file
211  // massTolerance : used to determine overlapping groups, and is +- the range bounds
212  void splitOverlapping(std::vector<MultiRange> &decomposedRanges, float massTolerance=0) const;
213 
214 
216  // such as parent ion, molecule formation data, representation etc.
217  // Does not copy group data, as this makes no sense, instead range group data is erased if present
218  void copyDataFromRange(const MultiRange &src,unsigned int srcRngId);
219 
220 #ifdef DEBUG
221  static bool testSplit(const AtomProbe::AbundanceData &abundance);
223 
224  static bool test();
225 #endif
226 };
227 
228 }
229 #endif
unsigned int atomicNumber
Number of protons (element number)
Definition: multiRange.h:37
bool operator<(const SIMPLE_SPECIES &a, const SIMPLE_SPECIES &b)
Definition: multiRange.cpp:117
Data holder for colour as float.
Definition: misc.h:26
Data storage and retrieval class for "ranging" a spectra, where overlapping ranges are permitted...
Definition: multiRange.h:63
bool operator==(const SIMPLE_SPECIES &oth) const
Definition: multiRange.cpp:112
Class to load abundance information for natural isotopes.
Definition: abundance.h:54
Structure that allows for the multirange data to be mapped into.
Definition: multiRange.h:47
This is a data holding class for POS file ions, from.
Definition: ionHit.h:36
Data storage and retrieval class for various range files.
Definition: ranges.h:95
unsigned int count
Number of copies of this isotope.
Definition: multiRange.h:40
std::vector< unsigned int > containedRangeIDs
Definition: multiRange.h:52