GWOLLUM 4.2.0
Tools for gravitational-wave analyses
Loading...
Searching...
No Matches
Segments.h
Go to the documentation of this file.
1
6#ifndef __Segments__
7#define __Segments__
8
9#include "ReadAscii.h"
10#include "Monitor.h"
11#include <TFile.h>
12#include <TChain.h>
13#include <TH1.h>
14
15using namespace std;
16
30
31 public:
32
41 Segments(void);
42
57 Segments(const string aSegmentFile,
58 const int aDefaultTag=1);
59
67 Segments(TTree *aSegmentTree);
68
78 Segments(const string aFilePattern, const string aDirectory, const string aTreeName);
79
88 Segments(const vector<double> aGpsStarts, const vector<double> aGpsEnds);
89
97 Segments(const double aGpsStart, const double aGpsEnd, const int aTag=1);
98
102 virtual ~Segments(void);
110 inline bool GetStatus(void){ return mon->GetStatus(); };
111
116 void Reset(void);
117
123 int GetSegmentStartBefore(const double aGps);
124
128 inline unsigned int GetN(void){ return seg_start.size(); };
129
134 unsigned int GetNWithTag(const int aTag);
135
139 inline vector<double> GetStarts(void){ return seg_start; };
140
144 inline vector<double> GetEnds(void){ return seg_end; };
145
149 inline vector<int> GetTags(void){ return seg_tag; };
150
156 inline double GetStart(const unsigned int aSegmentIndex=0){ return seg_start[aSegmentIndex]; };
157
163 inline double GetEnd(const unsigned int aSegmentIndex=0){ return seg_end[aSegmentIndex]; };
164
169 inline double GetFirst(){ return GetStart(0); };
170
175 inline double GetLast(){ return GetEnd(GetN()-1); };
176
182 inline double GetDuration(const unsigned int aSegmentIndex=0){ return seg_end[aSegmentIndex]-seg_start[aSegmentIndex]; };
183
189 inline int GetTag(const unsigned int aSegmentIndex=0){ return seg_tag[aSegmentIndex]; };
190
195 inline double GetLiveTime(void){ return seg_livetime; };
196
202 double GetLiveTime(const double aGpsStart, const double aGpsEnd);
203
210 inline double GetLiveTime(const unsigned int aSegmentIndex){ return GetDuration(aSegmentIndex); };
211
216 double GetLiveTimeWithTag(const int aTag);
217
230 TH1I* GetHistogram(const double aGpsStart, const double aGpsEnd, const double aOffset=0.0);
231
241 inline TH1I* GetHisto(const double aOffset=0.0){
242 return GetHistogram(seg_start[0],seg_end[seg_start.size()-1],aOffset);
243 };
244
254 TTree* GetTree(const double aGpsStart, const double aGpsEnd);
255
262 TTree* GetTree(void){
263 if(GetN()) return GetTree(seg_start[0], seg_end[seg_start.size()-1]);
264 return GetTree(0.0, 0.0);
265 };
266
273 bool IsInsideSegment(unsigned int &aSegmentIndex, const double aGps);
274
280 bool IsInsideSegment(const double aGps);
281
288 inline void SetTag(const unsigned int aSegmentIndex, const int aTag){
289 seg_tag[aSegmentIndex] = aTag;
290 };
291
296 void SetTags(const int aTag);
297
307 bool Append(const double aGpsStart, const double aGpsEnd, const int aTag=1);
308
318 bool Append(vector <double> aGpsStarts, vector <double> aGpsEnds, vector <double> aTags);
319
327 bool Append(Segments *aSeg);
328
336 bool Append(TTree *aSegmentTree);
337
346 bool AddSegment(const double aGpsStart, const double aGpsEnd, const int aTag=1);
347
354 bool AddSegments(Segments *aSeg);
355
364 bool ApplyPadding(const double aPadStart, const double aPadEnd);
365
371 bool TruncateAfter(const double aGps);
372
378 bool TruncateBefore(const double aGps);
379
386 bool MakeSecond(void);
387
396 bool Reverse(void);
397
406 bool Intersect(Segments *aSeg);
407
418 bool Intersect(const double aGpsStart, const double aGpsEnd, const int aTag=1);
419
433 void Dump(unsigned int aNcols=2, const string aTxtFileName="");
434
435 private:
436
438 vector<double> seg_start;
439 vector<double> seg_end;
440 vector<int> seg_tag;
442 vector <unsigned int> time_mark;
443 vector <unsigned int> n_mark;
444
453 inline int GetTagProduct(const int aTag1, const int aTag2){
454 if(aTag1<0||aTag2<0) return TMath::Min(aTag1,aTag2);
455 return TMath::Max(aTag1,aTag2);
456 };
457
470 bool CheckConsistency(const unsigned int aStartIndex);
471
478 void MergeOverlaps(const unsigned int aStartIndex);
479
484 void SetMarks(void);
485
491 int GetMark(const unsigned int aGps);
492
497 void ComputeLiveTime(void);
498
499
500 ClassDef(Segments,0)
501};
502
503#endif
504
505
Process monitoring.
Parser for text files with columns.
Monitor a GWOLLUM processing.
Definition Monitor.h:39
bool GetStatus(void)
Returns the object status.
Definition Monitor.h:109
Manage time segment lists.
Definition Segments.h:29
bool TruncateBefore(const double aGps)
Removes segments before a given time.
Definition Segments.cc:826
bool Intersect(Segments *aSeg)
Intersects with a segment list.
Definition Segments.cc:939
vector< double > GetEnds(void)
Returns the list of segment ends.
Definition Segments.h:144
bool IsInsideSegment(unsigned int &aSegmentIndex, const double aGps)
Checks whether a given GPS time is inside a segment of the list.
Definition Segments.cc:466
unsigned int GetN(void)
Returns the number of time segments in the list.
Definition Segments.h:128
bool CheckConsistency(const unsigned int aStartIndex)
Checks the sanity of the segment list.
Definition Segments.cc:1124
bool MakeSecond(void)
Rounds segment boundaries to integer seconds.
Definition Segments.cc:867
unsigned int GetNWithTag(const int aTag)
Returns the number of time segments in the list with a given tag.
Definition Segments.cc:303
vector< double > seg_start
List of segment starts.
Definition Segments.h:438
bool GetStatus(void)
Returns the class status.
Definition Segments.h:110
vector< int > seg_tag
List of segment tags.
Definition Segments.h:440
bool Reverse(void)
Reverses the segment list.
Definition Segments.cc:890
TH1I * GetHisto(const double aOffset=0.0)
Returns an histogram representing the segments.
Definition Segments.h:241
TTree * GetTree(void)
Returns a formatted TTree representing the segments.
Definition Segments.h:262
void ComputeLiveTime(void)
Computes the sgments livetime.
Definition Segments.cc:1244
double GetEnd(const unsigned int aSegmentIndex=0)
Returns the ending time of a given segment.
Definition Segments.h:163
void SetTag(const unsigned int aSegmentIndex, const int aTag)
Tags a given segment.
Definition Segments.h:288
TH1I * GetHistogram(const double aGpsStart, const double aGpsEnd, const double aOffset=0.0)
Returns an histogram representing the segments.
Definition Segments.cc:373
void MergeOverlaps(const unsigned int aStartIndex)
Merges segment overlaps.
Definition Segments.cc:1175
bool AddSegments(Segments *aSeg)
Adds a segment list to the list.
Definition Segments.cc:722
double GetLast()
Returns the ending time of the last segment.
Definition Segments.h:175
void SetMarks(void)
Sets time marks for internal lookup.
Definition Segments.cc:1200
bool AddSegment(const double aGpsStart, const double aGpsEnd, const int aTag=1)
Adds a segment to the list.
Definition Segments.cc:681
double GetLiveTimeWithTag(const int aTag)
Returns the livetime for segments with a given tag.
Definition Segments.cc:358
vector< double > seg_end
List of segment ends.
Definition Segments.h:439
bool TruncateAfter(const double aGps)
Removes segments after a given time.
Definition Segments.cc:786
bool Append(const double aGpsStart, const double aGpsEnd, const int aTag=1)
Appends a segment to the list.
Definition Segments.cc:511
bool ApplyPadding(const double aPadStart, const double aPadEnd)
Apply a padding to all segments.
Definition Segments.cc:746
int GetSegmentStartBefore(const double aGps)
Returns the segment starting just before a given GPS time.
Definition Segments.cc:280
int GetMark(const unsigned int aGps)
Returns the time mark for a given GPS time.
Definition Segments.cc:1226
double GetFirst()
Returns the starting time of the first segment.
Definition Segments.h:169
double GetStart(const unsigned int aSegmentIndex=0)
Returns the starting time of a given segment.
Definition Segments.h:156
vector< unsigned int > n_mark
Segment index marks.
Definition Segments.h:443
double GetDuration(const unsigned int aSegmentIndex=0)
Returns the duration of a given segment.
Definition Segments.h:182
vector< double > GetStarts(void)
Returns the list of segment starts.
Definition Segments.h:139
double GetLiveTime(const unsigned int aSegmentIndex)
Returns the livetime of one segment.
Definition Segments.h:210
Segments(void)
Default constructor of the Segments class.
Definition Segments.cc:21
int GetTag(const unsigned int aSegmentIndex=0)
Returns the tag of a given segment.
Definition Segments.h:189
void Reset(void)
Resets the segment list.
Definition Segments.cc:267
void Dump(unsigned int aNcols=2, const string aTxtFileName="")
Dumps the current segment list in the standard output or in a file.
Definition Segments.cc:1041
vector< int > GetTags(void)
Returns the list of segment tags.
Definition Segments.h:149
double GetLiveTime(void)
Returns the total livetime.
Definition Segments.h:195
Monitor * mon
Class monitor.
Definition Segments.h:437
double seg_livetime
Integrated livetime.
Definition Segments.h:441
vector< unsigned int > time_mark
Time marks.
Definition Segments.h:442
int GetTagProduct(const int aTag1, const int aTag2)
Computes the product of 2 tags.
Definition Segments.h:453
virtual ~Segments(void)
Destructor of the Segments class.
Definition Segments.cc:260
void SetTags(const int aTag)
Tags all segments.
Definition Segments.cc:503