Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
netevent.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 // Gabriele Vedovato, University of Padova
6 //////////////////////////////////////////////////////////
7 
8 
9 #ifndef netevent_h
10 #define netevent_h
11 
12 #include <TROOT.h>
13 #include <TChain.h>
14 #include <TFile.h>
15 #include "watfun.hh"
16 #include "watversion.hh"
17 #include "injection.hh"
18 #include "detector.hh"
19 #include "netcluster.hh"
20 #include "network.hh"
21 
22 //class detector;
23 //class netcluster;
24 //class network;
25 
26 /* Structure of WaveBurst network event */
27 
28 #define NETEVENT_INIT \
29  iFile(NULL),fChain(NULL),run(0),nevent(0),eventID(NULL),type(NULL),name(NULL),rate(NULL),\
30  volume(NULL),size(NULL),usize(0),gap(NULL),lag(NULL),slag(NULL),strain(NULL), \
31  phi(NULL),theta(NULL),psi(NULL),iota(NULL),bp(NULL),bx(NULL),time(NULL),gps(NULL), \
32  right(NULL),left(NULL),duration(NULL),start(NULL),stop(NULL),frequency(NULL), \
33  low(NULL),high(NULL),bandwidth(NULL),hrss(NULL),noise(NULL),erA(NULL),Psave(0), \
34  Psm(NULL),null(NULL),nill(NULL),netcc(NULL),neted(NULL),rho(NULL),gnet(0.),anet(0.), \
35  ecor(0.),norm(0.),ECOR(0.),penalty(0.),likelihood(0.),factor(0.),range(NULL), \
36  chirp(NULL),eBBH(NULL),Deff(NULL),mass(NULL),spin(NULL),snr(NULL),xSNR(NULL),sSNR(NULL), \
37  iSNR(NULL),oSNR(NULL),ioSNR(NULL),fP(NULL)
38 
39 class netevent {
40 
41 public :
42 
43  TFile *iFile; //!root input file cointainig the analyzed TTree
44 
45  TTree *fChain; //!pointer to the analyzed TTree or TChain
46  Int_t fCurrent; //!current Tree number in a TChain
47  Int_t ndim; //! number of detectors
48  Int_t Psave; //! max size used by allocate() for the probability maps
49 
50 // Declaration of leaves types
51 // for arrays: ifo1 - first index, ifo2 - second index, .....
52 
53  Int_t run; //! run ID
54  Int_t nevent; //! event count
55  Int_t* eventID; //! event ID
56  Int_t* type; //! event type: [0] - prod, [1]-sim
57  string* name; //! event name: "" - prod, mdc_name - sim
58  Int_t* rate; //! 1/rate - wavelet time resolution
59 
60  Int_t* volume; //! cluster volume
61  Int_t* size; //! cluster size (black pixels only)
62  Int_t usize; //! cluster union size
63 
64  Float_t* gap; //! time between consecutive events
65  Float_t* lag; //! time lag [sec]
66  Float_t* slag; //! time slag [sec]
67  Double_t* strain; //! sqrt(h+*h+ + hx*hx)
68  Float_t* phi; //! [0]-reconstructed, [1]-injected phi angle, [2]-RA
69  Float_t* theta; //! [0]-reconstructed, [1]-injected theta angle, [2]-DEC
70  Float_t* psi; //! [0]-reconstructed psi or phase of gc, [1]-injected psi angle
71  Float_t* iota; //! [0]-reconstructed iota angle, [1]-injected iota angle
72  Float_t* bp; //! beam pattern coefficients for hp
73  Float_t* bx; //! beam pattern coefficients for hx
74 
75  Double_t* time; //! average center_of_gravity time
76  Double_t* gps; //! segment start GPS time
77  Float_t* right; //! min cluster time relative to segment start
78  Float_t* left; //! max cluster time relative to segment start
79  Float_t* duration; //! cluster duration = stopW-startW
80  Double_t* start; //! GPS start time of the cluster
81  Double_t* stop; //! GPS stop time of the cluster
82 
83  Float_t* frequency; //! average center_of_snr frequency
84  Float_t* low; //! min frequency
85  Float_t* high; //! max frequency
86  Float_t* bandwidth; //! high-low
87  Double_t* hrss; //! hrss
88  Double_t* noise; //! noise rms
89  Float_t* erA; //! error angle
90  skymap* Psm; //! probability cc skymap
91  Float_t* null; //! un-biased null statistics
92  Float_t* nill; //! biased null statistics
93  Float_t* rho; //! effective correlated SNR
94  Float_t* netcc; //! network correlation coefficients: 0-net,1-pc,2-cc,3-net2
95  Float_t* neted; //! network energy disbalance: 0 - total, 1 - 00-phase, 2 - 90-phase
96  // 3 - L00-L90, 4 - total abs ED
97 
98  Float_t gnet; // network sensitivity
99  Float_t anet; // network alignment factor
100  Float_t inet; // network index
101  Float_t ecor; // correlated energy
102  Float_t norm; // norm Factor or ellipticity
103  Float_t ECOR; // effective correlated energy
104  Float_t penalty; // penalty factor
105  Float_t likelihood; // network likelihood
106 
107  Float_t factor; // Multiplicative amplitude factor - simulation only
108  Float_t* range; //! range to source: [0/1]-rec/inj
109  Float_t* chirp; //! chirp array: 0-injmass,1-recmass,2-merr,3-tmrgr,4-terr,5-chi2
110  Float_t* eBBH; //! eBBH array
111  Float_t* Deff; //! effective range for each detector
112  Float_t* mass; //! mass[2], binary mass parameters
113  Float_t* spin; //! spin[6], binary spin parameters
114 
115  Float_t* snr; //! energy/noise_variance
116  Float_t* xSNR; //! data-signal correlation Xk*Sk
117  Float_t* sSNR; //! energy of reconstructed responses Sk*Sk
118  Float_t* iSNR; //! injected snr waveform
119  Float_t* oSNR; //! reconstructed snr waveform
120  Float_t* ioSNR; //! injected reconstructed xcor waveform
121 
122  FILE* fP; //! dump file
123 
124  std::vector<detector*> ifoList; // detectors
125 
126 //List of branches
127 
128  TBranch *b_ndim; //!
129  TBranch *b_run; //!
130  TBranch *b_nevent; //!
131  TBranch *b_eventID; //!
132  TBranch *b_type; //!
133  TBranch *b_name; //!
134  TBranch *b_rate; //!
135 
136  TBranch *b_volume; //!
137  TBranch *b_size; //!
138  TBranch *b_usize; //!
139 
140  TBranch *b_gap; //!
141  TBranch *b_lag; //!
142  TBranch *b_slag; //!
143  TBranch *b_strain; //!
144  TBranch *b_phi; //!
145  TBranch *b_theta; //!
146  TBranch *b_psi; //!
147  TBranch *b_iota; //!
148  TBranch *b_bp; //!
149  TBranch *b_bx; //!
150 
151  TBranch *b_time; //!
152  TBranch *b_gps; //!
153  TBranch *b_right; //!
154  TBranch *b_left; //!
155  TBranch *b_duration; //!
156  TBranch *b_start; //!
157  TBranch *b_stop; //!
158 
159  TBranch *b_frequency; //!
160  TBranch *b_low; //!
161  TBranch *b_high; //!
162  TBranch *b_bandwidth; //!
163 
164  TBranch *b_hrss; //!
165  TBranch *b_noise; //!
166  TBranch *b_erA; //!
167  TBranch *b_Psm; //!
168  TBranch *b_null; //!
169  TBranch *b_nill; //!
170  TBranch *b_netcc; //!
171  TBranch *b_neted; //!
172  TBranch *b_rho; //!
173 
174  TBranch *b_gnet; //!
175  TBranch *b_anet; //!
176  TBranch *b_inet; //!
177  TBranch *b_ecor; //!
178  TBranch *b_norm; //!
179  TBranch *b_ECOR; //!
180  TBranch *b_penalty; //!
181  TBranch *b_likelihood;//!
182 
183  TBranch *b_factor; //!
184  TBranch *b_range; //!
185  TBranch *b_chirp; //!
186  TBranch *b_eBBH; //!
187  TBranch *b_Deff; //!
188  TBranch *b_mass; //!
189  TBranch *b_spin; //!
190 
191  TBranch *b_snr; //!
192  TBranch *b_xSNR; //!
193  TBranch *b_sSNR; //!
194  TBranch *b_iSNR; //!
195  TBranch *b_oSNR; //!
196  TBranch *b_ioSNR; //!
197 
199  { ndim=1; Psave=0; allocate(); init(); return; }
200 
201  netevent(int n, int Psave=0) : NETEVENT_INIT
202  { ndim=n; this->Psave=Psave; allocate(); init(); return; }
203 
204  netevent(const netevent& a) : NETEVENT_INIT
205  { ndim=a.ndim; Psave=a.Psave; allocate(); init(); *this = a; return; }
206 
207  netevent(TTree *tree, int n) : NETEVENT_INIT
208  { ndim=n; allocate(); init(); if(tree) Init(tree); return; }
209 
211  { TTree* tree=Init(fName,n); allocate(); init(); if(tree) Init(tree); return; }
212 
213  virtual ~netevent() {
214 // if (fChain) free(fChain->GetCurrentFile());
215  if (eventID) free(eventID); // event ID: 1/2 - prod/sim
216  if (type) free(type); // event type: 1/2 - prod/sim
217  if (name) delete name; // event name: "" - prod, mdc_name - sim
218  if (rate) free(rate); // 1/rate - wavelet time resolution
219 
220  if (volume) free(volume); // cluster volume
221  if (size) free(size); // cluster size (black pixels only)
222 
223  if (gap) free(gap); // time between consecutive events
224  if (lag) free(lag); // time lag [sec]
225  if (slag) free(slag); // time slag [sec]
226  if (strain) free(strain); // GW strain: 1/2 - prod/sim
227  if (phi) free(phi); // phi: 1/2 - prod/sim
228  if (theta) free(theta); // theta: 1/2 - prod/sim
229  if (psi) free(psi); // psi: 1/2 - prod/sim
230  if (iota) free(iota); // iota: 0/1 - prod/sim
231 
232  if (bp) free(bp); // beam pattern coefficients for hp
233  if (bx) free(bx); // beam pattern coefficients for hx
234 
235  if (time) free(time); // average center_of_snr time for prod and sim
236  if (right) free(right); // min cluster time
237  if (left) free(left); // max cluster time
238  if (duration) free(duration); // cluster duration = stopW-startW
239  if (start) free(start); // actual start GPS time
240  if (stop) free(stop); // actual stop GPS time
241 
242  if (frequency) free(frequency); // average center_of_snr frequency
243  if (low) free(low); // min frequency
244  if (high) free(high); // max frequency
245  if (bandwidth) free(bandwidth); // high-low
246  if (hrss) free(hrss); // log10(calibrated hrss)
247  if (noise) free(noise); // log10(calibrated noise rms)
248  if (erA) free(erA); // error angle
249  if (null) free(null); // un-biased null statistics
250  if (nill) free(nill); // biased null statistics
251  if (rho) free(rho); // effective correlated SNR
252  if (netcc) free(netcc); // correlation coefficients: 0-net,1-pc,2-cc,3-net2
253  if (neted) free(neted); // network energy disbalabce
254  if (range) free(range); // range array: 0-reconstructed, 1-injected
255  if (chirp) free(chirp); // chirp array: 0-injmass,1-recmass,2-merr,3-tmrgr,4-terr,5-chi2
256  if (eBBH) free(eBBH); // eBBH array: 0-rec ecc, 1-inj ecc
257  if (Deff) free(Deff); // sim effective distance for each detector
258  if (mass) free(mass); // mass[2], binary mass parameters
259  if (spin) free(spin); // spin[6], binary spin parameters
260  if (snr) free(snr); // energy/noise_variance
261  if (xSNR) free(xSNR); // x-s snr of the detector responses
262  if (sSNR) free(sSNR); // signal snr of the detector responses
263  if (iSNR) free(iSNR); // injected snr
264  if (oSNR) free(oSNR); // recontructed snr
265  if (ioSNR) free(ioSNR); // injected recontructed xcor
266 
267  if (Psm) delete Psm; // probability skymap
268 
269  if(iFile) delete iFile;
270  };
271 
272  virtual netevent& operator=(const netevent &);
273 
274  Int_t GetEntries();
275  Int_t GetEntry(Int_t);
276  void allocate();
277  void init();
278  TTree* Init(TString fName, int n);
279  void Init(TTree *);
280  Bool_t Notify();
281  TTree* setTree();
282  void setSLags(float* slag);
283  void output(TTree* = NULL, network* = NULL, double = 0., size_t = 0, int = -1);
284  void output2G(TTree*, network* , size_t, int, double);
285 
286 // void Loop();
287 // Int_t Cut(Int_t entry);
288 // Int_t LoadTree(Int_t entry);
289  void Show(Int_t entry = -1);
290 
291  inline void dopen(const char *fname, char* mode, bool header=true) {
292  if(fP != NULL) fclose(fP);
293  if((fP = fopen(fname, mode)) == NULL) {
294  cout << "netevent::Dump() error: cannot open file " << fname <<". \n";
295  return;
296  };
297  if(header) fprintf(fP,"# WAT Version : %s - GIT Revision : %s - Tag/Branch : %s",watversion('f'),watversion('r'),watversion('b'));
298  }
299 
300  inline void dclose() {
301  if(fP!=NULL) fclose(fP);
302  fP = NULL;
303  return;
304  }
305 
306  inline void Dump(TString analysis="2G") {
307  if(fP==NULL || ndim<1) return;
308  size_t i;
309  size_t I = ndim;
310 
311  fprintf(fP,"nevent: %d\n",nevent);
312  fprintf(fP,"ndim: %d\n",ndim);
313  fprintf(fP,"run: %d\n",run);
314  fprintf(fP,"name: %s\n",name->c_str());
315  fprintf(fP,"rho: %f\n",rho[analysis=="2G"?0:1]);
316  fprintf(fP,"netCC: %f\n",netcc[analysis=="2G"?0:0]);
317  fprintf(fP,"netED: %f\n",neted[0]/ecor);
318  fprintf(fP,"penalty: %f\n",penalty);
319  fprintf(fP,"gnet: %f\n",gnet);
320  fprintf(fP,"anet: %f\n",anet);
321  fprintf(fP,"inet: %f\n",inet);
322  fprintf(fP,"likelihood: %e\n",likelihood);
323  fprintf(fP,"ecor: %e\n",ecor);
324  fprintf(fP,"ECOR: %e\n",ECOR);
325  fprintf(fP,"factor: %f\n",factor);
326  fprintf(fP,"range: %f\n",range[0]);
327  fprintf(fP,"mchirp: %f\n",chirp[0]);
328  fprintf(fP,"norm: %f\n",norm);
329  fprintf(fP,"usize: %d\n",usize);
330 
331  fprintf(fP,"ifo: "); for(i=0; i<I; i++) fprintf(fP,"%s ",ifoList[i]->Name);fprintf(fP,"\n");
332  fprintf(fP,"eventID: "); for(i=0; i<2; i++) fprintf(fP,"%d ",eventID[i]); fprintf(fP,"\n");
333  fprintf(fP,"rho: "); for(i=0; i<2; i++) fprintf(fP,"%f ",rho[i]); fprintf(fP,"\n");
334  fprintf(fP,"type: "); for(i=0; i<2; i++) fprintf(fP,"%d ",type[i]); fprintf(fP,"\n");
335  fprintf(fP,"rate: "); for(i=0; i<I; i++) fprintf(fP,"%d ",rate[i]); fprintf(fP,"\n");
336  fprintf(fP,"volume: "); for(i=0; i<I; i++) fprintf(fP,"%d ",volume[i]); fprintf(fP,"\n");
337  fprintf(fP,"size: "); for(i=0; i<I; i++) fprintf(fP,"%d ",size[i]); fprintf(fP,"\n");
338  fprintf(fP,"lag: "); for(i=0; i<I; i++) fprintf(fP,"%f ",lag[i]); fprintf(fP,"\n");
339  fprintf(fP,"slag: "); for(i=0; i<I; i++) fprintf(fP,"%f ",slag[i]); fprintf(fP,"\n");
340  fprintf(fP,"phi: "); for(i=0; i<4; i++) fprintf(fP,"%f ",phi[i]); fprintf(fP,"\n");
341  fprintf(fP,"theta: "); for(i=0; i<4; i++) fprintf(fP,"%f ",theta[i]); fprintf(fP,"\n");
342  fprintf(fP,"psi: "); for(i=0; i<2; i++) fprintf(fP,"%f ",psi[i]); fprintf(fP,"\n");
343  fprintf(fP,"iota: "); for(i=0; i<2; i++) fprintf(fP,"%f ",iota[i]); fprintf(fP,"\n");
344  fprintf(fP,"bp: "); for(i=0; i<I; i++) fprintf(fP,"%7.4f ",bp[i]); fprintf(fP,"\n");
345  fprintf(fP,"inj_bp: "); for(i=I; i<2*I; i++) fprintf(fP,"%7.4f ",bp[i]); fprintf(fP,"\n");
346  fprintf(fP,"bx: "); for(i=0; i<I; i++) fprintf(fP,"%7.4f ",bx[i]); fprintf(fP,"\n");
347  fprintf(fP,"inj_bx: "); for(i=I; i<2*I; i++) fprintf(fP,"%7.4f ",bx[i]); fprintf(fP,"\n");
348  fprintf(fP,"chirp: "); for(i=0; i<6; i++) fprintf(fP,"%f ",chirp[i]); fprintf(fP,"\n");
349  fprintf(fP,"range: "); for(i=0; i<2; i++) fprintf(fP,"%f ",range[i]); fprintf(fP,"\n");
350  fprintf(fP,"Deff: "); for(i=0; i<I; i++) fprintf(fP,"%f ",Deff[i]); fprintf(fP,"\n");
351  fprintf(fP,"mass: "); for(i=0; i<2; i++) fprintf(fP,"%f ",mass[i]); fprintf(fP,"\n");
352  fprintf(fP,"spin: "); for(i=0; i<6; i++) fprintf(fP,"%f ",spin[i]); fprintf(fP,"\n");
353  fprintf(fP,"eBBH: "); for(i=0; i<4; i++) fprintf(fP,"%f ",eBBH[i]); fprintf(fP,"\n");
354  fprintf(fP,"null: "); for(i=0; i<I; i++) fprintf(fP,"%e ",null[i]); fprintf(fP,"\n");
355  fprintf(fP,"strain: "); for(i=0; i<2; i++) fprintf(fP,"%e ",strain[i]); fprintf(fP,"\n");
356  fprintf(fP,"hrss: "); for(i=0; i<I; i++) fprintf(fP,"%e ",hrss[i]); fprintf(fP,"\n");
357  fprintf(fP,"inj_hrss: "); for(i=I; i<2*I; i++) fprintf(fP,"%e ",hrss[i]); fprintf(fP,"\n");
358  fprintf(fP,"noise: "); for(i=0; i<I; i++) fprintf(fP,"%e ",noise[i]); fprintf(fP,"\n");
359 
360  // use seglen of detector at lag=0 to set segment value of the other detectors
361  int mdet=0; for(i=0; i<I; i++) if(lag[i]==0) mdet=i;
362  fprintf(fP,"segment: "); for(i=0; i<I; i++) {
363  double seglen = left[mdet]+right[mdet]+duration[1];
364  fprintf(fP,"%12.4f %12.4f ",gps[i],gps[i]+seglen);
365  }
366  fprintf(fP,"\n");
367 
368  fprintf(fP,"start: "); for(i=0; i<I; i++) fprintf(fP,"%12.4f ",start[i]); fprintf(fP,"\n");
369  fprintf(fP,"time: "); for(i=0; i<I; i++) fprintf(fP,"%12.4f ",time[i]); fprintf(fP,"\n");
370  fprintf(fP,"stop: "); for(i=0; i<I; i++) fprintf(fP,"%12.4f ",stop[i]); fprintf(fP,"\n");
371  fprintf(fP,"inj_time: "); for(i=I; i<2*I; i++) fprintf(fP,"%12.4f ",time[i]); fprintf(fP,"\n");
372  fprintf(fP,"left: "); for(i=0; i<I; i++) fprintf(fP,"%f ",left[i]); fprintf(fP,"\n");
373  fprintf(fP,"right: "); for(i=0; i<I; i++) fprintf(fP,"%f ",right[i]); fprintf(fP,"\n");
374  fprintf(fP,"duration: "); for(i=0; i<2; i++) fprintf(fP,"%f ",duration[i]); fprintf(fP,"\n");
375  fprintf(fP,"frequency: "); for(i=0; i<2; i++) fprintf(fP,"%f ",frequency[i]); fprintf(fP,"\n");
376  fprintf(fP,"low: "); for(i=0; i<1; i++) fprintf(fP,"%f ",low[i]); fprintf(fP,"\n");
377  fprintf(fP,"high: "); for(i=0; i<1; i++) fprintf(fP,"%f ",high[i]); fprintf(fP,"\n");
378  fprintf(fP,"bandwidth: "); for(i=0; i<2; i++) fprintf(fP,"%f ",bandwidth[i]); fprintf(fP,"\n");
379 
380  fprintf(fP,"snr: "); for(i=0; i<I; i++) fprintf(fP,"%e ",snr[i]); fprintf(fP,"\n");
381  fprintf(fP,"xSNR: "); for(i=0; i<I; i++) fprintf(fP,"%e ",xSNR[i]); fprintf(fP,"\n");
382  fprintf(fP,"sSNR: "); for(i=0; i<I; i++) fprintf(fP,"%e ",sSNR[i]); fprintf(fP,"\n");
383  fprintf(fP,"iSNR: "); for(i=0; i<I; i++) fprintf(fP,"%f ",iSNR[i]); fprintf(fP,"\n");
384  fprintf(fP,"oSNR: "); for(i=0; i<I; i++) fprintf(fP,"%f ",oSNR[i]); fprintf(fP,"\n");
385  fprintf(fP,"ioSNR: "); for(i=0; i<I; i++) fprintf(fP,"%f ",ioSNR[i]); fprintf(fP,"\n");
386 
387  fprintf(fP,"netcc: "); for(i=0; i<4; i++) fprintf(fP,"%f ",netcc[i]); fprintf(fP,"\n");
388  fprintf(fP,"neted: "); for(i=0; i<5; i++) fprintf(fP,"%f ",neted[i]); fprintf(fP,"\n");
389  fprintf(fP,"erA: "); for(i=0; i<11; i++) fprintf(fP,"%6.3f ",erA[i]); fprintf(fP,"\n");
390  };
391 
392  // used by THtml doc
393  ClassDef(netevent,1)
394 };
395 #endif
TTree * tree
Definition: TimeSortTree.C:20
Float_t anet
Definition: netevent.hh:99
TBranch * b_Psm
Definition: netevent.hh:167
TBranch * b_high
Definition: netevent.hh:161
Float_t * mass
effective range for each detector
Definition: netevent.hh:112
TBranch * b_ECOR
Definition: netevent.hh:179
Float_t * eBBH
chirp array: 0-injmass,1-recmass,2-merr,3-tmrgr,4-terr,5-chi2
Definition: netevent.hh:110
TBranch * b_hrss
Definition: netevent.hh:164
Int_t ndim
current Tree number in a TChain
Definition: netevent.hh:47
Float_t * rho
biased null statistics
Definition: netevent.hh:93
Float_t * high
min frequency
Definition: netevent.hh:85
Float_t * range
Definition: netevent.hh:108
TBranch * b_volume
Definition: netevent.hh:136
Float_t ecor
Definition: netevent.hh:101
void setSLags(float *slag)
Definition: netevent.cc:404
TBranch * b_range
Definition: netevent.hh:184
Bool_t Notify()
Definition: netevent.cc:304
Double_t * start
cluster duration = stopW-startW
Definition: netevent.hh:80
fprintf(stdout,"start=%f duration=%f rate=%f\n", x.start(), x.size()/x.rate(), x.rate())
Float_t * duration
max cluster time relative to segment start
Definition: netevent.hh:79
Float_t * low
average center_of_snr frequency
Definition: netevent.hh:84
TBranch * b_theta
Definition: netevent.hh:145
wavearray< double > a(hp.size())
TBranch * b_type
Definition: netevent.hh:132
int n
Definition: cwb_net.C:10
TTree * setTree()
Definition: netevent.cc:412
char * watversion(char c='s')
Definition: watversion.hh:21
TString("c")
TBranch * b_ndim
Definition: netevent.hh:128
Float_t inet
Definition: netevent.hh:100
Int_t * size
cluster volume
Definition: netevent.hh:61
Float_t * right
segment start GPS time
Definition: netevent.hh:77
TBranch * b_neted
Definition: netevent.hh:171
Int_t GetEntry(Int_t)
Definition: netevent.cc:387
Float_t ECOR
Definition: netevent.hh:103
TBranch * b_duration
Definition: netevent.hh:155
Float_t * gap
cluster union size
Definition: netevent.hh:64
TBranch * b_rho
Definition: netevent.hh:172
Float_t * left
min cluster time relative to segment start
Definition: netevent.hh:78
Int_t fCurrent
pointer to the analyzed TTree or TChain
Definition: netevent.hh:46
TBranch * b_norm
Definition: netevent.hh:178
Float_t * ioSNR
reconstructed snr waveform
Definition: netevent.hh:120
TTree * fChain
root input file cointainig the analyzed TTree
Definition: netevent.hh:45
TBranch * b_ioSNR
Definition: netevent.hh:196
TTree * Init(TString fName, int n)
Definition: netevent.cc:29
TBranch * b_snr
Definition: netevent.hh:191
Int_t * type
event ID
Definition: netevent.hh:56
i drho i
TBranch * b_ecor
Definition: netevent.hh:177
TBranch * b_inet
Definition: netevent.hh:176
TBranch * b_iota
Definition: netevent.hh:147
void dopen(const char *fname, char *mode, bool header=true)
Definition: netevent.hh:291
Int_t run
max size used by allocate() for the probability maps
Definition: netevent.hh:53
size_t mode
TBranch * b_psi
Definition: netevent.hh:146
TBranch * b_phi
Definition: netevent.hh:144
Float_t * frequency
GPS stop time of the cluster.
Definition: netevent.hh:83
TBranch * b_rate
Definition: netevent.hh:134
TBranch * b_netcc
Definition: netevent.hh:170
Float_t factor
Definition: netevent.hh:107
void output2G(TTree *, network *, size_t, int, double)
Definition: netevent.cc:687
Int_t nevent
run ID
Definition: netevent.hh:54
TBranch * b_low
Definition: netevent.hh:160
Int_t Psave
number of detectors
Definition: netevent.hh:48
TBranch * b_left
Definition: netevent.hh:154
Int_t * volume
1/rate - wavelet time resolution
Definition: netevent.hh:60
TBranch * b_spin
Definition: netevent.hh:189
TBranch * b_penalty
Definition: netevent.hh:180
Float_t * Deff
eBBH array
Definition: netevent.hh:111
Double_t * gps
average center_of_gravity time
Definition: netevent.hh:76
#define NETEVENT_INIT
Definition: netevent.hh:28
TBranch * b_eBBH
Definition: netevent.hh:186
TBranch * b_eventID
Definition: netevent.hh:131
TBranch * b_nevent
Definition: netevent.hh:130
TBranch * b_erA
Definition: netevent.hh:166
TBranch * b_sSNR
Definition: netevent.hh:193
Float_t penalty
Definition: netevent.hh:104
Double_t * stop
GPS start time of the cluster.
Definition: netevent.hh:81
TBranch * b_null
Definition: netevent.hh:168
Float_t gnet
network energy disbalance: 0 - total, 1 - 00-phase, 2 - 90-phase
Definition: netevent.hh:98
void Dump(TString analysis="2G")
Definition: netevent.hh:306
TBranch * b_time
Definition: netevent.hh:151
Float_t * lag
time between consecutive events
Definition: netevent.hh:65
char fname[1024]
Float_t * xSNR
energy/noise_variance
Definition: netevent.hh:116
TBranch * b_lag
Definition: netevent.hh:141
TBranch * b_name
Definition: netevent.hh:133
TFile * iFile
Definition: netevent.hh:43
std::vector< detector * > ifoList
dump file
Definition: netevent.hh:124
void Show(Int_t entry=-1)
Definition: netevent.cc:393
TBranch * b_likelihood
Definition: netevent.hh:181
Double_t * strain
time slag [sec]
Definition: netevent.hh:67
Float_t * netcc
effective correlated SNR
Definition: netevent.hh:94
Float_t * slag
time lag [sec]
Definition: netevent.hh:66
TBranch * b_gnet
Definition: netevent.hh:174
TBranch * b_anet
Definition: netevent.hh:175
void output(TTree *=NULL, network *=NULL, double=0., size_t=0, int=-1)
Definition: netevent.cc:1166
TBranch * b_bandwidth
Definition: netevent.hh:162
TBranch * b_usize
Definition: netevent.hh:138
Definition: skymap.hh:45
double * entry
Definition: cwb_setcuts.C:206
TBranch * b_bp
Definition: netevent.hh:148
TBranch * b_slag
Definition: netevent.hh:142
Float_t * iota
[0]-reconstructed psi or phase of gc, [1]-injected psi angle
Definition: netevent.hh:71
TBranch * b_strain
Definition: netevent.hh:143
TBranch * b_frequency
Definition: netevent.hh:159
Float_t * psi
[0]-reconstructed, [1]-injected theta angle, [2]-DEC
Definition: netevent.hh:70
Double_t * hrss
high-low
Definition: netevent.hh:87
TBranch * b_gap
Definition: netevent.hh:140
TBranch * b_start
Definition: netevent.hh:156
virtual ~netevent()
Definition: netevent.hh:213
TBranch * b_size
Definition: netevent.hh:137
TBranch * b_gps
Definition: netevent.hh:152
TBranch * b_oSNR
Definition: netevent.hh:195
FILE * fP
injected reconstructed xcor waveform
Definition: netevent.hh:122
Float_t * phi
sqrt(h+*h+ + hx*hx)
Definition: netevent.hh:68
Float_t norm
Definition: netevent.hh:102
TBranch * b_Deff
Definition: netevent.hh:187
TBranch * b_noise
Definition: netevent.hh:165
Float_t * theta
[0]-reconstructed, [1]-injected phi angle, [2]-RA
Definition: netevent.hh:69
Double_t * time
beam pattern coefficients for hx
Definition: netevent.hh:75
TBranch * b_right
Definition: netevent.hh:153
Int_t * rate
event name: "" - prod, mdc_name - sim
Definition: netevent.hh:58
void dclose()
Definition: netevent.hh:300
TBranch * b_iSNR
Definition: netevent.hh:194
Float_t * nill
un-biased null statistics
Definition: netevent.hh:92
Float_t * spin
mass[2], binary mass parameters
Definition: netevent.hh:113
TBranch * b_mass
Definition: netevent.hh:188
Float_t * sSNR
data-signal correlation Xk*Sk
Definition: netevent.hh:117
Float_t * erA
noise rms
Definition: netevent.hh:89
TBranch * b_stop
Definition: netevent.hh:157
Float_t * null
probability cc skymap
Definition: netevent.hh:91
Float_t likelihood
Definition: netevent.hh:105
Float_t * bx
beam pattern coefficients for hp
Definition: netevent.hh:73
TBranch * b_bx
Definition: netevent.hh:149
virtual netevent & operator=(const netevent &)
Definition: netevent.cc:587
Double_t * noise
hrss
Definition: netevent.hh:88
string * name
event type: [0] - prod, [1]-sim
Definition: netevent.hh:57
Float_t * snr
spin[6], binary spin parameters
Definition: netevent.hh:115
TBranch * b_factor
Definition: netevent.hh:183
Int_t GetEntries()
Definition: netevent.cc:381
Float_t * neted
network correlation coefficients: 0-net,1-pc,2-cc,3-net2
Definition: netevent.hh:95
fclose(ftrig)
TBranch * b_xSNR
Definition: netevent.hh:192
char fName[256]
Float_t * chirp
range to source: [0/1]-rec/inj
Definition: netevent.hh:109
TBranch * b_nill
Definition: netevent.hh:169
TBranch * b_chirp
Definition: netevent.hh:185
TBranch * b_run
Definition: netevent.hh:129
skymap * Psm
error angle
Definition: netevent.hh:90
netevent(int n, int Psave=0) this Psave
Definition: netevent.hh:202
Float_t * bp
[0]-reconstructed iota angle, [1]-injected iota angle
Definition: netevent.hh:72
Float_t * iSNR
energy of reconstructed responses Sk*Sk
Definition: netevent.hh:118
Int_t * eventID
event count
Definition: netevent.hh:55
Float_t * bandwidth
max frequency
Definition: netevent.hh:86
Int_t usize
cluster size (black pixels only)
Definition: netevent.hh:62
netevent(TTree *tree, int n allocate)()
Definition: netevent.hh:208
Float_t * oSNR
injected snr waveform
Definition: netevent.hh:119