GWOLLUM  4.2.0
Tools for gravitational-wave analyses
Segments Class Reference

Manage time segment lists. More...

#include <Segments.h>

Inheritance diagram for Segments:

Public Member Functions

bool AddSegment (const double aGpsStart, const double aGpsEnd, const int aTag=1)
 Adds a segment to the list. More...
 
bool AddSegments (Segments *aSeg)
 Adds a segment list to the list. More...
 
bool Append (const double aGpsStart, const double aGpsEnd, const int aTag=1)
 Appends a segment to the list. More...
 
bool Append (Segments *aSeg)
 Appends segments to the list. More...
 
bool Append (TTree *aSegmentTree)
 Appends segments to the list. More...
 
bool Append (vector< double > aGpsStarts, vector< double > aGpsEnds, vector< double > aTags)
 Appends segments to the list. More...
 
bool ApplyPadding (const double aPadStart, const double aPadEnd)
 Apply a padding to all segments. More...
 
void Dump (unsigned int aNcols=2, const string aTxtFileName="")
 Dumps the current segment list in the standard output or in a file. More...
 
double GetDuration (const unsigned int aSegmentIndex=0)
 Returns the duration of a given segment. More...
 
double GetEnd (const unsigned int aSegmentIndex=0)
 Returns the ending time of a given segment. More...
 
vector< double > GetEnds (void)
 Returns the list of segment ends. More...
 
double GetFirst ()
 Returns the starting time of the first segment. More...
 
TH1I * GetHisto (const double aOffset=0.0)
 Returns an histogram representing the segments. More...
 
TH1I * GetHistogram (const double aGpsStart, const double aGpsEnd, const double aOffset=0.0)
 Returns an histogram representing the segments. More...
 
double GetLast ()
 Returns the ending time of the last segment. More...
 
double GetLiveTime (const double aGpsStart, const double aGpsEnd)
 Returns the livetime between two GPS times. More...
 
double GetLiveTime (const unsigned int aSegmentIndex)
 Returns the livetime of one segment. More...
 
double GetLiveTime (void)
 Returns the total livetime. More...
 
double GetLiveTimeWithTag (const int aTag)
 Returns the livetime for segments with a given tag. More...
 
unsigned int GetN (void)
 Returns the number of time segments in the list. More...
 
unsigned int GetNWithTag (const int aTag)
 Returns the number of time segments in the list with a given tag. More...
 
int GetSegmentStartBefore (const double aGps)
 Returns the segment starting just before a given GPS time. More...
 
double GetStart (const unsigned int aSegmentIndex=0)
 Returns the starting time of a given segment. More...
 
vector< double > GetStarts (void)
 Returns the list of segment starts. More...
 
bool GetStatus (void) const
 Returns the status of the Segments object. More...
 
int GetTag (const unsigned int aSegmentIndex=0)
 Returns the tag of a given segment. More...
 
vector< int > GetTags (void)
 Returns the list of segment tags. More...
 
TTree * GetTree (const double aGpsStart, const double aGpsEnd)
 Returns a formatted TTree representing the segments. More...
 
TTree * GetTree (void)
 Returns a formatted TTree representing the segments. More...
 
bool Intersect (const double aGpsStart, const double aGpsEnd, const int aTag=1)
 Intersects with a segment. More...
 
bool Intersect (Segments *aSeg)
 Intersects with a segment list. More...
 
bool IsInsideSegment (const double aGps)
 Checks whether a given GPS time is inside a segment of the list. More...
 
bool IsInsideSegment (unsigned int &aSegmentIndex, const double aGps)
 Checks whether a given GPS time is inside a segment of the list. More...
 
bool MakeSecond (void)
 Rounds segment boundaries to integer seconds. More...
 
void Reset (void)
 Resets the segment list. More...
 
bool Reverse (void)
 Reverses the segment list. More...
 
void SetTag (const unsigned int aSegmentIndex, const int aTag)
 Tags a given segment. More...
 
void SetTags (const int aTag)
 Tags all segments. More...
 
