GWOLLUM 4.2.0
Tools for gravitational-wave analyses
Loading...
Searching...
No Matches
InjRea.h
Go to the documentation of this file.
1
6#ifndef __InjRea__
7#define __InjRea__
8
9#include <TChain.h>
10#include <TTreeIndex.h>
11#include <TH1.h>
12#include <TGraph.h>
13#include <TMath.h>
14#include "InjTyp.h"
15#include "Monitor.h"
16
17using namespace std;
18
24class InjRea{
25
26public:
27
40 InjRea(const string aPattern, const unsigned int aVerbose=0);
41
45 virtual ~InjRea();
53 inline Long64_t GetInjectionIndex(void){ return inj_index; };
54
61 Long64_t GetInjectionIndexAfter(const double aTime);
62
71 inline Long64_t GetFirstInjectionIndexEndAfter(const double aTime){ return GetInjectionIndexAfter(aTime - inj_durmax_p); };
72
81 inline Long64_t GetLastInjectionIndexStartBefore(const double aTime){ return GetInjectionIndexAfter(aTime + inj_durmax_m); };
82
89 int LoadInjection(const Long64_t aInjIndex);
90
95 double GetSineGaussh0plus(void);
96
101 double GetSineGaussh0cross(void);
102
106 inline unsigned int GetInjectionNameIndex(void){
107 unsigned int n = 0;
108 for(n=0; n<wave_names.size(); n++)
109 if(!(*wave_name).compare(wave_names[n])) break;
110 return n;
111 };
112
116 inline string GetInjectionName(void){ return *wave_name; };
117
121 inline unsigned int GetInjectionNameN(void){ return wave_names.size(); };
122
128 inline string GetInjectionName(const unsigned int aNameIndex){ return wave_names[aNameIndex]; };
129
133 inline injtype GetInjectionType(void){ return (injtype)wave_type; };
134
138 inline double GetInjectionTime(void){ return inj_time; };
139
144 inline double GetInjectionTimeMin(const unsigned int aNameIndex){
145 if(aNameIndex<wave_names.size()) return inj_timemin[aNameIndex];
146 return inj_timemin[wave_names.size()];
147 };
148
153 inline double GetInjectionTimeMax(const unsigned int aNameIndex){
154 if(aNameIndex<wave_names.size()) return inj_timemax[aNameIndex];
155 return inj_timemax[wave_names.size()];
156 };
157
161 inline double GetInjectionRa(void){ return inj_ra; };
162
166 inline double GetInjectionDec(void){ return inj_dec; };
167
171 inline double GetInjectionEccentricity(void){ return inj_ecc; };
172
176 inline double GetInjectionPolarization(void){ return inj_psi; };
177
181 inline double GetInjectionAmplitude(void){ return inj_amp; };
182
187 inline double GetInjectionAmplitudeMin(const unsigned int aNameIndex){
188 if(aNameIndex<wave_names.size()) return inj_ampmin[aNameIndex];
189 return inj_ampmin[wave_names.size()];
190 };
191
196 inline double GetInjectionAmplitudeMax(const unsigned int aNameIndex){
197 if(aNameIndex<wave_names.size()) return inj_ampmax[aNameIndex];
198 return inj_ampmax[wave_names.size()];
199 };
200
204 inline double GetInjectionSigma(void){ return inj_sigma; };
205
210 inline double GetInjectionSigmaMin(const unsigned int aNameIndex){
211 if(aNameIndex<wave_names.size()) return inj_sigmamin[aNameIndex];
212 return inj_sigmamin[wave_names.size()];
213 };
214
219 inline double GetInjectionSigmaMax(const unsigned int aNameIndex){
220 if(aNameIndex<wave_names.size()) return inj_sigmamax[aNameIndex];
221 return inj_sigmamax[wave_names.size()];
222 };
223
227 inline double GetInjectionFrequency(void){ return inj_f0; };
228
233 inline double GetInjectionFrequencyMin(const unsigned int aNameIndex){
234 if(aNameIndex<wave_names.size()) return inj_f0min[aNameIndex];
235 return inj_f0min[wave_names.size()];
236 };
237
242 inline double GetInjectionFrequencyMax(const unsigned int aNameIndex){
243 if(aNameIndex<wave_names.size()) return inj_f0max[aNameIndex];
244 return inj_f0max[wave_names.size()];
245 };
246
255 double GetInjectionTimeStart(void);
256
265 double GetInjectionTimeEnd(void);
266
270 inline bool GetInjectionTag(void){ return inj_tag[InjTree->GetReadEntry()]; };
271
276 inline bool GetInjectionTag(const Long64_t aInjIndex){ return inj_tag[sort_index->GetIndex()[aInjIndex]]; };
277
282 inline void SetInjectionTag(const bool aTag){ inj_tag[InjTree->GetReadEntry()]=aTag; };
283
289 inline void SetInjectionTag(const Long64_t aInjIndex, const bool aTag){ inj_tag[sort_index->GetIndex()[aInjIndex]]=aTag; };
290
295 inline Long64_t GetN(const unsigned int aNameIndex){
296 return wave_name_n[aNameIndex];
297 };
298
302 inline Long64_t GetN(void){ return InjTree->GetEntries(); };
303
312 TH1D* GetInjectionParamDist(const string aParamName, const unsigned int aNbins=1, const string aBinType="UNIFORM");
313
317 inline string GetInputFilePattern(void){ return pattern; };
318
319protected:
320
322 TGraph *wave_hplus;
323 TGraph *wave_hcross;
324 TChain *InjTree;
325
326private:
327
328 string pattern;
329
330 // INJECTION TREE
331 Long64_t inj_index;
332 TTreeIndex *sort_index;
333 double inj_ra;
334 double inj_dec;
335 double inj_psi;
336 double inj_time;
337 double inj_amp;
338 double inj_f0;
339 double inj_sigma;
340 double inj_ecc;
341
342 // MIN/MAX
343 double *inj_timemin;
344 double *inj_timemax;
345 double *inj_sigmamin;
346 double *inj_sigmamax;
347 double *inj_ampmin;
348 double *inj_ampmax;
349 double *inj_f0min;
350 double *inj_f0max;
351
352 bool *inj_tag;
355
356 // WAVEFORM TREE
357 TChain *WaveTree;
358 UInt_t wave_type;
359 string *wave_name;
360 vector<string> wave_names;
361 Long64_t *wave_name_n;
362
363 ClassDef(InjRea,0)
364};
365
366#endif
367
368
Manage injection types.
injtype
List of injection types.
Definition InjTyp.h:17
Process monitoring.
Read a set of injections.
Definition InjRea.h:24
double * inj_timemax
Time max.
Definition InjRea.h:344
double * inj_sigmamax
Sigma max.
Definition InjRea.h:346
virtual ~InjRea()
Destructor of the InjRea class.
Definition InjRea.cc:214
injtype GetInjectionType(void)
Returns the current injection type.
Definition InjRea.h:133
string GetInjectionName(const unsigned int aNameIndex)
Returns a given injection name in the list.
Definition InjRea.h:128
double inj_durmax_m
Maximum injection duration before inj_time.
Definition InjRea.h:353
bool GetInjectionTag(void)
Returns the current injection tag.
Definition InjRea.h:270
Long64_t GetN(const unsigned int aNameIndex)
Returns the number of injections of a given name.
Definition InjRea.h:295
double inj_f0
Injection frequency.
Definition InjRea.h:338
double GetInjectionTimeEnd(void)
Returns the GPS ending time of the injection [s] at the center of the Earth.
Definition InjRea.cc:344
double inj_amp
Injection amplitude.
Definition InjRea.h:337
double GetInjectionTimeMax(const unsigned int aNameIndex)
Returns the maximum injection GPS time.
Definition InjRea.h:153
double GetInjectionTime(void)
Returns the current injection GPS time [s].
Definition InjRea.h:138
Long64_t inj_index
Current injection index.
Definition InjRea.h:331
TTreeIndex * sort_index
TTree index.
Definition InjRea.h:332
string GetInjectionName(void)
Returns the current injection name.
Definition InjRea.h:116
void SetInjectionTag(const bool aTag)
Sets a new tag value to the current injection.
Definition InjRea.h:282
double GetInjectionPolarization(void)
Returns the current injection polarization [rad].
Definition InjRea.h:176
double * inj_f0min
Frequency min.
Definition InjRea.h:349
double GetInjectionSigma(void)
Returns the current injection sigma [s].
Definition InjRea.h:204
int LoadInjection(const Long64_t aInjIndex)
Loads a given injection.
Definition InjRea.cc:267
double inj_durmax_p
Maximum injection duration after inj_time.
Definition InjRea.h:354
double GetInjectionRa(void)
Returns the current injection right-ascension [rad].
Definition InjRea.h:161
TH1D * GetInjectionParamDist(const string aParamName, const unsigned int aNbins=1, const string aBinType="UNIFORM")
Returns the 1D distribution of an injection parameter.
Definition InjRea.cc:361
double * inj_ampmin
Amplitude min.
Definition InjRea.h:347
bool GetInjectionTag(const Long64_t aInjIndex)
Returns an injection tag.
Definition InjRea.h:276
double GetInjectionSigmaMin(const unsigned int aNameIndex)
Returns the minimum injection sigma [s].
Definition InjRea.h:210
TChain * InjTree
Injection tree.
Definition InjRea.h:324
double GetInjectionDec(void)
Returns the current injection declination [rad].
Definition InjRea.h:166
double GetInjectionFrequencyMax(const unsigned int aNameIndex)
Returns the maximum injection frequency [Hz].
Definition InjRea.h:242
TChain * WaveTree
Waveform tree.
Definition InjRea.h:357
Monitor * mon
Class monitor.
Definition InjRea.h:321
Long64_t GetFirstInjectionIndexEndAfter(const double aTime)
Returns the index of the first injection ending after a given time.
Definition InjRea.h:71
double GetInjectionSigmaMax(const unsigned int aNameIndex)
Returns the maximum injection sigma [s].
Definition InjRea.h:219
Long64_t * wave_name_n
Definition InjRea.h:361
double GetInjectionTimeStart(void)
Returns the GPS starting time of the injection [s] at the center of the Earth.
Definition InjRea.cc:327
double GetSineGaussh0plus(void)
Returns the SineGauss peak amplitude of current injection.
Definition InjRea.cc:285
Long64_t GetN(void)
Returns the number of injections.
Definition InjRea.h:302
TGraph * wave_hcross
Waveform .
Definition InjRea.h:323
double * inj_f0max
Frequency max.
Definition InjRea.h:350
double GetInjectionEccentricity(void)
Returns the current injection eccentricity.
Definition InjRea.h:171
Long64_t GetLastInjectionIndexStartBefore(const double aTime)
Returns the index of the last injection starting before a given time.
Definition InjRea.h:81
TGraph * wave_hplus
Waveform .
Definition InjRea.h:322
double * inj_sigmamin
Sigma min.
Definition InjRea.h:345
double * inj_ampmax
Amplitude max.
Definition InjRea.h:348
double GetSineGaussh0cross(void)
Returns the SineGauss peak amplitude of current injection.
Definition InjRea.cc:306
double GetInjectionAmplitudeMin(const unsigned int aNameIndex)
Returns the minimum injection amplitude.
Definition InjRea.h:187
string pattern
Input file pattern.
Definition InjRea.h:328
string * wave_name
Waveform name.
Definition InjRea.h:359
double * inj_timemin
Time min.
Definition InjRea.h:343
UInt_t wave_type
Waveform type: see InjTyp.h.
Definition InjRea.h:358
bool * inj_tag
Injection tags.
Definition InjRea.h:352
void SetInjectionTag(const Long64_t aInjIndex, const bool aTag)
Sets a new tag value to an injection.
Definition InjRea.h:289
Long64_t GetInjectionIndex(void)
Returns the index of the current injection.
Definition InjRea.h:53
double GetInjectionAmplitudeMax(const unsigned int aNameIndex)
Returns the maximum injection amplitude.
Definition InjRea.h:196
double inj_psi
Injection polarization angle.
Definition InjRea.h:335
unsigned int GetInjectionNameIndex(void)
Returns the current injection name index.
Definition InjRea.h:106
double GetInjectionFrequencyMin(const unsigned int aNameIndex)
Returns the minimum injection frequency [Hz].
Definition InjRea.h:233
double inj_dec
Injection declination.
Definition InjRea.h:334
vector< string > wave_names
List of waveform names.
Definition InjRea.h:360
double GetInjectionFrequency(void)
Returns the current injection frequency [Hz].
Definition InjRea.h:227
string GetInputFilePattern(void)
Returns the input file pattern provided in the constructor.
Definition InjRea.h:317
double inj_sigma
Injection sigma.
Definition InjRea.h:339
double inj_ra
Injection right ascension.
Definition InjRea.h:333
double GetInjectionTimeMin(const unsigned int aNameIndex)
Returns the minimum injection GPS time.
Definition InjRea.h:144
double inj_ecc
Injection eccentricity.
Definition InjRea.h:340
double GetInjectionAmplitude(void)
Returns the current injection amplitude.
Definition InjRea.h:181
Long64_t GetInjectionIndexAfter(const double aTime)
Returns the index of the first injection after a given time.
Definition InjRea.cc:237
double inj_time
Injection time.
Definition InjRea.h:336
unsigned int GetInjectionNameN(void)
Returns the number of injection names.
Definition InjRea.h:121
Monitor a GWOLLUM processing.
Definition Monitor.h:39