Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
injection.hh
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////
2 // class for WaveBurst network event
3 // used as ROOT macro
4 // Sergey Klimenko, University of Florida
5 //////////////////////////////////////////////////////////
6 
7 
8 #ifndef injection_h
9 #define injection_h
10 
11 #include <TROOT.h>
12 #include <TChain.h>
13 #include <TFile.h>
14 #include "wseries.hh"
15 #include "network.hh"
16 #include "detector.hh"
17 #include "netcluster.hh"
18 
19 //class detector;
20 //class netcluster;
21 //class network;
22 
23 /* Structure of WaveBurst network event */
24 
25 #define INJECTION_INIT \
26  iFile(NULL),run(0),nevent(0),eventID(0),type(0),name(NULL),gps(0.),strain(0), \
27  psi(NULL),iota(NULL),phi(NULL),theta(NULL),bp(NULL),bx(NULL),time(NULL), \
28  duration(NULL),frequency(NULL),bandwidth(NULL),hrss(NULL),snr(NULL), \
29  Deff(NULL),mass(NULL),spin(NULL),pwf(NULL)
30 
31 class injection {
32 
33 //private :
34 
35 // detector* pD[NIFO_MAX]; //!temporary detectors pointers
36 
37 public :
38 
39  TFile *iFile; //!root input file cointainig the mdc TTree
40 
41  TTree *fChain; //!pointer to the analyzed TTree or TChain
42  Int_t fCurrent; //!current Tree number in a TChain
43  Int_t ndim; //! number of detectors
44 
45 //Declaration of leaves types
46 // for arrays: ifo1 - first index, ifo2 - second index, .....
47  Int_t run; // run ID
48  Int_t nevent; // event count
49  Int_t eventID; // event ID
50  Int_t type; // injection type
51  string* name; //! injection name
52 
53  Float_t factor; // simulation factor
54  Float_t distance; // distance to source in Mpc
55  Float_t mchirp; // chirp mass in Mo
56 
57  Float_t rp0; // eBBH binary distance
58  Float_t e0; // eBBH eccentricity
59  Float_t redshift; // eBBH redshift
60 
61  Double_t gps; // start time of data segment
62  Double_t strain; // strain of injected simulated signals
63  Float_t* psi; //! source psi angle
64  Float_t* iota; //! source iota angle
65  Float_t* phi; //! source phi angle
66  Float_t* theta; //! source theta angle
67  Float_t* bp; //! beam pattern coefficients for hp
68  Float_t* bx; //! beam pattern coefficients for hx
69 
70  Double_t* time; //! injection gps time
71  Float_t* duration; //! estimated duration
72 
73  Float_t* frequency; //! average center_of_hrss frequency
74  Float_t* bandwidth; //! estimated bandwidth
75 
76  Double_t* hrss; //! injected hrss in the detectors
77  Float_t* snr; //! injected snr in the detectors
78  Float_t* Deff; //! detector specific effective distance
79  Float_t* mass; //! [m1,m2], binary mass parameters
80  Float_t* spin; //! [x1,y1,z1,x2,y2,z2] components of spin vector
81 
82  wavearray<double>** pwf; //! pointer to the reconstructed waveform
83 
84 //List of branches
85 
86  TBranch *b_ndim; //!
87  TBranch *b_run; //!
88  TBranch *b_nevent; //!
89  TBranch *b_eventID; //!
90  TBranch *b_type; //!
91  TBranch *b_name; //!
92 
93  TBranch *b_factor; //!
94  TBranch *b_distance; //!
95  TBranch *b_mchirp; //!
96  TBranch *b_rp0; //!
97  TBranch *b_e0; //!
98  TBranch *b_redshift; //!
99  TBranch *b_gps; //!
100  TBranch *b_strain; //!
101  TBranch *b_psi; //!
102  TBranch *b_iota; //!
103  TBranch *b_phi; //!
104  TBranch *b_theta; //!
105  TBranch *b_bp; //!
106  TBranch *b_bx; //!
107 
108  TBranch *b_time; //!
109  TBranch *b_duration; //!
110 
111  TBranch *b_frequency; //!
112  TBranch *b_bandwidth; //!
113 
114  TBranch *b_hrss; //!
115  TBranch *b_snr; //!
116  TBranch *b_Deff; //!
117  TBranch *b_mass; //!
118  TBranch *b_spin; //!
119 
121  { ndim=1; allocate(); init(); return; }
122 
123  injection(int n) : INJECTION_INIT
124  { ndim=n; allocate(); init(); return; }
125 
127  { ndim=a.ndim; allocate(); init(); *this = a; return; }
128 
129  injection(TTree *tree, int n) : INJECTION_INIT
130  { ndim=n; allocate(); init(); if(tree) Init(tree); }
131 
133  { TTree* tree=Init(fName,n); allocate(); init(); if(tree) Init(tree); return; }
134 
135  virtual ~injection() {
136 
137  if (name) delete name; // injection name
138 
139  if (psi) free(psi); // polarization angles
140  if (iota) free(iota); // polarization angles
141  if (phi) free(phi); // source phi angles
142  if (theta) free(theta); // source theta angles
143  if (bp) free(bp); // beam pattern coefficients for hp
144  if (bx) free(bx); // beam pattern coefficients for hx
145 
146  if (time) free(time); // average center_of_snr time
147  if (duration) free(duration); // injection duration
148 
149  if (frequency) free(frequency); // average center_of_hrss frequency
150  if (bandwidth) free(bandwidth); // injection bandwidth
151  if (hrss) free(hrss); // injected hrss
152  if (snr) free(snr); // injected snr
153  if (Deff) free(Deff); // effective distance
154  if (mass) free(mass); // source mass vector
155  if (spin) free(spin); // source spin vector
156 
157  if (pwf) free(pwf); // pointer to the reconstructed waveform
158 
159  if(iFile) {if(fChain) delete fChain; delete iFile;}
160  };
161 
162  virtual injection& operator=(const injection &);
163 
164  Int_t GetEntry(Int_t);
165  Int_t GetEntries();
166  void allocate();
167  void init();
168  TTree* Init(TString fName, int n);
169  void Init(TTree *);
170  Bool_t Notify();
171  TTree* setTree();
172  Bool_t fill_in(network*,int,bool=true);
173  void output(TTree*, network*, double, bool=true);
174 
175 // void Loop();
176 // Int_t Cut(Int_t entry);
177 // Int_t LoadTree(Int_t entry);
178  void Show(Int_t entry = -1);
179 
180  // used by THtml doc
181  ClassDef(injection,3)
182 };
183 
184 #endif
185 
Int_t fCurrent
pointer to the analyzed TTree or TChain
Definition: injection.hh:42
TTree * tree
Definition: TimeSortTree.C:20
TBranch * b_nevent
Definition: injection.hh:88
TBranch * b_hrss
Definition: injection.hh:114
TTree * Init(TString fName, int n)
Definition: injection.cc:11
Double_t * time
beam pattern coefficients for hx
Definition: injection.hh:70
Float_t rp0
Definition: injection.hh:57
string * name
Definition: injection.hh:51
TBranch * b_psi
Definition: injection.hh:101
Float_t e0
Definition: injection.hh:58
TBranch * b_bandwidth
Definition: injection.hh:112
Int_t GetEntry(Int_t)
Definition: injection.cc:178
TFile * iFile
Definition: injection.hh:39
TBranch * b_distance
Definition: injection.hh:94
wavearray< double > a(hp.size())
int n
Definition: cwb_net.C:10
Bool_t Notify()
Definition: injection.cc:141
TString("c")
void output(TTree *, network *, double, bool=true)
Definition: injection.cc:573
Float_t * snr
injected hrss in the detectors
Definition: injection.hh:77
virtual ~injection()
Definition: injection.hh:135
TBranch * b_Deff
Definition: injection.hh:116
TBranch * b_name
Definition: injection.hh:91
injection() allocate()
Float_t * iota
source psi angle
Definition: injection.hh:64
TBranch * b_duration
Definition: injection.hh:109
Float_t * duration
injection gps time
Definition: injection.hh:71
Float_t * bandwidth
average center_of_hrss frequency
Definition: injection.hh:74
TTree * setTree()
Definition: injection.cc:201
Double_t strain
Definition: injection.hh:62
TBranch * b_snr
Definition: injection.hh:115
TBranch * b_theta
Definition: injection.hh:104
TBranch * b_bp
Definition: injection.hh:105
Float_t * theta
source phi angle
Definition: injection.hh:66
Double_t gps
Definition: injection.hh:61
Int_t ndim
current Tree number in a TChain
Definition: injection.hh:43
TBranch * b_redshift
Definition: injection.hh:98
TBranch * b_e0
Definition: injection.hh:97
TBranch * b_iota
Definition: injection.hh:102
TBranch * b_spin
Definition: injection.hh:118
Float_t redshift
Definition: injection.hh:59
Float_t * psi
Definition: injection.hh:63
void Show(Int_t entry=-1)
Definition: injection.cc:189
Float_t factor
injection name
Definition: injection.hh:53
TBranch * b_type
Definition: injection.hh:90
Double_t * hrss
estimated bandwidth
Definition: injection.hh:76
virtual injection & operator=(const injection &)
Definition: injection.cc:274
Float_t mchirp
Definition: injection.hh:55
TBranch * b_frequency
Definition: injection.hh:111
TBranch * b_mass
Definition: injection.hh:117
#define INJECTION_INIT
Definition: injection.hh:25
Int_t GetEntries()
Definition: injection.cc:183
Float_t * mass
detector specific effective distance
Definition: injection.hh:79
TBranch * b_phi
Definition: injection.hh:103
TBranch * b_rp0
Definition: injection.hh:96
TBranch * b_factor
Definition: injection.hh:93
Int_t run
number of detectors
Definition: injection.hh:47
TBranch * b_eventID
Definition: injection.hh:89
double * entry
Definition: cwb_setcuts.C:206
Float_t * phi
source iota angle
Definition: injection.hh:65
Int_t type
Definition: injection.hh:50
TBranch * b_run
Definition: injection.hh:87
Float_t * bx
beam pattern coefficients for hp
Definition: injection.hh:68
injection(TTree *tree, int n allocate)()
Definition: injection.hh:130
Float_t distance
Definition: injection.hh:54
TBranch * b_strain
Definition: injection.hh:100
TBranch * b_gps
Definition: injection.hh:99
Float_t * bp
source theta angle
Definition: injection.hh:67
TBranch * b_bx
Definition: injection.hh:106
Int_t eventID
Definition: injection.hh:49
Float_t * frequency
estimated duration
Definition: injection.hh:73
Bool_t fill_in(network *, int, bool=true)
Definition: injection.cc:325
char fName[256]
TBranch * b_mchirp
Definition: injection.hh:95
wavearray< double > ** pwf
[x1,y1,z1,x2,y2,z2] components of spin vector
Definition: injection.hh:82
TBranch * b_ndim
pointer to the reconstructed waveform
Definition: injection.hh:86
Float_t * spin
[m1,m2], binary mass parameters
Definition: injection.hh:80
TTree * fChain
root input file cointainig the mdc TTree
Definition: injection.hh:41
Int_t nevent
Definition: injection.hh:48
Float_t * Deff
injected snr in the detectors
Definition: injection.hh:78
TBranch * b_time
Definition: injection.hh:108