bool TruncateAfter (const double aGps)
 Removes segments after a given time. More...
 
bool TruncateBefore (const double aGps)
 Removes segments before a given time. More...
 
Constructors and destructors
 Segments (void)
 Default constructor of the Segments class. More...
 
 Segments (const string aSegmentFile, const int aDefaultTag=1)
 Constructor of the Segments class. More...
 
 Segments (TTree *aSegmentTree)
 Constructor of the Segments class. More...
 
 Segments (const string aFilePattern, const string aDirectory, const string aTreeName)
 Constructor of the Segments class. More...
 
 Segments (const vector< double > aGpsStarts, const vector< double > aGpsEnds)
 Constructor of the Segments class. More...
 
 Segments (const double aGpsStart, const double aGpsEnd, const int aTag=1)
 Constructor of the Segments class. More...
 
virtual ~Segments (void)
 Destructor of the Segments class. More...
 

Private Member Functions

bool CheckConsistency (const unsigned int aStartIndex)
 Checks the sanity of the segment list. More...
 
void ComputeLiveTime (void)
 Computes the sgments livetime. More...
 
int GetMark (const unsigned int aGps)
 Returns the time mark for a given GPS time. More...
 
int GetTagProduct (const int aTag1, const int aTag2)
 Computes the product of 2 tags. More...
 
void MergeOverlaps (const unsigned int aStartIndex)
 Merges segment overlaps. More...
 
void SetMarks (void)
 Sets time marks for internal lookup. More...
 

Private Attributes

vector< unsigned int > n_mark
 Segment index marks. More...
 
vector< double > seg_end
 List of segment ends. More...
 
double seg_livetime
 Integrated livetime. More...
 
vector< double > seg_start
 List of segment starts. More...
 
vector< int > seg_tag
 List of segment tags. More...
 
bool status
 Status of the Segments object. More...
 
vector< unsigned int > time_mark
 Time marks. More...
 

Detailed Description

Manage time segment lists.

Every gravitational-wave analysis has to deal with segment lists at some point. This class is designed to manage segment lists in an easy way.

A segment is a time interval limited by 2 GPS times: [GPS_START] (included) and [GPS_END] (excluded). A segment list must respect the rules documented in CheckConsistency()

Segments can be individually tagged by an signed integer [TAG].

This class offers methods to transform segment lists: new segments, intersection, padding, truncation...

Constructor & Destructor Documentation

◆ Segments() [1/6]

Segments::Segments ( void  )

Default constructor of the Segments class.

An empty segment list is initialized.

◆ Segments() [2/6]

Segments::Segments ( const string  aSegmentFile,
const int  aDefaultTag = 1 
)

Constructor of the Segments class.

The original segment list is downloaded from a text file. Several formats are supported to extract the segments and the tags:

  • 2 columns: [GPS_START] [GPS_END]. In this case the segment tags are initialized to the default value.
  • 3 or 4 columns: (second column)=[GPS_START], (third column)=[GPS_END]. In this case the segment tags are initialized to the default value.
  • 5 columns: (second column)=[GPS_START], (third column)=[GPS_END], (fifth column)=[TAG].
  • any other format: the Segments object is corrupted.
Note
Input segments can overlap in time (overlaps will be internally merged).
Parameters
[in]aSegmentFilePath to the segment text file.
[in]aDefaultTagDefault tag value aaplied to the segments when the number of columns is less than 5.

◆ Segments() [3/6]

Segments::Segments ( TTree *  aSegmentTree)

Constructor of the Segments class.

The original segment list is downloaded from a TTree. 3 branches must be found: start/D, end/D, and tag/I.

Note
Input segments can overlap in time (overlaps will be internally merged).
Parameters
[in]aSegmentTreePointer to the segment TTree structure.

◆ Segments() [4/6]

Segments::Segments ( const string  aFilePattern,
const string  aDirectory,
const string  aTreeName 
)

Constructor of the Segments class.

The original segment list is downloaded from a TTree in a list of ROOT files. 3 branches must be found: start/D, end/D, and tag/I.

