 |
OpenMS
2.4.0
|
Go to the documentation of this file.
41 #include <OpenMS/OpenMSConfig.h>
89 QualityType getOverallQuality()
const;
93 void setOverallQuality(QualityType q);
96 QualityType getQuality(
Size index)
const;
98 void setQuality(
Size index, QualityType q);
107 const std::vector<ConvexHull2D>& getConvexHulls()
const;
110 std::vector<ConvexHull2D>& getConvexHulls();
112 void setConvexHulls(
const std::vector<ConvexHull2D>& hulls);
122 bool encloses(
double rt,
double mz)
const;
132 const std::vector<Feature>& getSubordinates()
const;
135 std::vector<Feature>& getSubordinates();
138 void setSubordinates(
const std::vector<Feature>& rhs);
152 template <
typename Type>
155 Size assignments = 0;
156 assignments += ((*this).*member_function)();
157 for (std::vector<Feature>::iterator iter = subordinates_.begin(); iter != subordinates_.end(); ++iter)
159 assignments += iter->applyMemberFunction(member_function);
165 template <
typename Type>
168 Size assignments = 0;
169 assignments += ((*this).*member_function)();
170 for (std::vector<Feature>::const_iterator iter = subordinates_.begin(); iter != subordinates_.end(); ++iter)
172 assignments += iter->applyMemberFunction(member_function);
180 QualityType qualities_[2];
QualityLess OverallQualityLess
Compare by quality.
Definition: Feature.h:101
bool convex_hulls_modified_
Flag that indicates if the overall convex hull needs to be recomputed (i.e. mass trace convex hulls w...
Definition: Feature.h:186
A basic LC-MS feature.
Definition: BaseFeature.h:55
std::vector< Feature > subordinates_
subordinate features (e.g. features that represent alternative explanations, usually with lower quali...
Definition: Feature.h:192
size_t Size
Size type e.g. used as variable which can hold result of size()
Definition: Types.h:127
Size applyMemberFunction(Size(Type::*member_function)() const) const
The "const" variant.
Definition: Feature.h:166
std::vector< ConvexHull2D > convex_hulls_
Array of convex hulls (one for each mass trace)
Definition: Feature.h:183
A 2-dimensional hull representation in [counter]clockwise direction - depending on axis labelling.
Definition: ConvexHull2D.h:72
Main OpenMS namespace.
Definition: FeatureDeconvolution.h:46
ConvexHull2D convex_hull_
Overall convex hull of the feature.
Definition: Feature.h:189
bool operator==(_Iterator< _Val, _Ref, _Ptr > const &, _Iterator< _Val, _Ref, _Ptr > const &)
Definition: KDTree.h:806
An LC-MS feature.
Definition: Feature.h:70
Compare by quality.
Definition: BaseFeature.h:108
Size applyMemberFunction(Size(Type::*member_function)())
Applies a member function of Type to the feature (including subordinates). The returned values are ac...
Definition: Feature.h:153