Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
CWB_Plugin_QLWveto.C
Go to the documentation of this file.
1 #define XIFO 4
2 
3 #pragma GCC system_header
4 
5 #include "cwb.hh"
6 #include "config.hh"
7 #include "network.hh"
8 #include "wavearray.hh"
9 #include "TString.h"
10 #include "TObjArray.h"
11 #include "TObjString.h"
12 #include "TRandom.h"
13 #include "TComplex.h"
14 #include "TMath.h"
15 #include "mdc.hh"
16 #include "watplot.hh"
17 #include "gwavearray.hh"
18 #include <vector>
19 
20 //#define PLOT_LIKELIHOOD
21 //#define PLOT_WHITENED_WAVEFORMS
22 
23 #define NTHR 1
24 #define ATHR 7.58859
25 
27 void GetLveto(netcluster* pwc, int cid, int nifo, float* Lveto);
28 void GetWveto(netcluster* pwc, int cid, int nifo, float* Wveto);
30  CWB::config* cfg, bool fft=false, bool strain=false);
31 
32 
33 void
35 //!MISCELLANEA
36 // Extract whitened reconstructed waveforms, and compute the Qveto, Lveto & Wveto parameters
37 
38  cout << endl;
39  cout << "-----> CWB_Plugin_QLWveto.C" << endl;
40  cout << "ifo " << ifo.Data() << endl;
41  cout << "type " << type << endl;
42  cout << endl;
43 
44  float Qveto[2*NIFO_MAX]; // Qveto
45  float Lveto[3]; // Lveto
46  float Wveto[2]; // Wveto
47 
48  if(type==CWB_PLUGIN_CONFIG) {
49  cfg->outPlugin=true; // disable built-in output root file
50  }
51 
52  if(type==CWB_PLUGIN_ILIKELIHOOD) {
53  NET->wfsave=true; // enable waveform save
54 
55  // search output root file in the system list
56  TFile* froot = NULL;
57  TList *files = (TList*)gROOT->GetListOfFiles();
58  TString outDump="";
59  netevent* EVT;
60  int nIFO = NET->ifoListSize(); // number of detectors
61  if (files) {
62  TIter next(files);
63  TSystemFile *file;
64  TString fname;
65  bool check=false;
66  while ((file=(TSystemFile*)next())) {
67  fname = file->GetName();
68  // set output root file as the current file
69  if(fname.Contains("wave_")) {
70  froot=(TFile*)file;froot->cd();
71  outDump=fname;
72  outDump.ReplaceAll(".root.tmp",".txt");
73  //cout << "output file name : " << fname << endl;
74  }
75  }
76  if(!froot) {
77  cout << "CWB_Plugin_QLWveto.C : Error - output root file not found" << endl;
78  gSystem->Exit(1);
79  }
80  } else {
81  cout << "CWB_Plugin_QLWveto.C : Error - output root file not found" << endl;
82  gSystem->Exit(1);
83  }
84 
85  TTree* net_tree = (TTree *) froot->Get("waveburst");
86  if(net_tree==NULL) {
87  EVT = new netevent(nIFO);
88  net_tree = EVT->setTree();
89  net_tree->Branch("Qveto",Qveto,TString::Format("Qveto[%i]/F",2*cfg->nIFO));
90  net_tree->Branch("Lveto",Lveto,TString::Format("Lveto[%i]/F",3));
91  net_tree->Branch("Wveto",Wveto,TString::Format("Wveto[%i]/F",2));
92  }
93  }
94 
95  if(type==CWB_PLUGIN_OLIKELIHOOD) {
96 
97  if(TString(cfg->analysis)!="2G") {
98  cout << "CWB_Plugin_QLWveto.C -> "
99  << "CWB_PLUGIN_OLIKELIHOOD implemented only for 2G" << endl;
100  gSystem->Exit(1);
101  }
102 
103  // import ifactor
104  int gIFACTOR=-1; IMPORT(int,gIFACTOR)
105  cout << "-----> CWB_Plugin_QLWveto.C -> "
106  << " gIFACTOR : " << gIFACTOR << endl;
107 
108  // import slagShift
109  float* gSLAGSHIFT=NULL; IMPORT(float*,gSLAGSHIFT)
110 
111  int nIFO = NET->ifoListSize(); // number of detectors
112  int K = NET->nLag; // number of time lag
113  netevent* EVT;
115  //double factor = cfg->simulation==3||cfg->simulation==4 ? -gIFACTOR : cfg->factors[gIFACTOR];
116  double factor = cfg->factors[gIFACTOR];
117  int rate = 0; // select all resolutions
118 
119  // search output root file in the system list
120  TFile* froot = NULL;
121  TList *files = (TList*)gROOT->GetListOfFiles();
122  TString outDump="";
123  if (files) {
124  TIter next(files);
125  TSystemFile *file;
126  TString fname;
127  bool check=false;
128  while ((file=(TSystemFile*)next())) {
129  fname = file->GetName();
130  // set output root file as the current file
131  if(fname.Contains("wave_")) {
132  froot=(TFile*)file;froot->cd();
133  outDump=fname;
134  outDump.ReplaceAll(".root.tmp",".txt");
135  //cout << "output file name : " << fname << endl;
136  }
137  }
138  if(!froot) {
139  cout << "CWB_Plugin_QLWveto.C : Error - output root file not found" << endl;
140  gSystem->Exit(1);
141  }
142  } else {
143  cout << "CWB_Plugin_QLWveto.C : Error - output root file not found" << endl;
144  gSystem->Exit(1);
145  }
146 
147  TTree* net_tree = (TTree *) froot->Get("waveburst");
148  if(net_tree!=NULL) {
149  EVT = new netevent(net_tree,nIFO);
150  net_tree->SetBranchAddress("Qveto",Qveto);
151  net_tree->SetBranchAddress("Lveto",Lveto);
152  net_tree->SetBranchAddress("Wveto",Wveto);
153  } else {
154  EVT = new netevent(nIFO);
155  net_tree = EVT->setTree();
156  net_tree->Branch("Qveto",Qveto,TString::Format("Qveto[%i]/F",2*cfg->nIFO));
157  net_tree->Branch("Lveto",Lveto,TString::Format("Lveto[%i]/F",3));
158  net_tree->Branch("Wveto",Wveto,TString::Format("Wveto[%i]/F",2));
159  }
160  EVT->setSLags(gSLAGSHIFT); // set slags into netevent
161 
162  for(int k=0; k<K; k++) { // loop over the lags
163 
164  id = NET->getwc(k)->get(const_cast<char*>("ID"), 0, 'L', rate);
165 
166  for(int j=0; j<(int)id.size(); j++) { // loop over cluster index
167 
168  int ID = size_t(id.data[j]+0.5);
169 
170  if(NET->getwc(k)->sCuts[ID-1]!=-1) continue; // skip rejected/processed clusters
171 
172  double ofactor=0;
173  if(cfg->simulation==4) ofactor=-factor;
174  else if(cfg->simulation==3) ofactor=-gIFACTOR;
175  else ofactor=factor;
176 
177  EVT->output2G(NULL,NET,ID,k,ofactor); // get reconstructed parameters
178 
179  wavearray<double>** pwfREC = new wavearray<double>*[nIFO];
180  detector* pd = NET->getifo(0);
181  int idSize = pd->RWFID.size();
182 
183  int wfIndex=-1;
184  for (int mm=0; mm<idSize; mm++) if (pd->RWFID[mm]==ID) wfIndex=mm;
185  if(wfIndex==-1) continue;
186 
187  netcluster* pwc = NET->getwc(k);
188  cout << endl << "----------------------------------------------------------------" << endl;
189  GetLveto(pwc, ID, nIFO, Lveto);
190  cout << "Lveto : " << "fmean : " << Lveto[0] << " frms : " << Lveto[1]
191  << " Energy Ratio : " << Lveto[2] << endl << endl;
192  GetWveto(pwc, ID, nIFO, Wveto);
193  cout << "Wveto : " << " Slope : " << Wveto[0] << " Correlation : " << Wveto[1] << endl << endl;
194 
195  // extract whitened reconstructed waveforms
196  for(int n=0; n<nIFO; n++) {
197 
198  pd = NET->getifo(n);
199 
200  pwfREC[n] = pd->RWFP[wfIndex];
201  wavearray<double>* wfREC = pwfREC[n]; // array of reconstructed waveforms
202 
203 #ifdef PLOT_WHITENED_WAVEFORMS
204  //PlotWaveform(NET->ifoName[n], wfREC, cfg, false, false);
205  PlotWaveform(NET->ifoName[n], wfREC, cfg, true, false);
206 #endif
207  // reconstructed whitened waveform
208  NET->getMRAwave(ID,k,'S',0,true);
209  Qveto[n] = GetQveto(&(pd->waveForm));
210  // whitened waveform
211  NET->getMRAwave(ID,k,'W',0,true);
212  Qveto[n+nIFO] = GetQveto(&(pd->waveBand));
213 
214  //Qveto[n] = GetQveto(wfREC);
215  cout << "Qveto : " << pd->Name << " Qveto[R] = " << Qveto[n]
216  << " Qveto[W] = " << Qveto[n+nIFO] << endl;
217  }
218  cout << "----------------------------------------------------------------" << endl;
219  delete [] pwfREC;
220 
221  std::vector<int> sCuts = NET->getwc(k)->sCuts; // save cCuts
222  // set sCuts=1 to the events which must be not copied with cps to pwc
223  for(int i=0; i<(int)sCuts.size(); i++) if(i!=ID-1) NET->getwc(k)->sCuts[i]=1;
224 
225  // ID can not be used to get the event, to get event use ID=1 (ex: for watplot)
226  NET->getwc(k)->sCuts = sCuts; // restore cCuts
227 
228  if(cfg->dump) EVT->dopen(outDump.Data(),const_cast<char*>("a"),false);
229  EVT->output2G(net_tree,NET,ID,k,ofactor); // get reconstructed parameters
230  if(cfg->dump) {
231  // add Qveto to dump file
232  fprintf(EVT->fP,"Qveto: ");
233  for(int i=0; i<2*nIFO; i++) fprintf(EVT->fP,"%f ",Qveto[i]);
234  fprintf(EVT->fP,"\n");
235  // add Lveto to dump file
236  fprintf(EVT->fP,"Lveto: ");
237  for(int i=0; i<3; i++) fprintf(EVT->fP,"%f ",Lveto[i]);
238  fprintf(EVT->fP,"\n");
239  // add Wveto to dump file
240  fprintf(EVT->fP,"Wveto: ");
241  for(int i=0; i<3; i++) fprintf(EVT->fP,"%f ",Wveto[i]);
242  fprintf(EVT->fP,"\n");
243  }
244  if(cfg->dump) EVT->dclose();
245  if(!cfg->cedDump) NET->getwc(k)->sCuts[ID-1]=1; // mark as processed
246  }
247  }
248 
249  jfile->cd();
250  if(EVT) delete EVT;
251  }
252  return;
253 }
254 
255 float
257 
258  wavearray<double> x = *wf;;
259 
260  // resample data by a factor 4
261  int xsize=x.size();
262  x.FFTW(1);
263  x.resize(4*x.size());
264  x.rate(4*x.rate());
265  for(int j=xsize;j<x.size();j++) x[j]=0;
266  x.FFTW(-1);
267 
268  // extract max/min values and save the absolute values in the array 'a'
269  wavearray<double> a(x.size());
270  int size=0;
271  double dt = 1./x.rate();
272  double prev=x[0];
273  double xmax=0;
274  for (int i=1;i<x.size();i++) {
275  if(fabs(x[i])>xmax) xmax=fabs(x[i]);
276  if(prev*x[i]<0) {
277  a[size]=xmax;
278  size++;
279  xmax=0;
280  }
281  prev=x[i];
282  }
283 
284  // find max value/index ans save on amax/imax
285  int imax=-1;
286  double amax=0;
287  for (int i=1;i<size;i++) {
288  if(a[i]>amax) {amax=a[i];imax=i;}
289  }
290 
291 /*
292  cout << endl;
293  cout << "a[imax-2] " << a[imax-2] << endl;
294  cout << "a[imax-1] " << a[imax-1] << endl;
295  cout << "a[imax] " << a[imax] << endl;
296  cout << "a[imax+1] " << a[imax+1] << endl;
297  cout << "a[imax+2] " << a[imax+2] << endl;
298  cout << endl;
299 */
300 
301  // compute Qveto
302  double ein=0; // energy of max values inside NTHR
303  double eout=0; // energy of max values outside NTHR
304  for (int i=0;i<size;i++) {
305  if(abs(imax-i)<=NTHR) {
306  ein+=a[i]*a[i];
307  //cout << i << " ein " << a[i] << " " << amax << endl;
308  } else {
309  if(a[i]>amax/ATHR) eout+=a[i]*a[i];
310  //if(a[i]>amax/ATHR) cout << i << " eout " << a[i] << " " << amax << endl;
311  }
312  }
313  float Qveto = ein>0 ? eout/ein : 0.;
314  //cout << "Qveto : " << Qveto << " ein : " << ein << " eout : " << eout << endl;
315 
316  return Qveto;
317 }
318 
319 void
320 GetLveto(netcluster* pwc, int cid, int nifo, float* Lveto) {
321 //
322 // input
323 // pwc : pointer to netcluster object
324 // cid : cluster id
325 // nifo : number of detectors
326 // output
327 // Lveto[0] : line frequency
328 // Lveto[1] : line bandwitdh
329 // Lveto[2] : line enery ratio (line_energy / total_energy)
330 //
331 
332  Lveto[0] = Lveto[1] = Lveto[2] = 0;
333 
334  std::vector<int>* vint = &(pwc->cList[cid-1]); // pixel list
335  int V = vint->size(); // cluster size
336  if(!V) return;
337 
338  // ------------------------------------------------------------------
339  // Find max pixel parameters
340  // ------------------------------------------------------------------
341 
342  double likeMax=0; // maximum pixel's energy
343  double likeTot=0; // total cluster energy
344  double freqMax; // frequency of the pixel with max energy
345  double dfreqMax; // df of the pixel with max energy
346  for(int n=0; n<V; n++) {
347  netpixel* pix = pwc->getPixel(cid,n);
348  if(pix->layers%2==0) {
349  cout << "CWB_Plugin_QLWveto.C - Error : is enabled only for WDM (2G)" << endl;
350  exit(1);
351  }
352  if(!pix->core) continue; // select only the principal components pixels
353 
354  double likePix=0;
355  for(int m=0; m<nifo; m++) {
356  likePix += pow(pix->getdata('S',m),2); // like whitened reconstructed signal 00
357  likePix += pow(pix->getdata('P',m),2); // like whitened reconstructed signal 90
358  }
359 
360  double freq = pix->frequency*pix->rate/2.;
361  double df = pix->rate/2.;
362 
363  likeTot+=likePix;
364  if(likePix>likeMax) {likeMax=likePix;freqMax=freq;dfreqMax=df;}
365  }
366  //cout << "likeMax : " << likeMax << " likeTot : " << likeTot
367  // << " freqMax : " << freqMax << " dfreqMax : " << dfreqMax << endl;
368 
369  // ------------------------------------------------------------------
370  // Compute Lveto parameters
371  // ------------------------------------------------------------------
372 
373  double fmean=0; // line mean frequency
374  double frms=0; // line bandwidth
375  double likeLin=0; // line energy
376  for(int n=0; n<V; n++) {
377  netpixel* pix = pwc->getPixel(cid,n);
378  if(!pix->core) continue; // select only the principal components pixels
379 
380  double likePix=0;
381  for(int m=0; m<nifo; m++) {
382  likePix += pow(pix->getdata('S',m),2); // like whitened reconstructed signal 00
383  likePix += pow(pix->getdata('P',m),2); // like whitened reconstructed signal 90
384  }
385 
386  // the estimation is done for all pixels
387  // with freq in the range [freqMax-dfreqMax, freqMax+dfreqMax]
388  double freq = pix->frequency*pix->rate/2.;
389  if(fabs(freq-freqMax)<=dfreqMax) {
390  likeLin += likePix;
391  fmean += likePix*freq;
392  frms += likePix*freq*freq;
393  }
394  }
395 
396  fmean = fmean/likeLin;
397  frms = frms/likeLin-fmean*fmean;
398  frms = frms>0 ? sqrt(frms) : 0.;
399 
400  if(frms<dfreqMax/2.) frms=dfreqMax/2.;
401 
402  // ------------------------------------------------------------------
403  // Save Lveto parameters
404  // ------------------------------------------------------------------
405 
406  Lveto[0] = fmean; // line mean frequency
407  Lveto[1] = frms; // line bandwidth
408  Lveto[2] = likeTot>0. ? likeLin/likeTot : 0.; // energy ratio energy inside_line/total
409 
410  // ------------------------------------------------------------------
411  // plot time-frequency energy
412  // ------------------------------------------------------------------
413 
414 #if defined PLOT_LIKELIHOOD
415  watplot WTS(const_cast<char*>("wts"));
416  WTS.plot(pwc, cid, nifo, 'L', 0, const_cast<char*>("COLZ"));
417  WTS.canvas->Print("l_tfmap_scalogram.png");
418 #endif
419 
420  return;
421 }
422 
423 void
424 GetWveto(netcluster* pwc, int cid, int nifo, float* Wveto) {
425 //
426 // input
427 // pwc : pointer to netcluster object
428 // cid : cluster id
429 // nifo : number of detectors
430 // output
431 // Wveto[0] : whistle slope
432 // Wveto[1] : whistle correlation
433 //
434 
435  Wveto[0] = Wveto[1] = Wveto[2] = Wveto[3] = Wveto[4] = 0;
436 
437  std::vector<int>* vint = &(pwc->cList[cid-1]); // pixel list
438  int V = vint->size(); // cluster size
439  if(!V) return;
440 
441  std::vector<double> x;
442  std::vector<double> y;
443  std::vector<double> w;
444 
445  // extract pixels
446  for(int j=0; j<V; j++) {
447  netpixel* pix = pwc->getPixel(cid,j);
448  if(pix->layers%2==0) {
449  cout << "CWB_Plugin_QLWveto.C - Error : is enabled only for WDM (2G)" << endl;
450  exit(1);
451  }
452  if(pix->likelihood<1. || pix->frequency==0) continue;
453  //if(!pix->core) continue; // select only the principal components pixels
454 
455  double time = int(double(pix->time)/pix->layers)/pix->rate; // time in seconds from the start
456  double freq = pix->frequency*pix->rate/2.;
457 
458  x.push_back(time);
459  y.push_back(freq);
460  w.push_back(pix->likelihood);
461  }
462  int size = x.size();
463  if(size<5) return;
464 
465  double xcm, ycm, qxx, qyy, qxy, ew;
466  xcm = ycm = qxx = qyy = qxy = ew = 0;
467 
468  for(int i=0; i<size; ++i) {
469  xcm += x[i]*w[i];
470  ycm += y[i]*w[i];
471  ew += w[i];
472  }
473  xcm /= ew;
474  ycm /= ew;
475 
476  for(int i=0; i<size; ++i) {
477  qxx += (x[i] - xcm)*(x[i] - xcm)*w[i];
478  qyy += (y[i] - ycm)*(y[i] - ycm)*w[i];
479  qxy += (x[i] - xcm)*(y[i] - ycm)*w[i];
480  }
481 
482  double beta = qxy/qxx; // slope
483  double alpha = ycm-beta*xcm; // intercept
484  double corr = qxy/sqrt(qxx*qyy); // correlation
485  double duration = sqrt(qxx/ew); // duration
486  double bandwidth = sqrt(qyy/ew); // bandwidth
487 
488  //printf("alpha : %lf , beta : %lf , corr : %lf, dur : %lf , bw : %lf \n",
489  // alpha, beta, corr, duration, bandwidth);
490 
491  Wveto[0] = beta; // slope
492  Wveto[1] = fabs(corr); // correlation
493 
494  return;
495 }
496 
497 void
499  CWB::config* cfg, bool fft, bool strain) {
500 
501  watplot PTS(const_cast<char*>("ptswrc"),200,20,800,500);
502 
503  //cout << "Print " << fname << endl;
504  double tmin = wfREC->start();
505  wfREC->start(wfREC->start()-tmin);
506  if(fft) {
507  PTS.plot(wfREC, const_cast<char*>("ALP"), 1, 0, 0, true, cfg->fLow, cfg->fHigh);
508  } else {
509  PTS.plot(wfREC, const_cast<char*>("ALP"), 1, 0, 0);
510  }
511  PTS.graph[0]->SetLineWidth(1);
512  wfREC->start(wfREC->start()+tmin);
513 
514  char label[64]="";
515  if(fft) sprintf(label,"%s","fft");
516  else sprintf(label,"%s","time");
517  if(strain) sprintf(label,"%s_%s",label,"strain");
518  else sprintf(label,"%s_%s",label,"white");
519 
520  char fname[1024];
521  //sprintf(fname, "%s_wf_%s_rec_gps_%d.root",ifo.Data(),label,int(tmin));
522  sprintf(fname, "%s_wf_%s_rec_gps_%d.png",ifo.Data(),label,int(tmin));
523  PTS.canvas->Print(fname);
524  cout << "write : " << fname << endl;
525  //PTS.canvas->Write(REPLACE(fname,dirCED,gtype));
526 }
std::vector< char * > ifoName
Definition: network.hh:591
detector * getifo(size_t n)
param: detector index
Definition: network.hh:418
CWB::config * cfg
Definition: TestCWB_Plugin.C:5
char analysis[8]
Definition: config.hh:99
TString outDump
virtual size_t size() const
Definition: wavearray.hh:127
size_t nLag
Definition: network.hh:555
double duration
void setSLags(float *slag)
Definition: netevent.cc:404
fprintf(stdout,"start=%f duration=%f rate=%f\n", x.start(), x.size()/x.rate(), x.rate())
virtual void rate(double r)
Definition: wavearray.hh:123
float GetQveto(wavearray< double > *wf)
wavearray< double > a(hp.size())
size_t frequency
Definition: netpixel.hh:93
float likelihood
Definition: netpixel.hh:96
int n
Definition: cwb_net.C:10
cout<< "skymap size : "<< L<< endl;for(int l=0;l< L;l++) sm.set(l, l);sm > const_cast< char * >("skymap.dat")
double beta
TTree * setTree()
Definition: netevent.cc:412
TString("c")
int ID
Definition: TestMDC.C:70
cout<< endl;cout<< "ts size = "<< ts.size()<< " ts rate = "<< ts.rate()<< endl;tf.Forward(ts, wdm);int levels=tf.getLevel();cout<< "tf size = "<< tf.size()<< endl;double dF=tf.resolution();double dT=1./(2 *dF);cout<< "rate(hz) : "<< RATE<< "\t layers : "<< nLAYERS<< "\t dF(hz) : "<< dF<< "\t dT(ms) : "<< dT *1000.<< endl;int itime=TIME_PIXEL_INDEX;int ifreq=FREQ_PIXEL_INDEX;int index=(levels+1)*itime+ifreq;double time=itime *dT;double freq=(ifreq >0)?ifreq *dF:dF/4;cout<< endl;cout<< "PIXEL TIME = "<< time<< " sec "<< endl;cout<< "PIXEL FREQ = "<< freq<< " Hz "<< endl;cout<< endl;wavearray< double > x
wavearray< double > get(char *name, size_t index=0, char atype='R', int type=1, bool=true)
param: string with parameter name param: index in the amplitude array, which define detector param: c...
Definition: netcluster.cc:2188
std::vector< wavearray< double > * > RWFP
Definition: detector.hh:364
double fLow
Definition: config.hh:122
netpixel pix(nifo)
netcluster * pwc
Definition: cwb_job_obj.C:20
std::vector< TGraph * > graph
Definition: watplot.hh:176
bool cedDump
Definition: config.hh:279
waveform wf
Long_t size
WSeries< double > waveBand
Definition: detector.hh:338
size_t layers
Definition: netpixel.hh:94
int m
Definition: cwb_net.C:10
std::vector< vector_int > cList
Definition: netcluster.hh:379
virtual void start(double s)
Definition: wavearray.hh:119
int j
Definition: cwb_net.C:10
bool getMRAwave(size_t ID, size_t lag, char atype='S', int mode=0, bool tof=false)
Definition: network.cc:3635
i drho i
void plot(wavearray< double > &, char *=NULL, int=1, double=0., double=0., bool=false, float=0., float=0., bool=false, float=0., bool=false)
Definition: watplot.cc:132
bool outPlugin
Definition: config.hh:351
bool core
Definition: netpixel.hh:102
char ifo[NIFO_MAX][8]
void dopen(const char *fname, char *mode, bool header=true)
Definition: netevent.hh:291
size_t ifoListSize()
Definition: network.hh:413
wavearray< double > w
Definition: Test1.C:27
#define nIFO
TCanvas * canvas
Definition: watplot.hh:174
double factor
#define ATHR
double time[6]
Definition: cbc_plots.C:435
jfile
Definition: cwb_job_obj.C:25
void output2G(TTree *, network *, size_t, int, double)
Definition: netevent.cc:687
TTree * net_tree
int simulation
Definition: config.hh:181
void GetLveto(netcluster *pwc, int cid, int nifo, float *Lveto)
#define IMPORT(TYPE, VAR)
Definition: cwb.hh:51
watplot * WTS
Definition: ChirpMass.C:115
i() int(T_cor *100))
network NET
Definition: cwb_dump_inj.C:12
TString label
Definition: MergeTrees.C:21
const int NIFO_MAX
Definition: wat.hh:4
int gIFACTOR
#define NTHR
TIter next(twave->GetListOfBranches())
char fname[1024]
void GetWveto(netcluster *pwc, int cid, int nifo, float *Wveto)
std::vector< int > RWFID
Definition: detector.hh:363
int k
size_t time
Definition: netpixel.hh:92
TFile * froot
netpixel * getPixel(size_t n, size_t i)
Definition: netcluster.hh:395
double dt
virtual void FFTW(int=1)
Definition: wavearray.cc:878
int nifo
FILE * fP
injected reconstructed xcor waveform
Definition: netevent.hh:122
bool wfsave
Definition: network.hh:582
netevent EVT(itree, nifo)
std::vector< int > sCuts
Definition: netcluster.hh:374
WSeries< double > waveForm
Definition: detector.hh:337
void dclose()
Definition: netevent.hh:300
double fHigh
Definition: config.hh:123
void CWB_Plugin(TFile *jfile, CWB::config *cfg, network *NET, WSeries< double > *x, TString ifo, int type)
COHERENCE.
char Name[16]
Definition: detector.hh:309
double fabs(const Complex &x)
Definition: numpy.cc:37
string file
Definition: cwb_online.py:385
double df
netcluster * getwc(size_t n)
param: delay index
Definition: network.hh:421
sprintf(tfres,"(1/%g)x(%g) (sec)x(Hz)", 2 *df, df)
int nIFO
Definition: config.hh:102
Long_t id
double factors[FACTORS_MAX]
Definition: config.hh:184
list files
Definition: cwb_online.py:529
float rate
Definition: netpixel.hh:95
double getdata(char type='R', size_t n=0)
Definition: netpixel.hh:56
virtual void resize(unsigned int)
Definition: wavearray.cc:445
int check
wavearray< double > y
Definition: Test10.C:31
void PlotWaveform(TString ifo, wavearray< double > *wfREC, CWB::config *cfg, bool fft=false, bool strain=false)
exit(0)
bool dump
Definition: config.hh:277