Note
Input segments can overlap in time (overlaps will be internally merged).
Parameters
[in]aFilePatternInput file pattern.
[in]aDirectoryROOT sub-directory. Use "" if none.
[in]aTreeNameSegment tree name.

◆ Segments() [5/6]

Segments::Segments ( const vector< double >  aGpsStarts,
const vector< double >  aGpsEnds 
)

Constructor of the Segments class.

The original segment list is downloaded from vectors.

Note
Input segments can overlap in time (overlaps will be internally merged).
All segments are given a tag 1.
Parameters
[in]aGpsStartsVector of GPS starts.
[in]aGpsEndsVector of GPS ends.

◆ Segments() [6/6]

Segments::Segments ( const double  aGpsStart,
const double  aGpsEnd,
const int  aTag = 1 
)

Constructor of the Segments class.

The original segment list contains one single segment.

Parameters
[in]aGpsStartSegment GPS start.
[in]aGpsEndSegment GPS end.
[in]aTagSegment tag.

◆ ~Segments()

Segments::~Segments ( void  )
virtual

Destructor of the Segments class.

Member Function Documentation

◆ AddSegment()

bool Segments::AddSegment ( const double  aGpsStart,
const double  aGpsEnd,
const int  aTag = 1 
)

Adds a segment to the list.

The Segments object is automatically re-organized to remove overlaps.

Returns
true if the segment was correctly added, false otherwise.
Parameters
[in]aGpsStartSegment GPS start.
[in]aGpsEndSegment GPS end.
[in]aTagSegment tag.

◆ AddSegments()

bool Segments::AddSegments ( Segments aSeg)

Adds a segment list to the list.

The Segments object is automatically re-organized to remove overlaps.

Returns
true if the segments were correctly added, false otherwise.
Parameters
[in]aSegSegment list to add.

◆ Append() [1/4]

bool Segments::Append ( const double  aGpsStart,
const double  aGpsEnd,
const int  aTag = 1 
)

Appends a segment to the list.

This function appends a time segment at the end of the current Segments object. It only works if the new segment starts after the last segment start (an overlap is possible).

Returns
true if the segment was correctly added, false otherwise.
Parameters
[in]aGpsStartSegment GPS start.
[in]aGpsEndSegment GPS end.
[in]aTagSegment tag.

◆ Append() [2/4]

bool Segments::Append ( Segments aSeg)

Appends segments to the list.

This function appends a list of segments at the end of the current Segments object. It only works if the first segment starts after the last segment start (an overlap is possible).

Returns
true if the segments were correctly added, false otherwise.
Parameters
[in]aSegSegment list to append.

◆ Append() [3/4]

bool Segments::Append ( TTree *  aSegmentTree)

Appends segments to the list.

The Segments object is automatically re-organized to remove overlaps.

Note
Three branches must be found: start/D, end/D, and tag/I.
Returns
true if the segments were correctly added, false otherwise.
Parameters
[in]aSegmentTreeSegment tree with 3 branches: "start", "end", and "tag".

◆ Append() [4/4]

bool Segments::Append ( vector< double >  aGpsStarts,
vector< double >  aGpsEnds,
vector< double >  aTags 
)

Appends segments to the list.

This function appends a list of segments at the end of the current Segments object. It only works if the first segment starts after the last segment start (an overlap is possible).

Returns
true if the segments were correctly added, false otherwise.
Parameters
[in]aGpsStartsSegment GPS starts.
[in]aGpsEndsSegment GPS ends.
[in]aTagsSegment tags.

◆ ApplyPadding()

bool Segments::ApplyPadding ( const double  aPadStart,
const double  aPadEnd 
)

Apply a padding to all segments.

This function shifts the segment start and end by a given padding value. The Segments object is automatically re-organized to remove overlaps. Some segments are removed if the duration becomes null or negative.

Returns
true if the padding was correctly added, false otherwise.
Parameters
[in]aPadStartPadding value added to the segment starts.
[in]aPadEndPadding value added to the segment ends.

