Handle time segment lists. More...
#include <Segments.h>
Public Member Functions | |
bool | AddSegment (const double gps_start, const double gps_end, const int tag=-1) |
Adds a segment to the list. | |
bool | Append (Segments *Seg) |
Appends segments to the list. | |
bool | ApplyPadding (const double pad_low, const double pad_up) |
Apply a padding to segments. | |
bool | Dump (int ncol=2) |
Dumps the current segment list in the standard output. | |
int | GetClosestSegment (const double gps, vector< double > &seg) |
Returns the closest segment of a given GPS time. | |
double | GetDuration (const int s=0) |
Returns the duration of segment number s. | |
double | GetEnd (const int s=0) |
Returns the ending time of segment number s. | |
vector< double > | GetEnds (void) |
Returns segment ends. | |
TH1I * | GetHisto (const double gps_start, const double gps_end, const double offset=0.0) |
Returns an histogram representing the segments. | |
TH1I * | GetHisto (const double offset=0.0) |
Returns an histogram representing the segments. | |
double | GetLiveTime (const int s) |
Returns the livetime value of segment number s. | |
double | GetLiveTime (const double gps_start, const double gps_end) |
Returns the livetime value. | |
double | GetLiveTime (void) |
Returns the total livetime value. | |
double | GetLiveTimeWithTag (const int tag) |
Returns the livetime value for segments with a given tag. | |
int | GetNsegments (void) |
Returns the number of disjoint time segments in the list. | |
int | GetNSegmentsWithTag (const int tag) |
Returns the number of tagged segments in the list. | |
double | GetStart (const int s=0) |
Returns the starting time of segment number s. | |
vector< double > | GetStarts (void) |
Returns segment starts. | |
bool | GetStatus (void) const |
Returns the status of the Segments object. | |
int | GetTag (const int s) |
Returns the tag value of segment number s. | |
vector< int > | GetTags (void) |
Returns segment tags. | |
bool | GetTimeSeries (bool *timeseries, const int sampling_rate, const int gps_start, const int gps_end) |
Returns a time-series representing the segments. | |
TTree * | GetTree (const double gps_start, const double gps_end) |
Returns a formatted TTree representing the segments. | |
TTree * | GetTree (void) |
Returns a formatted TTree representing the segments. | |
bool | Intersect (const double gps_start, const double gps_end, const int tag=-1) |
Intersects with a segment. | |
bool | Intersect (Segments *Seg) |
Intersects two segment lists. | |
bool | IsInsideSegment (const double gps, int &seg_n) |
Checks whether a given GPS time is inside a segment of the list. | |
bool | IsInsideSegment (const double gps) |
Checks whether a given GPS time is inside a segment of the list. | |
void | Reset (void) |
Resets the segment list. | |
bool | Reverse (void) |
Reverses the segment list. | |
bool | SetTag (const int s, const int tag) |
Tags segment number s. | |
bool | SetTags (const int tag) |
Tags all segments. | |
double | SmartShift (void) |
Apply a 'smart' time shift. | |
bool | TruncateAfter (const double gps) |
Remove segments after a given time. | |
bool | Write (const string txtfilename, const int ncol=2) |
Writes segments in a text file. | |
bool | WriteROOT (const string rootfilename) |
Writes segments in a ROOT file. | |
Constructors and destructors | |
Segments (const double segment_start, const double segment_end, const int segment_tag=-1) | |
Constructor of the Segments class. | |
Segments (const vector< double > segment_starts, const vector< double > segment_ends) | |
Constructor of the Segments class. | |
Segments (TTree *segment_tree) | |
Constructor of the Segments class. | |
Segments (const string segment_file, const int ncol=2, const int default_tag=-1) | |
Constructor of the Segments class. | |
Segments (void) | |
Constructor of the Segments class. | |
virtual | ~Segments (void) |
Destructor of the Segments class. | |
Protected Attributes | |
double | livetime_tot |
integrated livetime | |
vector< double > | seg_end |
list of segment ends | |
vector< double > | seg_start |
list of segment starts | |
vector< int > | seg_tag |
list of segment tags | |
bool | status |
status of the Segments object | |
Friends | |
class | ReadTriggerSegments |
Handle time segment lists.
Every GW analysis has to deal with segment lists at some point. The class Segments was built to handle segment lists in an easy way. A segment is a GPS interval and the Segments class was designed to handle lists of segments. Once the Segments class is constructed, many functions offer operations on the segments. Convention for segments is described in the GWOLLUM convention.
Segments::Segments | ( | void | ) |
Segments::Segments | ( | const string | segment_file, | |
const int | ncol = 2 , |
|||
const int | default_tag = -1 | |||
) |
Constructor of the Segments class.
A segment list contained in a text file is loaded Note that input segments can overlap in time (overlaps will be internally merged).
segment_file | path to the segment text file. Several formats are supported and are descibed in the GWOLLUM convention for segments. | |
ncol | number of columns in the segment file. It defines the format of the input list. See the GWOLLUM convention for supported formats. If ncol<0, the format is automatically detected. | |
default_tag | a default tag value can be specified to be attributed to every segments. When the segment file has a 5-column format, this default value is not used unless it is >=0. In that case, only segments with tag=default_tag are loaded. |
Segments::Segments | ( | TTree * | segment_tree | ) |
Constructor of the Segments class.
A segment list contained in a TTree structure is loaded. Note that input segments can overlap in time (overlaps will be internally merged).
segment_tree | pointer to the segment TTree structure. See the GWOLLUM convention for the TTree container. |
Segments::Segments | ( | const vector< double > | segment_starts, | |
const vector< double > | segment_ends | |||
) |
Constructor of the Segments class.
A segment list is given as lists of segment starts and ends. The input segments must follow the GWOLLUM convention for segments. Note that input segments can overlap in time (overlaps will be internally merged).
segment_starts | list of segment starts (sorted by increasing GPS values). | |
segment_ends | list of segment ends |
Segments::Segments | ( | const double | segment_start, | |
const double | segment_end, | |||
const int | segment_tag = -1 | |||
) |
Constructor of the Segments class.
Constructor for a single segment. The input segment must follow the GWOLLUM convention for segments.
segment_start | segment start | |
segment_end | segment end | |
segment_tag | segment tag |
Segments::~Segments | ( | void | ) | [virtual] |
Destructor of the Segments class.
bool Segments::AddSegment | ( | const double | gps_start, | |
const double | gps_end, | |||
const int | tag = -1 | |||
) |
Adds a segment to the list.
This function adds the segment [gps_start; gps_end[ to the Segments object. The Segments object is automatically re-organized to remove the overlaps. Optionally the new segment can be tagged.
true is returned in case of success.
gps_start | segment start | |
gps_end | segment end | |
tag | segment tag |
bool Segments::Append | ( | Segments * | Seg | ) |
Appends segments to the list.
This function appends the segments in 'Seg' to the current segment object. This only works if the first segment of 'Seg' starts after the last segment (start) of the current segment object.
true is returned in case of success.
Seg | segment object to append |
bool Segments::ApplyPadding | ( | const double | pad_low, | |
const double | pad_up | |||
) |
Apply a padding to segments.
This function shifts the segments start and end by a given value: [gps_start; gps_end[ [gps_start+pad_low; gps_end+pad_up[
The Segments object is automatically re-organized to remove overlaps. Some segments are removed if the duration becomes null or negative.
true is returned in case of success.
pad_low | padding added to the segment lower boundary | |
pad_up | padding added to the segment upper boundary |
bool Segments::Dump | ( | int | ncol = 2 |
) |
Dumps the current segment list in the standard output.
Several formats are possible and are given by the 'ncol' argument. See the GWOLLUM convention for segments. In addition:
true is returned in case of success.
ncol | number of columns in the segment file. It defines the format of the input list. See the GWOLLUM convention for supported formats |
int Segments::GetClosestSegment | ( | const double | gps, | |
vector< double > & | seg | |||
) |
Returns the closest segment of a given GPS time.
This function fills the vector 'seg' with the segment which is found to be the closest in time of the GPS value gps. As a first priority the segment which contains 'gps' is returned. If 'gps' is not inside a segment, 'seg' is filled with the segment whose lower/upper bound is the closest of 'gps'.
In case of success, the segment number is returned. In case of failure, -1 is returned and 'seg' is empty.
gps | GPS time of interest | |
seg | returned segment as a vector (GPS_START GPS_END) |
double Segments::GetDuration | ( | const int | s = 0 |
) |
Returns the duration of segment number s.
Returns -1 if not found.
s | segment number |
double Segments::GetEnd | ( | const int | s = 0 |
) |
Returns the ending time of segment number s.
Returns -1 if not found.
s | segment number |
vector<double> Segments::GetEnds | ( | void | ) | [inline] |
Returns segment ends.
This function returns a vector which contains all the segment ends.
TH1I * Segments::GetHisto | ( | const double | gps_start, | |
const double | gps_end, | |||
const double | offset = 0.0 | |||
) |
Returns an histogram representing the segments.
Only segments between 'gps_start' and 'gps_end' are represented. The GPS time is represented on the X axis:
gps_start | starting GPS time of the histogram. | |
gps_end | ending GPS time of the histogram. | |
offset | time offset applied to segment GPS times. |
TH1I * Segments::GetHisto | ( | const double | offset = 0.0 |
) |
Returns an histogram representing the segments.
The GPS time is represented on the X axis:
offset | time offset applied to segment GPS times. |
double Segments::GetLiveTime | ( | const int | s | ) |
Returns the livetime value of segment number s.
Returns -1 if not found.
s | segment number |
double Segments::GetLiveTime | ( | const double | gps_start, | |
const double | gps_end | |||
) |
Returns the livetime value.
Only segments between 'gps_start' and 'gps_end' are considered.
gps_start | starting GPS time | |
gps_end | ending GPS time |
double Segments::GetLiveTime | ( | void | ) | [inline] |
Returns the total livetime value.
Every segments in the list are considered.
double Segments::GetLiveTimeWithTag | ( | const int | tag | ) |
Returns the livetime value for segments with a given tag.
tag | tag value |
int Segments::GetNsegments | ( | void | ) | [inline] |
Returns the number of disjoint time segments in the list.
int Segments::GetNSegmentsWithTag | ( | const int | tag | ) |
Returns the number of tagged segments in the list.
tag | tag value of segments to consider. |
double Segments::GetStart | ( | const int | s = 0 |
) |
Returns the starting time of segment number s.
Returns -1 if not found.
s | segment number |
vector<double> Segments::GetStarts | ( | void | ) | [inline] |
Returns segment starts.
This function returns a vector which contains all the segment starts.
bool Segments::GetStatus | ( | void | ) | const [inline] |
Returns the status of the Segments object.
int Segments::GetTag | ( | const int | s | ) |
Returns the tag value of segment number s.
Returns -1 if not found.
s | segment number |
vector<int> Segments::GetTags | ( | void | ) | [inline] |
Returns segment tags.
This function returns a vector which contains all the segment tags.
bool Segments::GetTimeSeries | ( | bool * | timeseries, | |
const int | sampling_rate, | |||
const int | gps_start, | |||
const int | gps_end | |||
) |
Returns a time-series representing the segments.
The returned vector 'timeseries' is filled with true when inside a segment, false otherwise. The memory pointed by timeseries must be allocated before calling this function. false is returned if this function fails.
IMPORTANT: the array pointed by timeseries must be the right size: (gps_end-gps_start)*sampling_rate. No check is performed by this function!
timeseries | pointer to the time-series array | |
sampling_rate | sampling rate [Hz] | |
gps_start | starting GPS time | |
gps_end | ending GPS time |
TTree * Segments::GetTree | ( | const double | gps_start, | |
const double | gps_end | |||
) |
Returns a formatted TTree representing the segments.
This tree is formatted as described in the GWOLLUM convention for segments. Only segments between 'gps_start' and 'gps_end' are copied in the tree.
gps_start | starting GPS time | |
gps_end | ending GPS time |
TTree * Segments::GetTree | ( | void | ) |
Returns a formatted TTree representing the segments.
This tree is formatted as described in the GWOLLUM convention for segments.
bool Segments::Intersect | ( | const double | gps_start, | |
const double | gps_end, | |||
const int | tag = -1 | |||
) |
Intersects with a segment.
This function intersects the current segment list with a segment given in argument. This function is destructive and there is no way to go back to the original segment list. Optionally, the input segment can be tagged.
true is returned in case of success.
gps_start | segment start | |
gps_end | segment end | |
tag | segment tag |
bool Segments::Intersect | ( | Segments * | Seg | ) |
Intersects two segment lists.
This function intersects the current segment list with a second segment list given in argument. This function is destructive and there is no way to go back to the original segment list.
true is returned in case of success.
Seg | external Segments objects (not modified) |
bool Segments::IsInsideSegment | ( | const double | gps, | |
int & | seg_n | |||
) |
Checks whether a given GPS time is inside a segment of the list.
This function returns true if the GPS value gps is found inside a time segment.
gps | GPS time of interest | |
seg_n | if true is returned, the segment number of the segment containing gps is returned. |
bool Segments::IsInsideSegment | ( | const double | gps | ) |
Checks whether a given GPS time is inside a segment of the list.
This function returns true if the GPS value gps is found inside a time segment.
gps | GPS time of interest |
void Segments::Reset | ( | void | ) |
bool Segments::Reverse | ( | void | ) |
Reverses the segment list.
This means that the gaps between segments are now considered as segments.
BEFORE: |-------| |-------| |-| |----| AFTER: ----| |------| |--| |---| |----
Of course the first and last segment are finite. They start at 700,000,000 and stop at 2,000,000,000 respectively.
true is returned in case of success.
bool Segments::SetTag | ( | const int | s, | |
const int | tag | |||
) |
Tags segment number s.
true is returned in case of success.
s | segment number | |
tag | new tag value |
bool Segments::SetTags | ( | const int | tag | ) |
Tags all segments.
true is returned in case of success.
tag | new tag value for all segments |
double Segments::SmartShift | ( | void | ) |
Apply a 'smart' time shift.
Each segment is time-shited such as the new shifted start is positioned where the segment used to end:
BEFORE: |----| |-| |----| |----| AFTER: |----| |-| |-|
Notice that if the gap is not long enough for the shifted segment, the segment end is adjusted not to overlap the next segment. This way, the resulting segments never overlap the original segments. As a consequence, the resulting livetime is always smaller or equal to the original livetime. The final livetime is returned in case of success, -1.0 if failure.
bool Segments::TruncateAfter | ( | const double | gps | ) |
Remove segments after a given time.
gps | GPS time after which to truncate |
bool Segments::Write | ( | const string | txtfilename, | |
const int | ncol = 2 | |||
) |
Writes segments in a text file.
Several formats are supported and are defined by the 'ncol' argument. See the GWOLLUM convention for segments.
true is returned in case of success.
txtfilename | file name | |
ncol | number of columns in the segment file. See the GWOLLUM convention for supported formats. |
bool Segments::WriteROOT | ( | const string | rootfilename | ) |
Writes segments in a ROOT file.
The segments are saved in a TTree formatted as described in the GWOLLUM convention for segments.
true is returned in case of success.
rootfilename | output file name |
friend class ReadTriggerSegments [friend] |
double Segments::livetime_tot [protected] |
integrated livetime
vector<double> Segments::seg_end [protected] |
list of segment ends
vector<double> Segments::seg_start [protected] |
list of segment starts
vector<int> Segments::seg_tag [protected] |
list of segment tags
bool Segments::status [protected] |
status of the Segments object