OpenMS
ConvexHull2D Class Reference

#include <OpenMS/DATASTRUCTURES/ConvexHull2D.h>

Collaboration diagram for ConvexHull2D:
[legend]

Public Types

typedef DPosition< 2 > PointType
 
typedef std::vector< PointTypePointArrayType
 
typedef PointArrayType::size_type SizeType
 
typedef PointArrayType::const_iterator PointArrayTypeConstIterator
 
typedef std::map< PointType::CoordinateType, DBoundingBox< 1 > > HullPointType
 

Public Member Functions

 ConvexHull2D ()
 default constructor More...
 
 ConvexHull2D (const ConvexHull2D &)=default
 Copy constructor. More...
 
 ConvexHull2D (ConvexHull2D &&)=default
 Move constructor. More...
 
ConvexHull2Doperator= (const ConvexHull2D &rhs)
 assignment operator More...
 
ConvexHull2Doperator= (ConvexHull2D &&) &=default
 move assignment operator More...
 
bool operator== (const ConvexHull2D &rhs) const
 equality operator More...
 
void clear ()
 removes all points More...
 
const PointArrayTypegetHullPoints () const
 accessor for the outer points More...
 
void setHullPoints (const PointArrayType &points)
 accessor for the outer(!) points (no checking is performed if this is actually a convex hull) More...
 
DBoundingBox< 2 > getBoundingBox () const
 returns the bounding box of the feature hull points More...
 
bool addPoint (const PointType &point)
 
void addPoints (const PointArrayType &points)
 
Size compress ()
 Allows to reduce the disk/memory footprint of a hull. More...
 
void expandToBoundingBox ()
 
bool encloses (const PointType &point) const
 returns if the point lies in the feature hull More...
 

Protected Attributes

HullPointType map_points_
 internal structure maintaining the hull and enabling queries to encloses() More...
 
PointArrayType outer_points_
 just the list of points of the outer hull (derived from map_points_ or given by user) More...
 

Detailed Description

@brief A 2-dimensional hull representation in [counter]clockwise direction - depending on axis labelling.

The current implementation does not guarantee to produce convex hulls.
It can still store 'old' convex hulls from featureXML without problems, but does not support the enclose() query in this case,
and you will get an exception. As an alternative, you can use my_hull.getBoundingBox().encloses(), which yields similar results,
and will always work.

If you are creating new hull from peaks (e.g. during FeatureFinding), the generated hulls of a feature are defined as
a range in m/z dimension for each RT scan (thus might be non-convex). This has the advantage that one can clearly see
where points range within each scan (although missing points within this range are still not shown).
When hulls are created like this, the encloses() function is supported, and works as expected, i.e.
for the shape defined by this hull (view it in TOPPView) it answers whether the point is inside the shape.

However, once you store the hull in featureXML and load it again, the encloses() function is not supported any longer, because the old convex hulls did not save min&max for each scan. (to support encloses() at least for the new hulls, one would need to check if there exists a min&max value for each scan --> then the query would be valid and the inner representation can be filled. Old featureXML's are not supported in any case.)

The outer hullpoints can be queried by getHullPoints().

Improvement:
For chromatograms we could postprocess the input and remove points in intermediate RT scans, which are currently reported but make the number of points rather large.

Member Typedef Documentation

◆ HullPointType

◆ PointArrayType

typedef std::vector<PointType> PointArrayType

◆ PointArrayTypeConstIterator

typedef PointArrayType::const_iterator PointArrayTypeConstIterator

◆ PointType

typedef DPosition<2> PointType

◆ SizeType

typedef PointArrayType::size_type SizeType

Constructor & Destructor Documentation

◆ ConvexHull2D() [1/3]

default constructor

◆ ConvexHull2D() [2/3]

ConvexHull2D ( const ConvexHull2D )
default

Copy constructor.

◆ ConvexHull2D() [3/3]

ConvexHull2D ( ConvexHull2D &&  )
default

Move constructor.

Member Function Documentation

◆ addPoint()

bool addPoint ( const PointType point)

adds a point to the hull if it is not already contained. Returns if the point was added. this will trigger recomputation of the outer hull points (thus points set with setHullPoints() will be lost)

◆ addPoints()

void addPoints ( const PointArrayType points)

adds points to the hull if it is not already contained. this will trigger recomputation of the outer hull points (thus points set with setHullPoints() will be lost)

◆ clear()

void clear ( )

removes all points

◆ compress()

Size compress ( )

Allows to reduce the disk/memory footprint of a hull.

Removes points from the hull which lie on a straight line and do not contribute to the hulls shape. Should be called before saving to disk.

Example: Consider a series of 3 scans with the same dimension in m/z. After calling compress, the points from the second scan will be removed, since they do not contribute to the convex hull.

Returns
Number of removed scans

◆ encloses()

bool encloses ( const PointType point) const

returns if the point lies in the feature hull

This function is only supported if the hull is created using addPoint() or addPoints(), but not when using setHullPoints(). If you require the latter functionality, then augment this function.

Exceptions
Exception::NotImplementedif only hull points (outer_points_), but no internal structure (map_points_) is given

◆ expandToBoundingBox()

void expandToBoundingBox ( )
    @brief Expand a convex hull to its bounding box.

    This reduces the size of a convex hull to four points, its

bounding box, thus reducing size when storing the information. Note that this leads to an enclosed area that can be significantly larger than the original convex hull.

◆ getBoundingBox()

DBoundingBox<2> getBoundingBox ( ) const

returns the bounding box of the feature hull points

◆ getHullPoints()

const PointArrayType& getHullPoints ( ) const

accessor for the outer points

◆ operator=() [1/2]

ConvexHull2D& operator= ( const ConvexHull2D rhs)

assignment operator

◆ operator=() [2/2]

ConvexHull2D& operator= ( ConvexHull2D &&  ) &
default

move assignment operator

◆ operator==()

bool operator== ( const ConvexHull2D rhs) const

equality operator

◆ setHullPoints()

void setHullPoints ( const PointArrayType points)

accessor for the outer(!) points (no checking is performed if this is actually a convex hull)

Referenced by MRMTransitionGroupPicker::pickFragmentChromatograms(), and MRMTransitionGroupPicker::pickPrecursorChromatograms().

Member Data Documentation

◆ map_points_

HullPointType map_points_
protected

internal structure maintaining the hull and enabling queries to encloses()

◆ outer_points_

PointArrayType outer_points_
mutableprotected

just the list of points of the outer hull (derived from map_points_ or given by user)