◆ CheckConsistency()

bool Segments::CheckConsistency ( const unsigned int  aStartIndex)
private

Checks the sanity of the segment list.

List of checks:

  • The list of starts, ends and tags must be the same size.
  • All segments must start after SEGMENTS_START
  • All segments must end after SEGMENTS_END
  • All segments must have a strictly positive duration.
  • Segments must be sorted in increasing [GPS_START].
Parameters
[in]aStartIndexSegment index where to start the checks.
Returns
false if the segment list is not sane, true otherwise.

◆ ComputeLiveTime()

void Segments::ComputeLiveTime ( void  )
private

Computes the sgments livetime.

Note
This function assumes the segments consistency.

◆ Dump()

void Segments::Dump ( unsigned int  aNcols = 2,
const string  aTxtFileName = "" 
)

Dumps the current segment list in the standard output or in a file.

Several formats are possible:

  • if 'aNcols=0' The list of tags is printed
  • if 'aNcols=1' The list of duration is printed
  • if 'aNcols=2' The list of segment starts and ends is printed
  • if 'aNcols=3' The list of segment starts and ends is printed, as well as the segment index
  • if 'aNcols=4' The list of segment starts and ends is printed, as well as the segment index and duration
  • if 'aNcols>=5' The list of segment starts and ends is printed, as well as the segment index, duration and tag
Parameters
[in]aNcolsnumber of columns in the output.
[in]aTxtFileNameText file name. Use "" for the standard output.

◆ GetDuration()

double Segments::GetDuration ( const unsigned int  aSegmentIndex = 0)
inline

Returns the duration of a given segment.

Parameters
[in]aSegmentIndexSegment index.
Precondition
The segment index must be valid (less than the number of segments returned by GetN()).

◆ GetEnd()

double Segments::GetEnd ( const unsigned int  aSegmentIndex = 0)
inline

Returns the ending time of a given segment.

Parameters
[in]aSegmentIndexSegment index.
Precondition
The segment index must be valid (less than the number of segments returned by GetN()).

◆ GetEnds()

vector<double> Segments::GetEnds ( void  )
inline

Returns the list of segment ends.

◆ GetFirst()

double Segments::GetFirst ( )
inline

Returns the starting time of the first segment.

Precondition
There should be at least one segment in the list or else this function crashes.

◆ GetHisto()

TH1I* Segments::GetHisto ( const double  aOffset = 0.0)
inline

Returns an histogram representing the segments.

The GPS time is represented on the X axis:

  • 1 is used when inside a segment.
  • 0 is used when outside a segment.
    Parameters
    [in]aOffsetTime offset applied to GPS times.
    Note
    The user is in charge of deleting this histogram.
    Returns
    A NULL pointer if the Segments object is corrupted.

◆ GetHistogram()

TH1I * Segments::GetHistogram ( const double  aGpsStart,
const double  aGpsEnd,
const double  aOffset = 0.0 
)

Returns an histogram representing the segments.

Only segments between 'aGpsStart' and 'aGpsEnd' are represented. The GPS time is represented on the X axis:

  • 1 is used when inside a segment.
  • 0 is used when outside a segment.
    Parameters
    [in]aGpsStartStarting GPS time of the histogram.
    [in]aGpsEndEnding GPS time of the histogram.
    [in]aOffsetTime offset applied to GPS times.
    Note
    The user is in charge of deleting this histogram.
    Returns
    A NULL pointer if the Segments object is corrupted.

◆ GetLast()

double Segments::GetLast ( )
inline

Returns the ending time of the last segment.

Precondition
There should be at least one segment in the list or else this function crashes.

◆ GetLiveTime() [1/3]

double Segments::GetLiveTime ( const double  aGpsStart,
const double  aGpsEnd 
)

Returns the livetime between two GPS times.

Parameters
[in]aGpsStartStarting GPS time.
[in]aGpsEndEnding GPS time.

◆ GetLiveTime() [2/3]

double Segments::GetLiveTime ( const unsigned int  aSegmentIndex)
inline

