Logo coherent WaveBurst  
Library Reference Guide
Logo
injection.hh
Go to the documentation of this file.
1 /*
2 # Copyright (C) 2019 Sergey Klimenko, Gabriele Vedovato
3 #
4 # This program is free software: you can redistribute it and/or modify
5 # it under the terms of the GNU General Public License as published by
6 # the Free Software Foundation, either version 3 of the License, or
7 # (at your option) any later version.
8 #
9 # This program is distributed in the hope that it will be useful,
10 # but WITHOUT ANY WARRANTY; without even the implied warranty of
11 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 # GNU General Public License for more details.
13 #
14 # You should have received a copy of the GNU General Public License
15 # along with this program. If not, see <https://www.gnu.org/licenses/>.
16 */
17 
18 
19 //////////////////////////////////////////////////////////
20 // class for WaveBurst network event
21 // used as ROOT macro
22 // Sergey Klimenko, University of Florida
23 //////////////////////////////////////////////////////////
24 
25 
26 #ifndef injection_h
27 #define injection_h
28 
29 #include <TROOT.h>
30 #include <TChain.h>
31 #include <TFile.h>
32 #include "wseries.hh"
33 #include "network.hh"
34 #include "detector.hh"
35 #include "netcluster.hh"
36 
37 //class detector;
38 //class netcluster;
39 //class network;
40 
41 /* Structure of WaveBurst network event */
42 
43 #define INJECTION_INIT \
44  iFile(NULL),run(0),nevent(0),eventID(0),type(0),name(NULL),log(NULL),gps(0.),strain(0), \
45  psi(NULL),iota(NULL),phi(NULL),theta(NULL),bp(NULL),bx(NULL),time(NULL), \
46  duration(NULL),frequency(NULL),bandwidth(NULL),hrss(NULL),snr(NULL), \
47  Deff(NULL),mass(NULL),spin(NULL),pwf(NULL)
48 
49 class injection {
50 
51 //private :
52 
53 // detector* pD[NIFO_MAX]; //!temporary detectors pointers
54 
55 public :
56 
57  TFile *iFile; //!root input file cointainig the mdc TTree
58 
59  TTree *fChain; //!pointer to the analyzed TTree or TChain
60  Int_t fCurrent; //!current Tree number in a TChain
61  Int_t ndim; //! number of detectors
62 
63 //Declaration of leaves types
64 // for arrays: ifo1 - first index, ifo2 - second index, .....
65  Int_t run; // run ID
66  Int_t nevent; // event count
67  Int_t eventID; // event ID
68  Int_t type; // injection type
69  string* name; //! injection name
70  string* log; //! injection log
71 
72  Float_t factor; // simulation factor
73  Float_t distance; // distance to source in Mpc
74  Float_t mchirp; // chirp mass in Mo
75 
76  Float_t rp0; // eBBH binary distance
77  Float_t e0; // eBBH eccentricity
78  Float_t redshift; // eBBH redshift
79 
80  Double_t gps; // start time of data segment
81  Double_t strain; // strain of injected simulated signals
82  Float_t* psi; //! source psi angle
83  Float_t* iota; //! source iota angle
84  Float_t* phi; //! source phi angle
85  Float_t* theta; //! source theta angle
86  Float_t* bp; //! beam pattern coefficients for hp
87  Float_t* bx; //! beam pattern coefficients for hx
88 
89  Double_t* time; //! injection gps time
90  Float_t* duration; //! estimated duration
91 
92  Float_t* frequency; //! average center_of_hrss frequency
93  Float_t* bandwidth; //! estimated bandwidth
94 
95  Double_t* hrss; //! injected hrss in the detectors
96  Float_t* snr; //! injected snr in the detectors
97  Float_t* Deff; //! detector specific effective distance
98  Float_t* mass; //! [m1,m2], binary mass parameters
99  Float_t* spin; //! [x1,y1,z1,x2,y2,z2] components of spin vector
100 
101  wavearray<double>** pwf; //! pointer to the reconstructed waveform
102 
103 //List of branches
104 
105  TBranch *b_ndim; //!
106  TBranch *b_run; //!
107  TBranch *b_nevent; //!
108  TBranch *b_eventID; //!
109  TBranch *b_type; //!
110  TBranch *b_name; //!
111  TBranch *b_log; //!
112 
113  TBranch *b_factor; //!
114  TBranch *b_distance; //!
115  TBranch *b_mchirp; //!
116  TBranch *b_rp0; //!
117  TBranch *b_e0; //!
118  TBranch *b_redshift; //!
119  TBranch *b_gps; //!
120  TBranch *b_strain; //!
121  TBranch *b_psi; //!
122  TBranch *b_iota; //!
123  TBranch *b_phi; //!
124  TBranch *b_theta; //!
125  TBranch *b_bp; //!
126  TBranch *b_bx; //!
127 
128  TBranch *b_time; //!
129  TBranch *b_duration; //!
130 
131  TBranch *b_frequency; //!
132  TBranch *b_bandwidth; //!
133 
134  TBranch *b_hrss; //!
135  TBranch *b_snr; //!
136  TBranch *b_Deff; //!
137  TBranch *b_mass; //!
138  TBranch *b_spin; //!
139 
141  { ndim=1; allocate(); init(); return; }
142 
143  injection(int n) : INJECTION_INIT
144  { ndim=n; allocate(); init(); return; }
145 
147  { ndim=a.ndim; allocate(); init(); *this = a; return; }
148 
149  injection(TTree *tree, int n) : INJECTION_INIT
150  { ndim=n; allocate(); init(); if(tree) Init(tree); }
151 
153  { TTree* tree=Init(fName,n); allocate(); init(); if(tree) Init(tree); return; }
154 
155  virtual ~injection() {
156 
157  if (name) delete name; // injection name
158  if (log) delete log; // injection log
159 
160  if (psi) free(psi); // polarization angles
161  if (iota) free(iota); // polarization angles
162  if (phi) free(phi); // source phi angles
163  if (theta) free(theta); // source theta angles
164  if (bp) free(bp); // beam pattern coefficients for hp
165  if (bx) free(bx); // beam pattern coefficients for hx
166 
167  if (time) free(time); // average center_of_snr time
168  if (duration) free(duration); // injection duration
169 
170  if (frequency) free(frequency); // average center_of_hrss frequency
171  if (bandwidth) free(bandwidth); // injection bandwidth
172  if (hrss) free(hrss); // injected hrss
173  if (snr) free(snr); // injected snr
174  if (Deff) free(Deff); // effective distance
175  if (mass) free(mass); // source mass vector
176  if (spin) free(spin); // source spin vector
177 
178  if (pwf) free(pwf); // pointer to the reconstructed waveform
179 
180  if(iFile) {if(fChain) delete fChain; delete iFile;}
181  };
182 
183  virtual injection& operator=(const injection &);
184 
185  Int_t GetEntry(Int_t);
186  Int_t GetEntries();
187  void allocate();
188  void init();
189  TTree* Init(TString fName, int n);
190  void Init(TTree *);
191  Bool_t Notify();
192  TTree* setTree();
193  Bool_t fill_in(network*,int,bool=true);
194  void output(TTree*, network*, double, bool=true);
195 
196 // void Loop();
197 // Int_t Cut(Int_t entry);
198 // Int_t LoadTree(Int_t entry);
199  void Show(Int_t entry = -1);
200 
201  // used by THtml doc
202  ClassDef(injection,4)
203 };
204 
205 #endif
206 
Int_t fCurrent
pointer to the analyzed TTree or TChain
Definition: injection.hh:60
TTree * tree
Definition: TimeSortTree.C:20
TBranch * b_nevent
Definition: injection.hh:107
TBranch * b_hrss
Definition: injection.hh:134
TTree * Init(TString fName, int n)
Definition: injection.cc:29
Double_t * time
beam pattern coefficients for hx
Definition: injection.hh:89
Float_t rp0
Definition: injection.hh:76
string * name
Definition: injection.hh:69
TBranch * b_psi
Definition: injection.hh:121
Float_t e0
Definition: injection.hh:77
TBranch * b_bandwidth
Definition: injection.hh:132
Int_t GetEntry(Int_t)
Definition: injection.cc:200
TFile * iFile
Definition: injection.hh:57
void init()
Definition: injection.cc:156
TBranch * b_distance
Definition: injection.hh:114
wavearray< double > a(hp.size())
int n
Definition: cwb_net.C:28
Bool_t Notify()
Definition: injection.cc:162
TString("c")
void output(TTree *, network *, double, bool=true)
Definition: injection.cc:602
Float_t * snr
injected hrss in the detectors
Definition: injection.hh:96
TBranch * b_Deff
Definition: injection.hh:136
TBranch * b_name
Definition: injection.hh:110
Float_t * iota
source psi angle
Definition: injection.hh:83
TBranch * b_duration
Definition: injection.hh:129
Float_t * duration
injection gps time
Definition: injection.hh:90
Float_t * bandwidth
average center_of_hrss frequency
Definition: injection.hh:93
TTree * setTree()
Definition: injection.cc:223
Double_t strain
Definition: injection.hh:81
TBranch * b_snr
Definition: injection.hh:135
TBranch * b_theta
Definition: injection.hh:124
TBranch * b_bp
Definition: injection.hh:125
Float_t * theta
source phi angle
Definition: injection.hh:85
Double_t gps
Definition: injection.hh:80
Int_t ndim
current Tree number in a TChain
Definition: injection.hh:61
TBranch * b_redshift
Definition: injection.hh:118
TBranch * b_e0
Definition: injection.hh:117
TBranch * b_iota
Definition: injection.hh:122
TBranch * b_spin
Definition: injection.hh:138
Float_t redshift
Definition: injection.hh:78
Float_t * psi
Definition: injection.hh:82
void Show(Int_t entry=-1)
Definition: injection.cc:211
Float_t factor
injection log
Definition: injection.hh:72
TBranch * b_type
Definition: injection.hh:109
Double_t * hrss
estimated bandwidth
Definition: injection.hh:95
virtual injection & operator=(const injection &)
Definition: injection.cc:297
Float_t mchirp
Definition: injection.hh:74
TBranch * b_frequency
Definition: injection.hh:131
TBranch * b_mass
Definition: injection.hh:137
#define INJECTION_INIT
Definition: injection.hh:43
Int_t GetEntries()
Definition: injection.cc:205
Float_t * mass
detector specific effective distance
Definition: injection.hh:98
TBranch * b_phi
Definition: injection.hh:123
TBranch * b_rp0
Definition: injection.hh:116
TBranch * b_factor
Definition: injection.hh:113
Int_t run
number of detectors
Definition: injection.hh:65
TBranch * b_eventID
Definition: injection.hh:108
double * entry
Definition: cwb_setcuts.C:224
Float_t * phi
source iota angle
Definition: injection.hh:84
Int_t type
Definition: injection.hh:68
void allocate()
Definition: injection.cc:109
TBranch * b_run
Definition: injection.hh:106
Float_t * bx
beam pattern coefficients for hp
Definition: injection.hh:87
TBranch * b_log
Definition: injection.hh:111
Float_t distance
Definition: injection.hh:73
TBranch * b_strain
Definition: injection.hh:120
string * log
injection name
Definition: injection.hh:70
TBranch * b_gps
Definition: injection.hh:119
Float_t * bp
source theta angle
Definition: injection.hh:86
TBranch * b_bx
Definition: injection.hh:126
Int_t eventID
Definition: injection.hh:67
Float_t * frequency
estimated duration
Definition: injection.hh:92
Bool_t fill_in(network *, int, bool=true)
Definition: injection.cc:349
char fName[256]
TBranch * b_mchirp
Definition: injection.hh:115
wavearray< double > ** pwf
[x1,y1,z1,x2,y2,z2] components of spin vector
Definition: injection.hh:101
TBranch * b_ndim
pointer to the reconstructed waveform
Definition: injection.hh:105
Float_t * spin
[m1,m2], binary mass parameters
Definition: injection.hh:99
TTree * fChain
root input file cointainig the mdc TTree
Definition: injection.hh:59
Int_t nevent
Definition: injection.hh:66
Float_t * Deff
injected snr in the detectors
Definition: injection.hh:97
TBranch * b_time
Definition: injection.hh:128