GWOLLUM  4.2.0
Tools for gravitational-wave analyses
Segments.h
Go to the documentation of this file.
1 
6 #ifndef __Segments__
7 #define __Segments__
8 
9 #include "ReadAscii.h"
10 #include <TFile.h>
11 #include <TChain.h>
12 #include <TH1.h>
13 
14 using namespace std;
15 
28 class Segments {
29 
30  public:
31 
40  Segments(void);
41 
56  Segments(const string aSegmentFile,
57  const int aDefaultTag=1);
58 
66  Segments(TTree *aSegmentTree);
67 
77  Segments(const string aFilePattern, const string aDirectory, const string aTreeName);
78 
87  Segments(const vector<double> aGpsStarts, const vector<double> aGpsEnds);
88 
96  Segments(const double aGpsStart, const double aGpsEnd, const int aTag=1);
97 
101  virtual ~Segments(void);
110  void Reset(void);
111 
115  inline bool GetStatus(void) const { return status; };
116 
122  int GetSegmentStartBefore(const double aGps);
123 
127  inline unsigned int GetN(void){ return seg_start.size(); };
128 
133  unsigned int GetNWithTag(const int aTag);
134 
138  inline vector<double> GetStarts(void){ return seg_start; };
139 
143  inline vector<double> GetEnds(void){ return seg_end; };
144 
148  inline vector<int> GetTags(void){ return seg_tag; };
149 
155  inline double GetStart(const unsigned int aSegmentIndex=0){ return seg_start[aSegmentIndex]; };
156 
162  inline double GetEnd(const unsigned int aSegmentIndex=0){ return seg_end[aSegmentIndex]; };
163 
168  inline double GetFirst(){ return GetStart(0); };
169 
174  inline double GetLast(){ return GetEnd(GetN()-1); };
175 
181  inline double GetDuration(const unsigned int aSegmentIndex=0){ return seg_end[aSegmentIndex]-seg_start[aSegmentIndex]; };
182 
188  inline int GetTag(const unsigned int aSegmentIndex=0){ return seg_tag[aSegmentIndex]; };
189 
194  inline double GetLiveTime(void){ return seg_livetime; };
195 
201  double GetLiveTime(const double aGpsStart, const double aGpsEnd);
202 
209  inline double GetLiveTime(const unsigned int aSegmentIndex){ return GetDuration(aSegmentIndex); };
210 
215  double GetLiveTimeWithTag(const int aTag);
216 
229  TH1I* GetHistogram(const double aGpsStart, const double aGpsEnd, const double aOffset=0.0);
230 
240  inline TH1I* GetHisto(const double aOffset=0.0){
241  return GetHistogram(seg_start[0],seg_end[seg_start.size()-1],aOffset);
242  };
243 
253  TTree* GetTree(const double aGpsStart, const double aGpsEnd);
254 
261  TTree* GetTree(void){
262  if(GetN()) return GetTree(seg_start[0], seg_end[seg_start.size()-1]);
263  return GetTree(0.0, 0.0);
264  };
265 
272  bool IsInsideSegment(unsigned int &aSegmentIndex, const double aGps);
273 
279  bool IsInsideSegment(const double aGps);
280 
287  inline void SetTag(const unsigned int aSegmentIndex, const int aTag){
288  seg_tag[aSegmentIndex] = aTag;
289  };
290 
295  void SetTags(const int aTag);
296 
306  bool Append(const double aGpsStart, const double aGpsEnd, const int aTag=1);
307 
317  bool Append(vector <double> aGpsStarts, vector <double> aGpsEnds, vector <double> aTags);
318 
326  bool Append(Segments *aSeg);
327 
335  bool Append(TTree *aSegmentTree);
336 
345  bool AddSegment(const double aGpsStart, const double aGpsEnd, const int aTag=1);
346 
353  bool AddSegments(Segments *aSeg);
354 
363  bool ApplyPadding(const double aPadStart, const double aPadEnd);
364 
370  bool TruncateAfter(const double aGps);
371 
377  bool TruncateBefore(const double aGps);
378 
385  bool MakeSecond(void);
386 
395  bool Reverse(void);
396 
405  bool Intersect(Segments *aSeg);
406 
417  bool Intersect(const double aGpsStart, const double aGpsEnd, const int aTag=1);
418 
432  void Dump(unsigned int aNcols=2, const string aTxtFileName="");
433 
434  private:
435 
436  bool status;
437  vector<double> seg_start;
438  vector<double> seg_end;
439  vector<int> seg_tag;
440  double seg_livetime;
441  vector <unsigned int> time_mark;
442  vector <unsigned int> n_mark;
443 
452  inline int GetTagProduct(const int aTag1, const int aTag2){
453  if(aTag1<0||aTag2<0) return TMath::Min(aTag1,aTag2);
454  return TMath::Max(aTag1,aTag2);
455  };
456 
469  bool CheckConsistency(const unsigned int aStartIndex);
470 
477  void MergeOverlaps(const unsigned int aStartIndex);
478 
483  void SetMarks(void);
484 
490  int GetMark(const unsigned int aGps);
491 
496  void ComputeLiveTime(void);
497 
498 
499  ClassDef(Segments,0)
500 };
501 
502 #endif
503 
504 
Parser for text files with columns.
Manage time segment lists.
Definition: Segments.h:28
unsigned int GetN(void)
Returns the number of time segments in the list.
Definition: Segments.h:127
vector< double > seg_start
List of segment starts.
Definition: Segments.h:437
vector< int > seg_tag
List of segment tags.
Definition: Segments.h:439
double GetEnd(const unsigned int aSegmentIndex=0)
Returns the ending time of a given segment.
Definition: Segments.h:162
void SetTag(const unsigned int aSegmentIndex, const int aTag)
Tags a given segment.
Definition: Segments.h:287
vector< double > GetStarts(void)
Returns the list of segment starts.
Definition: Segments.h:138
double GetLast()
Returns the ending time of the last segment.
Definition: Segments.h:174
vector< double > seg_end
List of segment ends.
Definition: Segments.h:438
TTree * GetTree(void)
Returns a formatted TTree representing the segments.
Definition: Segments.h:261
bool status
Status of the Segments object.
Definition: Segments.h:436
double GetFirst()
Returns the starting time of the first segment.
Definition: Segments.h:168
double GetStart(const unsigned int aSegmentIndex=0)
Returns the starting time of a given segment.
Definition: Segments.h:155
vector< unsigned int > n_mark
Segment index marks.
Definition: Segments.h:442
double GetDuration(const unsigned int aSegmentIndex=0)
Returns the duration of a given segment.
Definition: Segments.h:181
double GetLiveTime(const unsigned int aSegmentIndex)
Returns the livetime of one segment.
Definition: Segments.h:209
TH1I * GetHisto(const double aOffset=0.0)
Returns an histogram representing the segments.
Definition: Segments.h:240
Segments(void)
Default constructor of the Segments class.
int GetTag(const unsigned int aSegmentIndex=0)
Returns the tag of a given segment.
Definition: Segments.h:188
bool GetStatus(void) const
Returns the status of the Segments object.
Definition: Segments.h:115
double GetLiveTime(void)
Returns the total livetime.
Definition: Segments.h:194
vector< int > GetTags(void)
Returns the list of segment tags.
Definition: Segments.h:148
double seg_livetime
Integrated livetime.
Definition: Segments.h:440
vector< unsigned int > time_mark
Time marks.
Definition: Segments.h:441
int GetTagProduct(const int aTag1, const int aTag2)
Computes the product of 2 tags.
Definition: Segments.h:452
vector< double > GetEnds(void)
Returns the list of segment ends.
Definition: Segments.h:143