Returns the livetime of one segment.

Parameters
[in]aSegmentIndexSegment index.
Precondition
The segment index must be valid (less than the number of segments returned by GetN()).
See also
GetDuration().

◆ GetLiveTime() [3/3]

double Segments::GetLiveTime ( void  )
inline

Returns the total livetime.

Every segments in the list are considered.

◆ GetLiveTimeWithTag()

double Segments::GetLiveTimeWithTag ( const int  aTag)

Returns the livetime for segments with a given tag.

Parameters
[in]aTagTag value of segments to consider.

◆ GetMark()

int Segments::GetMark ( const unsigned int  aGps)
private

Returns the time mark for a given GPS time.

Returns
-1 is returned if the GPS time is before the first segment.
Parameters
[in]aGpsGPS time of interest.

◆ GetN()

unsigned int Segments::GetN ( void  )
inline

Returns the number of time segments in the list.

◆ GetNWithTag()

unsigned int Segments::GetNWithTag ( const int  aTag)

Returns the number of time segments in the list with a given tag.

Parameters
[in]aTagTag value of segments to consider.

◆ GetSegmentStartBefore()

int Segments::GetSegmentStartBefore ( const double  aGps)

Returns the segment starting just before a given GPS time.

Returns
The segment index or -1 if no segment starts before the GPS time.
Parameters
[in]aGpsGPS time of interest

◆ GetStart()

double Segments::GetStart ( const unsigned int  aSegmentIndex = 0)
inline

Returns the starting time of a given segment.

Parameters
[in]aSegmentIndexSegment index.
Precondition
The segment index must be valid (less than the number of segments returned by GetN()).

◆ GetStarts()

vector<double> Segments::GetStarts ( void  )
inline

Returns the list of segment starts.

◆ GetStatus()

bool Segments::GetStatus ( void  ) const
inline

Returns the status of the Segments object.

◆ GetTag()

int Segments::GetTag ( const unsigned int  aSegmentIndex = 0)
inline

Returns the tag of a given segment.

Parameters
[in]aSegmentIndexSegment index.
Precondition
The segment index must be valid (less than the number of segments returned by GetN()).

◆ GetTagProduct()

int Segments::GetTagProduct ( const int  aTag1,
const int  aTag2 
)
inlineprivate

Computes the product of 2 tags.

  • If one of the 2 tags is strictly negative, the lowest value is returned.
  • If both tags are positive, the greatest value is returned
    Parameters
    [in]aTag1First tag.
    [in]aTag2Second tag.

◆ GetTags()

vector<int> Segments::GetTags ( void  )
inline

Returns the list of segment tags.

◆ GetTree() [1/2]

TTree * Segments::GetTree ( const double  aGpsStart,
const double  aGpsEnd 
)

Returns a formatted TTree representing the segments.

This tree has 3 branches: "start", "end", and "tag". Only segments between 'aGpsStart' and 'aGpsEnd' are copied in the tree.

Parameters
[in]aGpsStartStarting GPS time.
[in]aGpsEndEnding GPS time.
Note
The user is in charge of deleting this TTree.
Returns
A NULL pointer if the Segments object is corrupted.

◆ GetTree() [2/2]

TTree* Segments::GetTree ( void  )
inline

Returns a formatted TTree representing the segments.

This tree has 3 branches: "start", "end", and "tag".

Note
The user is in charge of deleting this TTree.
Returns
A NULL pointer if the Segments object is corrupted.

◆ Intersect() [1/2]

bool Segments::Intersect ( const double  aGpsStart,
const double  aGpsEnd,
const int  aTag = 1 
)

Intersects with a segment.

This function intersects the current segment list with a segment given in argument.

Warning
This function is destructive and there is no way to go back to the original segment list.
Returns
true if the segments were correctly intersected, false otherwise.
Note
The segment tags are obtained using GetTagProduct() for the overlapping segments.
Parameters
[in]aGpsStartSegment GPS start.
[in]aGpsEndSegment GPS end.
[in]aTagSegment tag.

◆ Intersect() [2/2]

bool Segments::Intersect ( Segments aSeg)

Intersects with a segment list.

This function intersects the current segment list with a second segment list given in argument.

Warning
This function is destructive and there is no way to go back to the original segment list.
Returns
true if the segments were correctly intersected, false otherwise.
Note
The segment tags are obtained using GetTagProduct() for the overlapping segments.
Parameters
[in]aSegExternal Segments objects (not modified)

◆ IsInsideSegment() [1/2]

bool Segments::IsInsideSegment ( const double  aGps)

Checks whether a given GPS time is inside a segment of the list.

Returns
This function returns true if the GPS value is found inside a time segment. It returns false otherwise
Parameters
[in]aGpsGPS time of interest.

◆ IsInsideSegment() [2/2]

bool Segments::IsInsideSegment ( unsigned int &  aSegmentIndex,
const double  aGps 
)

Checks whether a given GPS time is inside a segment of the list.

Returns
This function returns true if the GPS value is found inside a time segment. It returns false otherwise
Parameters
[out]aSegmentIndexSegment index of segment containing the GPS time. This is irrelevant if false is returned.
[in]aGpsGPS time of interest.

◆ MakeSecond()

bool Segments::MakeSecond ( void  )

Rounds segment boundaries to integer seconds.

Segments starts are rounded to the previous integer and segments ends are rounded to the next integer.

Note
The Segments object is automatically re-organized to remove overlaps.
Returns
true if the segment was correctly modified, false otherwise.

◆ MergeOverlaps()

void Segments::MergeOverlaps ( const unsigned int  aStartIndex)
private

Merges segment overlaps.

When merging 2 segments, the tag with the highest value is kept.

Note
This function assumes the segments consistency.
Parameters
[in]aStartIndexSegment index where to start the merge.

◆ Reset()

void Segments::Reset ( void  )

Resets the segment list.

As a result, the segment list is empty and the status is back to OK.

◆ Reverse()

bool Segments::Reverse ( void  )

Reverses the segment list.

It means that the gaps between segments are now considered as segments. The first segment starts at :: SEGMENTS_START. The last segment ends at :: SEGMENTS_END.

Returns
true if the segments were correctly reversed, false otherwise.
Note
All tags are lost in the reversing process. The new segments are tagged with 1.

◆ SetMarks()

void Segments::SetMarks ( void  )
private

Sets time marks for internal lookup.

Note
This function assumes the segments consistency.

◆ SetTag()

void Segments::SetTag ( const unsigned int  aSegmentIndex,
const int  aTag 
)
inline

Tags a given segment.

Parameters
[in]aSegmentIndexSegment index.
[in]aTagNew tag value.
Precondition
The segment index must be valid (less than the number of segments returned by GetN()).

◆ SetTags()

void Segments::SetTags ( const int  aTag)

Tags all segments.

Parameters
[in]aTagNew tag value.

◆ TruncateAfter()

bool Segments::TruncateAfter ( const double  aGps)

Removes segments after a given time.

Note
if the GPS time is inside a segment, the segment is truncated.
Parameters
[in]aGpsGPS time after which to truncate.

◆ TruncateBefore()

bool Segments::TruncateBefore ( const double  aGps)

Removes segments before a given time.

Note
if the GPS time is inside a segment, the segment is truncated.
Parameters
[in]aGpsGPS time before which to truncate.

Member Data Documentation

◆ n_mark

vector<unsigned int> Segments::n_mark
private

Segment index marks.

◆ seg_end

vector<double> Segments::seg_end
private

List of segment ends.

◆ seg_livetime

double Segments::seg_livetime
private

Integrated livetime.

◆ seg_start

vector<double> Segments::seg_start
private

List of segment starts.

◆ seg_tag

vector<int> Segments::seg_tag
private

List of segment tags.

◆ status

bool Segments::status
private

Status of the Segments object.

◆ time_mark

vector<unsigned int> Segments::time_mark
private

Time marks.


The documentation for this class was generated from the following files: