Logo coherent WaveBurst  
Library Reference Guide
Logo
CWB_Plugin_ChirpMass.C
Go to the documentation of this file.
1 /*
2 # Copyright (C) 2019 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 #define XIFO 4
20 
21 #pragma GCC system_header
22 
23 #define _USE_EBBH
24 
25 #include "cwb.hh"
26 #include "config.hh"
27 #include "network.hh"
28 #include "wavearray.hh"
29 #include "TString.h"
30 #include "TObjArray.h"
31 #include "TObjString.h"
32 #include "TRandom.h"
33 #include "TComplex.h"
34 #include "TGraphAsymmErrors.h"
35 #include "TMath.h"
36 #include "mdc.hh"
37 #include "cwb2G.hh"
38 #include "watplot.hh"
39 #include "gwavearray.hh"
40 #include "network.hh"
41 #include <fstream>
42 #include <vector>
43 
44 
45 // ---------------------------------------------------------------------------------
46 // DEFINES
47 // ---------------------------------------------------------------------------------
48 
49 #define NCHIRP_MAX 4
50 
51 #define CHI2_THR -2.5 // if negative the chirp pixels are tagged with pix->likelihood=0 after the selection
52 #define TMERGER_CUT 0.1
53 #define ZMAX_THR 0.2
54 
55 #define DUMP_WAVE false // dump rec/inj/wht waveforms to output root file
56 #define DUMP_CLUSTER false // dump pixel cluster to output root file
57 
58 // ---------------------------------------------------------------------------------
59 // USER CONFIG OPTIONS
60 // ---------------------------------------------------------------------------------
61 
62 struct uoptions {
63  float chi2_thr;
64  float tmerger_cut;
65  float zmax_thr;
66 
67  bool dump_wave;
68  bool dump_cluster;
69 };
70 
71 // ---------------------------------------------------------------------------------
72 // FUNCTIONS
73 // ---------------------------------------------------------------------------------
74 
76 
77 void SetOutputFile(network* NET, netevent* &EVT, CWB::config* cfg, bool dump_ebbh);
78 void DumpOutputFile(network* NET, netevent* &EVT, CWB::config* cfg, int ID, int k, int factor);
79 
80 std::vector<netpixel> GetCluster(network* NET, CWB::config* cfg, int lag, int id);
81 
82 void SetEventWindow(CWB::config* cfg, double gps);
83 
85 double GetTimeBoundaries(wavearray<double> x, double P, double& bT, double& eT);
86 double GetFrequencyBoundaries(wavearray<double> x, double P, double& bF, double& eF);
88 
89 std::vector<netpixel> DoPCA(network* NET, CWB::config* cfg, int lag, int id);
90 void ResetPCA(network* NET, CWB::config* cfg, netcluster* pwc, std::vector<netpixel>* vPIX, int ID);
91 
92 std::vector<wavearray<double> > GetInjWaveform(network* NET, netevent* EVT, int id, double factor);
93 std::vector<wavearray<double> > GetRecWaveform(network* NET, netevent* EVT, int id);
94 
95 wavearray<double> GetWaveform(int ifoId, network* NET, int lag, int id, char type, bool shift=true);
96 
100 
102  wavearray<double>* wf2, wavearray<double>* wf3, wavearray<double>* wref, bool fft=false,TString pdir="", double P=0.99);
103 
104 void GetEccentricityIndex(network* NET, int lag, int id);
105 
106 void ReadUserOptions(CWB::config* cfg);
107 void PrintUserOptions();
108 
109 // ----------------------------------------------------
110 // ROOT Output EBBH Parameters
111 // ----------------------------------------------------
112 
113 struct _EBBH { // structure for output estimated parameters
114 
115  int nch; // number of analyzed chirp
116  float ei; // eccentricity index
117  float ch_mass[NCHIRP_MAX]; // chirp mass
118  float ch_tmerger[NCHIRP_MAX]; // merger time
119  float ch_merr[NCHIRP_MAX]; // mchirp error
120  float ch_mchi2[NCHIRP_MAX]; // mchirp chi2
121  float ch_energy[NCHIRP_MAX]; // chirp energy
122 
123  wavearray<double>* wDAT[NIFO_MAX]; // whitened data (in the vREC time range)
124  wavearray<double>* wINJ[NIFO_MAX]; // whiten injected waveforms
125  wavearray<double>* sINJ[NIFO_MAX]; // strain injected waveforms
126  wavearray<double>* wREC[NIFO_MAX]; // whiten injected waveforms
127  wavearray<double>* sREC[NIFO_MAX]; // strain injected waveforms
128 
129  double segGPS; // segment start time
130  float cRATE; // netcluste::rate param
131 
132  std::vector<netpixel> *pCLUSTER; // netcluster
133 };
134 
135 // ---------------------------------------------------------------------------------
136 // Global Variables
137 // ---------------------------------------------------------------------------------
138 
139 uoptions gOPT; // global User Options
140 
141 _EBBH gEBBH; // global output for estimated parameters
142 TTree* gTREE; // output tree file name
143 TString gOUTPUT; // output root file name
145 
146 wavearray<double> gHOT[NIFO_MAX]; // HoT time series used to save whitened data
147 
148 double gINRATE; // input data sampling rate
149 int gRATEANA; // analysis data rate
150 
151 wavearray<double> gwDAT[NIFO_MAX]; // whitened data (in the vREC time range)
152 wavearray<double> gwINJ[NIFO_MAX]; // whiten injected waveforms
153 wavearray<double> gsINJ[NIFO_MAX]; // strain injected waveforms
154 wavearray<double> gwREC[NIFO_MAX]; // whiten injected waveforms
155 wavearray<double> gsREC[NIFO_MAX]; // strain injected waveforms
156 
157 double gSEGGPS; // segment start time
158 float gCRATE; // netcluster::rate
159 std::vector<netpixel> gCLUSTER; // netcluster
160 
161 void
162 CWB_Plugin(TFile* jfile, CWB::config* cfg, network* NET, WSeries<double>* x, TString ifo, int type) {
163 //!MISCELLANEA
164 // Plugin used for eBBH Parameter Estimation
165 
166  if(type==CWB_PLUGIN_CONFIG) {
167  cwb2G* gCWB2G; IMPORT(cwb2G*,gCWB2G)
168 
169  if(gCWB2G->istage!=CWB_STAGE_FULL) {cout<< "CWB_Plugin_ChirpMass - Error : EBBH can be executed only in CWB_STAGE_FULL mode" << endl;exit(1);}
170 
171  ReadUserOptions(cfg); // user config options : read from parPlugin
172 
173  if(cfg->pattern<=0) {
174  cout << "CWB_Plugin_ChirpMass Error : EBBH enable only with pattern>0 !!!" << endl;
175  exit(1);
176  }
177  cfg->outPlugin=true; // disable built-in output root file
178  gINRATE=cfg->inRate; // input data sampling rate
179  gRATEANA=gCWB2G->rateANA; // analysis data rate
180 
181  gEBBH.ei = 1.e-20;
182  gEBBH.nch = NCHIRP_MAX;
183  for(int i=0;i<NCHIRP_MAX;i++) {
184  gEBBH.ch_mass[i] = 0;
185  gEBBH.ch_tmerger[i] = 0;
186  gEBBH.ch_merr[i] = 0;
187  gEBBH.ch_mchi2[i] = 0;
188  gEBBH.ch_energy[i] = 0;
189  }
190  }
191 
192  if(type==CWB_PLUGIN_NETWORK) {
193  PrintUserOptions(); // print config options
194  }
195 
197  // save whitened HoT
198  int ifoID =0; for(int n=0;n<cfg->nIFO;n++) if(ifo==NET->getifo(n)->Name) {ifoID=n;break;}
199  gHOT[ifoID] = *x;
200  }
201 
202  if(type==CWB_PLUGIN_OLIKELIHOOD) { // AFTER EVENT RECONSTRUCTION
203 
204  // import trials
205  int gLRETRY=-1; IMPORT(int,gLRETRY)
206 
207  // import ifactor
208  int gIFACTOR=-1; IMPORT(int,gIFACTOR)
209  double factor = cfg->simulation==3||cfg->simulation==4 ? -gIFACTOR : cfg->factors[gIFACTOR];
210  double ofactor=0;
211  if(cfg->simulation==4) ofactor=-gIFACTOR;
212  else if(cfg->simulation==3) ofactor=-gIFACTOR;
213  else ofactor=factor;
214 
215  int nIFO = NET->ifoListSize(); // number of detectors
216  int K = NET->nLag; // number of time lag
217  int rate = 0; // select all resolutions
218  netevent* EVT;
220 
221  SetOutputFile(NET, EVT, cfg, false); // set output root file
222 
223  for(int k=0; k<K; k++) { // loop over the lags
224 
225  id = NET->getwc(k)->get(const_cast<char*>("ID"), 0, 'L', rate);
226 
227  for(int j=0; j<(int)id.size(); j++) { // loop over cluster index
228 
229  int ID = size_t(id.data[j]+0.5);
230 
231  if(NET->getwc(k)->sCuts[ID-1]!=-1) continue; // skip rejected/processed clusters
232 
233  EVT->output2G(NULL,NET,ID,k,ofactor); // get reconstructed parameters
234 
235  // print event parameters
236  cout << endl;
237  cout << "event parameters : ID -> " << ID << endl;
238  for(int n=0;n<nIFO;n++) printf("rec_time %s : %.4f\n",NET->ifoName[n], EVT->time[n]);
239  cout << "rec_theta : " << EVT->theta[0] << " rec_phi : " << EVT->phi[0] << endl;
240  cout << "SNRnet : " << sqrt(EVT->likelihood) << " netcc[0] : " << EVT->netcc[0]
241  << " rho[0] : " << EVT->rho[0] << " size : " << EVT->size[0] << endl;
242 
243  gSEGGPS = EVT->gps[0]; // segment start time
244 
245  // get waveforms
246  if(cfg->simulation) {
247  std::vector<wavearray<double> > vINJ; // injected
248  vINJ = GetInjWaveform(NET, EVT, ID, factor); // get injected waveform
249  if(vINJ.size()!=nIFO) {
250  cout << "CWB_Plugin_ChirpMass Error : Injection Waveform Not Found !!!" << endl;
251  exit(1);
252  }
253  for(int n=0;n<nIFO;n++) gwINJ[n] = vINJ[n];
254  }
255  for(int n=0;n<nIFO;n++) {
256  gwREC[n] = GetWaveform(n, NET, k, ID, 'S'); // get whiten reconstructed waveform
257  gsREC[n] = GetWaveform(n, NET, k, ID, 's'); // get strain reconstructed waveform
258  gwDAT[n] = GetAlignedWaveform(&gHOT[n], &gwREC[n]); // get whitened data (in the vREC time range)
259  }
260  gCRATE = NET->getwc(k)->rate;
261  gCLUSTER = GetCluster(NET, cfg, k, ID); // get netcluster
262 
263  GetEccentricityIndex(NET, k, ID); // get eccentricity index
264 
265  SetOutputFile(NET, EVT, cfg, true); // set output root file
266  DumpOutputFile(NET, EVT, cfg, ID, k, ofactor); // dump event to output root file
267 
268  for(int n=0;n<nIFO;n++) {
269  detector* pD = NET->getifo(n);
270  if(!cfg->simulation) ClearWaveforms(pD); // release waveform memory
271  }
272  }
273  }
274 
275  jfile->cd();
276  if(EVT) delete EVT;
277 
278  while(!gCLUSTER.empty()) gCLUSTER.pop_back();
279  gCLUSTER.clear(); std::vector<netpixel>().swap(gCLUSTER);
280  }
281 
282  return;
283 }
284 
285 std::vector<netpixel> DoPCA(network* NET, CWB::config* cfg, int lag, int id) {
286 
287  double ee;
288 
289  size_t nIFO = NET->ifoList.size();
290 
291  float En = 2*NET->acor*NET->acor*nIFO; // network energy threshold in the sky loop
292 
293  int size = NET->a_00.size();
294  int f_ = NIFO/4;
295 
296  netpixel* pix;
297  netcluster* pwc = NET->getwc(lag);
298  std::vector<netpixel> vPIX;
299 
300  std::vector<int> pI = NET->wdmMRA.getXTalk(pwc, id, false); // buffer for pixel IDs
301  int V = pI.size();
302  if(V>cfg->BATCH) return vPIX; // attach TD amp to pixels < V
303 
304  wavearray<float> xi(size); xi=0; // PC 00 array
305  wavearray<float> XI(size); XI=0; // PC 90 array
306 
307  __m128* _xi = (__m128*) xi.data; // set pointer to PC 00 array
308  __m128* _XI = (__m128*) XI.data; // set pointer to PC 90 array
309 
310  __m128* _aa = (__m128*) NET->a_00.data; // set pointer to 00 array
311  __m128* _AA = (__m128*) NET->a_90.data; // set pointer to 90 array
312 
313  int nPC = 0;
314  for(int j=0; j<V; j++) {
315  int jf = j*f_; // source sse pointer increment
316  _sse_zero_ps(_xi+jf); // zero MRA amplitudes
317  _sse_zero_ps(_XI+jf); // zero MRA amplitudes
318  ee = _sse_abs_ps(_aa+jf,_AA+jf); // total pixel energy / quadrature
319  if(ee>En) nPC++; else ee=0.; // count core pixels
320  NET->rNRG.data[j] = ee; // init residual energy array
321  NET->pNRG.data[j] = NET->rNRG.data[j]; // init residual energy array
322  }
323 
324  nPC = NET->_sse_mra_ps(xi.data,XI.data,En,nPC); // get principal components
325 
326  for(int j=0; j<V; j++) { // loop over principal components
327  pix = pwc->getPixel(id,pI[j]);
328  vPIX.push_back(*pix); // save original pixels
329  pix->core = false;
330  ee = NET->pNRG.data[j]; // total pixel energy
331  if(ee<En) continue;
332  pix->core = true;
333  for(int i=0; i<nIFO; i++) {
334  pix->setdata(double(xi.data[j*NIFO+i]),'S',i); // store 00 whitened response PC
335  pix->setdata(double(XI.data[j*NIFO+i]),'P',i); // store 90 whitened response PC
336  }
337  }
338 
339  return vPIX;
340 }
341 
342 void ResetPCA(network* NET, CWB::config* cfg, netcluster* pwc, std::vector<netpixel>* vPIX, int ID) {
343 
344  std::vector<int> pI = NET->wdmMRA.getXTalk(pwc, ID, false); // buffer for pixel IDs
345  int V = pI.size();
346  if(V>cfg->BATCH) return; // attach TD amp to pixels < V
347  for(int j=0; j<V; j++) {
348  netpixel* pix = pwc->getPixel(ID,pI[j]);
349  *pix = (*vPIX)[j];
350  }
351 
352  while(!vPIX->empty()) vPIX->pop_back();
353  vPIX->clear(); std::vector<netpixel>().swap(*vPIX);
354 }
355 
356 std::vector<wavearray<double> > GetInjWaveform(network* NET, netevent* EVT, int ID, double factor) {
357 
358  std::vector<wavearray<double> > xINJ; // injected
359 
360  int nIFO = NET->ifoListSize(); // number of detectors
361 
362  double recTime = EVT->time[0]; // rec event time det=0
363 
364  injection INJ(nIFO);
365  // get inj ID
366  int M = NET->mdc__IDSize();
367  double injTime = 1.e12;
368  int injID = -1;
369  for(int m=0; m<M; m++) {
370  int mdcID = NET->getmdc__ID(m);
371  double T = fabs(recTime - NET->getmdcTime(mdcID));
372  if(T<injTime && INJ.fill_in(NET,mdcID)) {
373  injTime = T;
374  injID = mdcID;
375  }
376  }
377 
378  if(INJ.fill_in(NET,injID)) { // get inj parameters
379 
380  wavearray<double>** pwfINJ = new wavearray<double>*[nIFO];
381 
382  // extract whitened injected & reconstructed waveforms
383  for(int n=0; n<nIFO; n++) {
384 
385  detector* pd = NET->getifo(n);
386 
387  pwfINJ[n] = INJ.pwf[n];
388  if (pwfINJ[n]==NULL) {
389  cout << "CWB_Plugin_ChirpMass.C : Error : Injected waveform not saved !!! : detector "
390  << NET->ifoName[n] << endl;
391  continue;
392  }
393 
394  double rFactor = 1.;
395  rFactor *= factor>0 ? factor : 1;
396  wavearray<double>* wf = pwfINJ[n];
397  *wf*=rFactor; // injected wf is multiplied by the custom factor
398  xINJ.push_back(*wf);
399  *wf*=1./rFactor; // restore injected amplitude
400  }
401  delete [] pwfINJ;
402  }
403 
404  return xINJ;
405 }
406 
407 std::vector<wavearray<double> > GetRecWaveform(network* NET, netevent* EVT, int ID) {
408 
409  std::vector<wavearray<double> > xREC; // reconstructed
410 
411  int nIFO = NET->ifoListSize(); // number of detectors
412 
413  wavearray<double>** pwfREC = new wavearray<double>*[nIFO];
414  for(int n=0; n<nIFO; n++) {
415 
416  detector* pd = NET->getifo(n);
417  int idSize = pd->RWFID.size();
418  int wfIndex=-1;
419  for(int mm=0; mm<idSize; mm++) if (pd->RWFID[mm]==ID) wfIndex=mm;
420 
421  if(wfIndex<0) {
422  cout << "CWB_Plugin_ChirpMass.C : Error : Reconstructed waveform not saved !!! : ID -> "
423  << ID << " : detector " << NET->ifoName[n] << endl;
424  continue;
425  }
426  if(wfIndex>=0) pwfREC[n] = pd->RWFP[wfIndex];
427 
428  wavearray<double>* wf = pwfREC[n];
429  xREC.push_back(*wf);
430  }
431  delete [] pwfREC;
432 
433  return xREC;
434 }
435 
436 wavearray<double> GetWaveform(int ifoId, network* NET, int lag, int id, char type, bool shift) {
437 
438  if(type!='S' && type!='s' && type!='W' && type!='w' && type!='N' && type!='n') {
439  cout << "GetWaveform : not valid type : Abort" << endl; exit(1);
440  }
441 
442  netcluster* pwc = NET->getwc(lag);
443  detector* pd = NET->getifo(ifoId);
445 
446  if(type=='S' || type=='s') {
447  NET->getMRAwave(id,lag,type,0,shift); // reconstruct whitened/strain shifted pd->waveForm
448  wf = &pd->waveForm;
449  wf->start(pwc->start+pd->waveForm.start());
450  }
451 
452  if(type=='W' || type=='w') {
453  NET->getMRAwave(id,lag,type,0,shift); // reconstruct+noise whitened shifted pd->waveBand
454  wf = &pd->waveBand;
455  wf->start(pwc->start+pd->waveBand.start());
456  }
457 
458  if(type=='N' || type=='n') {
459  NET->getMRAwave(id,lag,'W',0,shift); // reconstruct+noise whitened shifted pd->waveBand
460  NET->getMRAwave(id,lag,'S',0,shift); // reconstruct whitened shifted pd->waveForm
461  pd->waveNull = pd->waveBand;
462  pd->waveNull-= pd->waveForm;
463  wf = &pd->waveNull;
464  wf->start(pwc->start+pd->waveNull.start());
465  }
466 
467  return *wf;
468 }
469 
471  wavearray<double>* wf2, wavearray<double>* wf3, wavearray<double>* wref, bool fft, TString pdir, double P) {
472 
473  watplot PTS(const_cast<char*>("ptspe"),200,20,800,500);
474 
475  //cout << "Print " << ofname << endl;
476  double tmin=1.e20;
477  if(wref==NULL) return;
478  if(wf1==NULL) return;
479  else tmin=wf1->start();
480  if(wf2!=NULL) if(wf2->start()<tmin) tmin=wf2->start();
481  if(wf3!=NULL) if(wf3->start()<tmin) tmin=wf3->start();
482 
483  wf1->start(wf1->start()-tmin);
484  if(wf2!=NULL) if(wf2!=wf1) wf2->start(wf2->start()-tmin);
485  if(wf3!=NULL) if(wf3!=wf1 && wf3!=wf2) wf3->start(wf3->start()-tmin);
486  if(wref!=wf1 && wref!=wf2 && wref!=wf3) wref->start(wref->start()-tmin);
487 
488  if(fft) {
489  PTS.plot(wf1, const_cast<char*>("AL"), kRed, 0, 0, true, cfg->fLow, cfg->fHigh);
490  } else {
491  double bT, eT;
492  GetTimeBoundaries(*wref, P, bT, eT);
493  PTS.plot(wf1, const_cast<char*>("AL"), kRed, bT, eT);
494  PTS.graph[0]->GetXaxis()->SetRangeUser(bT, eT);
495  }
496  PTS.graph[0]->SetLineWidth(1);
497  PTS.graph[0]->SetTitle(title);
498 
499  TString xtitle = TString::Format("Time (sec) : GPS OFFSET = %.3f",tmin);
500  PTS.graph[0]->GetXaxis()->SetTitle(xtitle);
501 
502  if(wf2!=NULL) {
503  if(fft) {
504  if(wf2!=NULL) PTS.plot(wf2, const_cast<char*>("SAME"), kBlue, 0, 0, true, cfg->fLow, cfg->fHigh);
505  } else {
506  if(wf2!=NULL) PTS.plot(wf2, const_cast<char*>("SAME"), kBlue, 0, 0);
507  }
508  if(wf2!=NULL) PTS.graph[1]->SetLineWidth(1);
509  }
510 
511  if(wf3!=NULL) {
512  if(fft) {
513  if(wf3!=NULL) PTS.plot(wf3, const_cast<char*>("SAME"), kGreen+2, 0, 0, true, cfg->fLow, cfg->fHigh);
514  } else {
515  if(wf3!=NULL) PTS.plot(wf3, const_cast<char*>("SAME"), kGreen+2, 0, 0);
516  }
517  if(wf3!=NULL) PTS.graph[2]->SetLineWidth(1);
518  }
519 
520  wf1->start(wf1->start()+tmin);
521  if(wf2!=NULL) if(wf2!=wf1) wf2->start(wf2->start()+tmin);
522  if(wf3!=NULL) if(wf3!=wf1 && wf3!=wf2) wf3->start(wf3->start()+tmin);
523  if(wref!=wf1 && wref!=wf2 && wref!=wf3) wref->start(wref->start()+tmin);
524 
525  if(fft) {
526  PTS.canvas->SetLogx();
527  PTS.canvas->SetLogy();
528  }
529 
530  if(ofname!="") {
531  ofname = TString(pdir)+TString("/")+ofname;
532  PTS.canvas->Print(ofname);
533  cout << "write : " << ofname << endl;
534  }
535 }
536 
538 
539  double R=wf1->rate();
540 
541  double b_wf1 = wf1->start();
542  double e_wf1 = wf1->start()+wf1->size()/R;
543  double b_wf2 = wf2->start();
544  double e_wf2 = wf2->start()+wf2->size()/R;
545 
546  int o_wf1 = b_wf1>b_wf2 ? 0 : int((b_wf2-b_wf1)*R+0.5);
547  int o_wf2 = b_wf1<b_wf2 ? 0 : int((b_wf1-b_wf2)*R+0.5);
548 
549  double startXCOR = b_wf1>b_wf2 ? b_wf1 : b_wf2;
550  double endXCOR = e_wf1<e_wf2 ? e_wf1 : e_wf2;
551  int sizeXCOR = int((endXCOR-startXCOR)*R+0.5);
552 
553  wavearray<double> wfdif(sizeXCOR);
554  wfdif=0.;
555  wfdif.rate(R);
556  wfdif.start(b_wf1+double(o_wf1)/R);
557 
558  for(int i=0;i<sizeXCOR;i++) wfdif[i] = wf1->data[i+o_wf1] - wf2->data[i+o_wf2];
559 
560  return wfdif;
561 }
562 
563 void SetOutputFile(network* NET, netevent* &EVT, CWB::config* cfg, bool dump_ebbh) {
564 
565  // import slagShift
566  float* gSLAGSHIFT=NULL; IMPORT(float*,gSLAGSHIFT)
567 
568  int nIFO = NET->ifoListSize(); // number of detectors
569 
570  // search output root file in the system list
571  TFile* froot = NULL;
572  TList *files = (TList*)gROOT->GetListOfFiles();
573  gOUTPUT="";
574  if (files) {
575  TIter next(files);
576  TSystemFile *file;
577  TString fname;
578  bool check=false;
579  while ((file=(TSystemFile*)next())) {
580  fname = file->GetName();
581  // set output root file as the current file
582  if(fname.Contains("wave_")) {
583  froot=(TFile*)file;froot->cd();
584  gOUTPUT=fname;
585  gOUTPUT.ReplaceAll(".root.tmp",".txt");
586  //cout << "output file name : " << fname << endl;
587  }
588  }
589  if(!froot) {
590  cout << "CWB_Plugin_ChirpMass.C : Error - output root file not found" << endl;
591  gSystem->Exit(1);
592  }
593  } else {
594  cout << "CWB_Plugin_ChirpMass.C : Error - output root file not found" << endl;
595  gSystem->Exit(1);
596  }
597 
598  gTREE = (TTree *) froot->Get("waveburst");
599  if(gTREE!=NULL) {
600  EVT = new netevent(gTREE,nIFO);
601  } else {
602  EVT = new netevent(nIFO);
603  gTREE = EVT->setTree();
604  }
605  EVT->setSLags(gSLAGSHIFT); // set slags into netevent
606  EVT->Psave=cfg->Psave;
607 
608  static bool ebbh_tree_init = false;
609  if(dump_ebbh) {
610 
611  gEBBH.segGPS = gSEGGPS;
612  if(gOPT.dump_wave) {
613  if(cfg->simulation) for(int n=0;n<nIFO;n++) gEBBH.wINJ[n] = &gwINJ[n];
614  for(int n=0;n<nIFO;n++) gEBBH.wDAT[n] = &gwDAT[n];
615  for(int n=0;n<nIFO;n++) gEBBH.wREC[n] = &gwREC[n];
616  for(int n=0;n<nIFO;n++) gEBBH.sREC[n] = &gsREC[n];
617  }
618  if(gOPT.dump_cluster) {
619  gEBBH.cRATE = gCRATE;
620  gEBBH.pCLUSTER = &gCLUSTER;
621  }
622 
623  if(!ebbh_tree_init) {
624 
625  ebbh_tree_init = true;
626 
627  gTREE->Branch("ebbh_ei", &gEBBH.ei,"ebbh_ei/F");
628  gTREE->Branch("ebbh_nch",&gEBBH.nch,"ebbh_nch/I");
629 
630  gTREE->Branch("ebbh_ch_mass", gEBBH.ch_mass,TString::Format("ebbh_ch_mass[%i]/F",NCHIRP_MAX));
631  gTREE->Branch("ebbh_ch_tmerger",gEBBH.ch_tmerger,TString::Format("ebbh_ch_tmerger[%i]/F",NCHIRP_MAX));
632  gTREE->Branch("ebbh_ch_merr", gEBBH.ch_merr,TString::Format("ebbh_ch_merr[%i]/F",NCHIRP_MAX));
633  gTREE->Branch("ebbh_ch_mchi2", gEBBH.ch_mchi2,TString::Format("ebbh_ch_mchi2[%i]/F",NCHIRP_MAX));
634  gTREE->Branch("ebbh_ch_energy", gEBBH.ch_energy,TString::Format("ebbh_ch_energy[%i]/F",NCHIRP_MAX));
635 
636  gTREE->Branch("ebbh_seggps",&gEBBH.segGPS,"ebbh_seggps/D");
637  if(gOPT.dump_wave) {
638  for(int n=0;n<nIFO;n++) {
639  if(cfg->simulation) gTREE->Branch(TString::Format("ebbh_wINJ_%d",n).Data(),"wavearray<double>",&gEBBH.wINJ[n],32000,0);
640  gTREE->Branch(TString::Format("ebbh_wDAT_%d",n).Data(),"wavearray<double>",&gEBBH.wDAT[n],32000,0);
641  gTREE->Branch(TString::Format("ebbh_wREC_%d",n).Data(),"wavearray<double>",&gEBBH.wREC[n],32000,0);
642  gTREE->Branch(TString::Format("ebbh_sREC_%d",n).Data(),"wavearray<double>",&gEBBH.sREC[n],32000,0);
643  }
644  }
645  if(gOPT.dump_cluster) {
646  gTREE->Branch("ebbh_crate",&gEBBH.cRATE,"ebbh_crate/F");
647  gTREE->Branch("ebbh_cluster","vector<netpixel>",&gEBBH.pCLUSTER,32000,0);
648  }
649 
650  } else {
651 
652  gTREE->SetBranchAddress("ebbh_ei", &gEBBH.ei);
653  gTREE->SetBranchAddress("ebbh_nch",&gEBBH.nch);
654 
655  gTREE->SetBranchAddress("ebbh_ch_mass", gEBBH.ch_mass);
656  gTREE->SetBranchAddress("ebbh_ch_tmerger",gEBBH.ch_tmerger);
657  gTREE->SetBranchAddress("ebbh_ch_merr", gEBBH.ch_merr);
658  gTREE->SetBranchAddress("ebbh_ch_mchi2", gEBBH.ch_mchi2);
659  gTREE->SetBranchAddress("ebbh_ch_energy", gEBBH.ch_energy);
660 
661  gTREE->SetBranchAddress("ebbh_seggps",&gEBBH.segGPS);
662  if(gOPT.dump_wave) {
663  for(int n=0;n<nIFO;n++) {
664  if(cfg->simulation) gTREE->SetBranchAddress(TString::Format("ebbh_wINJ_%d",n).Data(),&gEBBH.wINJ[n]);
665  gTREE->SetBranchAddress(TString::Format("ebbh_wDAT_%d",n).Data(),&gEBBH.wDAT[n]);
666  gTREE->SetBranchAddress(TString::Format("ebbh_wREC_%d",n).Data(),&gEBBH.wREC[n]);
667  gTREE->SetBranchAddress(TString::Format("ebbh_sREC_%d",n).Data(),&gEBBH.sREC[n]);
668  }
669  }
670  if(gOPT.dump_cluster) {
671  gTREE->SetBranchAddress("ebbh_crate",&gEBBH.cRATE);
672  gTREE->SetBranchAddress("ebbh_cluster",&gEBBH.pCLUSTER);
673  }
674  }
675  }
676 }
677 
678 void DumpOutputFile(network* NET, netevent* &EVT, CWB::config* cfg, int ID, int k, int factor) {
679 
680  if(cfg->dump) EVT->dopen(gOUTPUT.Data(),const_cast<char*>("a"),false);
681  EVT->output2G(gTREE,NET,ID,k,factor); // get reconstructed parameters
682  if(cfg->dump) EVT->dclose();
683  if(!cfg->cedDump) NET->getwc(k)->sCuts[ID-1]=1; // mark as processed
684 }
685 
687 
688  wavearray<double> wf = *wf2;
689  wf=0;
690 
691  if(wf1==NULL) return wf;
692  if(wf1->size()==0) return wf;
693 
694  double R=wf1->rate();
695 
696  double b_wf1 = wf1->start();
697  double e_wf1 = wf1->start()+wf1->size()/R;
698  double b_wf2 = wf2->start();
699  double e_wf2 = wf2->start()+wf2->size()/R;
700 
701  int o_wf1 = b_wf1>b_wf2 ? 0 : int((b_wf2-b_wf1)*R+0.5);
702  int o_wf2 = b_wf1<b_wf2 ? 0 : int((b_wf1-b_wf2)*R+0.5);
703 
704  double startXCOR = b_wf1>b_wf2 ? b_wf1 : b_wf2;
705  double endXCOR = e_wf1<e_wf2 ? e_wf1 : e_wf2;
706  int sizeXCOR = int((endXCOR-startXCOR)*R+0.5);
707 
708  for(int i=0;i<sizeXCOR;i++) wf[i+o_wf2] = wf1->data[i+o_wf1];
709 
710  return wf;
711 }
712 
714 
715  double a;
716  double E=0.,T=0.;
717  int size=(int)x.size();
718  double rate=x.rate();
719  for(int j=0;j<size;j++) {
720  a = x[j];
721  T += a*a*j/rate; // central time
722  E += a*a; // energy
723  }
724  T = E>0 ? T/E : 0.5*size/rate;
725 
726  return T;
727 }
728 
729 double GetTimeBoundaries(wavearray<double> x, double P, double& bT, double& eT) {
730 
731  if(P<0) P=0;
732  if(P>1) P=1;
733 
734  int N = x.size();
735 
736  double E = 0; // signal energy
737  double avr = 0; // average
738  for(int i=0;i<N;i++) {avr+=i*x[i]*x[i]; E+=x[i]*x[i];}
739  int M=int(avr/E); // central index
740 
741  // search range which contains percentage P of the total energy E
742  int jB=0;
743  int jE=N-1;
744  double a,b;
745  double sum = ((M>=0)&&(M<N)) ? x[M]*x[M] : 0.;
746  for(int j=1; j<N; j++) {
747  a = ((M-j>=0)&&(M-j<N)) ? x[M-j] : 0.;
748  b = ((M+j>=0)&&(M+j<N)) ? x[M+j] : 0.;
749  if(a) jB=M-j;
750  if(b) jE=M+j;
751  sum += a*a+b*b;
752  if(sum/E > P) break;
753  }
754 
755  bT = x.start()+jB/x.rate();
756  eT = x.start()+jE/x.rate();
757 
758  return eT-bT;
759 }
760 
762 
763  double a;
764  double E=0.,F=0.;
765  int size=(int)x.size();
766  double rate=x.rate();
767  x.FFTW(1);
768  double dF=(rate/(double)size)/2.;
769  for(int j=0;j<size/2;j+=2) {
770  a = x[j]*x[j]+x[j+1]*x[j+1];
771  F += a*j*dF; // central frequency
772  E += a; // energy
773  }
774  F = E>0 ? F/E : 0.5*rate;
775 
776  return F;
777 }
778 
779 double GetFrequencyBoundaries(wavearray<double> x, double P, double& bF, double& eF) {
780 
781  if(P<0) P=0;
782  if(P>1) P=1;
783 
784  int N = x.size();
785 
786  x.FFTW(1);
787 
788  double E = 0; // signal energy
789  double avr = 0; // average
790  for(int i=0;i<N;i++) {avr+=i*x[i]*x[i]; E+=x[i]*x[i];}
791  int M=int(avr/E); // central index
792 
793  // search range which contains percentage P of the total energy E
794  int jB=0;
795  int jE=N-1;
796  double a,b;
797  double sum = ((M>=0)&&(M<N)) ? x[M]*x[M] : 0.;
798  for(int j=1; j<N; j++) {
799  a = ((M-j>=0)&&(M-j<N)) ? x[M-j] : 0.;
800  b = ((M+j>=0)&&(M+j<N)) ? x[M+j] : 0.;
801  if(a) jB=M-j;
802  if(b) jE=M+j;
803  sum += a*a+b*b;
804  if(sum/E > P) break;
805  }
806 
807  double dF=(x.rate()/(double)x.size())/2.;
808 
809  bF = jB*dF;
810  eF = jE*dF;
811 
812  return eF-bF;
813 }
814 
815 void
817 
818  int n;
819 
820  n = ifo->IWFP.size();
821  for (int i=0;i<n;i++) {
823  delete wf;
824  }
825  ifo->IWFP.clear();
826  ifo->IWFID.clear();
827 
828  n = ifo->RWFP.size();
829  for (int i=0;i<n;i++) {
831  delete wf;
832  }
833  ifo->RWFP.clear();
834  ifo->RWFID.clear();
835 }
836 
838 
839  wavearray<double> wf = *wf1;
840 
841  if(wf1==NULL) return wf;
842  if(wf1->size()==0) return wf;
843 
844  double R=wf1->rate();
845 
846  double b_wf1 = wf1->start();
847  double e_wf1 = wf1->start()+wf1->size()/R;
848  double b_wf2 = wf2->start();
849  double e_wf2 = wf2->start()+wf2->size()/R;
850 
851  int o_wf1 = b_wf1>b_wf2 ? 0 : int((b_wf2-b_wf1)*R+0.5);
852  int o_wf2 = b_wf1<b_wf2 ? 0 : int((b_wf1-b_wf2)*R+0.5);
853 
854  double startXCOR = b_wf1>b_wf2 ? b_wf1 : b_wf2;
855  double endXCOR = e_wf1<e_wf2 ? e_wf1 : e_wf2;
856  int sizeXCOR = int((endXCOR-startXCOR)*R+0.5);
857 
858  for(int i=0;i<sizeXCOR;i++) wf[i+o_wf1] += wf2->data[i+o_wf2];
859 
860  return wf;
861 }
862 
863 void SetEventWindow(CWB::config* cfg, double gps) {
864 
865  if(gps<0) return;
866 
867  // dq file list
868  // {ifo, dqcat_file, dqcat[0/1/2], shift[sec], inverse[false/true], 4columns[true/false]}
869 
870  for(int n=0; n<cfg->nIFO; n++) {
871 
872  strcpy(cfg->DQF[cfg->nDQF].ifo, cfg->ifo[n]);
873  sprintf(cfg->DQF[cfg->nDQF].file, "%s/%s_%s.gps_%d",cfg->tmp_dir,cfg->ifo[n],cfg->data_label,int(gps));
874  cfg->DQF[cfg->nDQF].cat = CWB_CAT2;
875  cfg->DQF[cfg->nDQF].shift = 0.;
876  cfg->DQF[cfg->nDQF].invert = false;
877  cfg->DQF[cfg->nDQF].c4 = true;
878  cfg->nDQF++;
879 
880  cout << cfg->DQF[cfg->nDQF-1].file << endl;
881 
882  ofstream out;
883  out.open(cfg->DQF[cfg->nDQF-1].file,ios::out);
884  cout << "Write file : " << cfg->DQF[cfg->nDQF-1].file << endl;
885  if (!out.good()) {cout << "Error Opening File : " << cfg->DQF[cfg->nDQF-1].file << endl;exit(1);}
886  out.precision(14);
887  int istart = int(gps)-cfg->iwindow;
888  int istop = int(gps)+cfg->iwindow;
889  out << "1 " << istart << " " << istop << " " << 2*cfg->iwindow << endl;
890  out.close();
891  }
892 }
893 
894 std::vector<netpixel> GetCluster(network* NET, CWB::config* cfg, int lag, int id) {
895 
896  netpixel* pix;
897  netcluster* pwc = NET->getwc(lag);
898  std::vector<netpixel> vPIX;
899 
900  std::vector<int> pI = NET->wdmMRA.getXTalk(pwc, id, false); // buffer for pixel IDs
901 
902  int V = pI.size();
903  for(int j=0; j<V; j++) { // loop over principal components
904  pix = pwc->getPixel(id,pI[j]);
905  vPIX.push_back(*pix); // save original pixels
906  }
907 
908  return vPIX;
909 }
910 
911 void GetEccentricityIndex(network* NET, int lag, int id) {
912 
913  netcluster* pwc = NET->getwc(lag);
914 
915  double chi2_thr = gOPT.chi2_thr;
916  double tmerger_cut = gOPT.tmerger_cut;
917  double zmax_thr = gOPT.zmax_thr;
918 
919  clusterdata CD = pwc->cData[id-1]; // save loudest chirp
920 
921  for(int i=0;i<NCHIRP_MAX;i++) {
922 
923  double echirp = pwc->mchirp(id,chi2_thr,tmerger_cut,zmax_thr);
924 
925  clusterdata* pcd = &(pwc->cData[id-1]);
926  TF1* fit = &(pcd->fit);
927 
928  gEBBH.ch_mass[i] = pcd->mchirp;
929  gEBBH.ch_tmerger[i] = -fit->GetParameter(1)/fit->GetParameter(0);
930  gEBBH.ch_merr[i] = pcd->mchirperr;
931  gEBBH.ch_mchi2[i] = pcd->chi2chirp;
932  gEBBH.ch_energy[i] = echirp;
933  }
934 
935  // compute the eccentrycity index
936  int imax=0;
937  float emax=0;
938  // find chirp with max energy
939  for(int i=0;i<NCHIRP_MAX;i++) if(gEBBH.ch_energy[i]>emax) {emax=gEBBH.ch_energy[i];imax=i;}
940  int imax2=imax;
941  // find chirp with highest chirp mass and energy > (max chirp energy)/4
942  for(int i=0;i<NCHIRP_MAX;i++) {
943  if(i==imax) continue;
944  if(gEBBH.ch_mass[i]>0) {
945  if(gEBBH.ch_mass[i]>gEBBH.ch_mass[imax]) if(gEBBH.ch_energy[i]>gEBBH.ch_energy[imax]/4.) imax2=i;
946  }
947  }
948  imax=imax2;
949  gEBBH.ei=0;
950  for(int i=0;i<NCHIRP_MAX;i++) {
951  if(gEBBH.ch_mass[i]>0) {
952  if(gEBBH.ch_mass[i]<gEBBH.ch_mass[imax]) gEBBH.ei+=gEBBH.ch_energy[i];
953  if(gEBBH.ch_mass[i]>gEBBH.ch_mass[imax]) gEBBH.ei-=gEBBH.ch_energy[i];
954  }
955  }
956  gEBBH.ei/=gEBBH.ch_energy[imax];
957  gEBBH.ei*=100;
958 
959  cout << endl;
960  cout << "---------------------------------------------------" << endl;
961  for(int i=0;i<NCHIRP_MAX;i++) cout << "mchirp " << i << " -> " << gEBBH.ch_mass[i] << "\t\t echirp -> " << gEBBH.ch_energy[i] << endl;
962  cout << "Eccentricity Index : " << gEBBH.ei << endl;
963  cout << "---------------------------------------------------" << endl;
964  cout << endl;
965 
966  pwc->cData[id-1] = CD; // restore loudest chirp
967 }
968 
970 
971  TString options = cfg->parPlugin;
972 
973  gOPT.chi2_thr = CHI2_THR;
974  gOPT.tmerger_cut = TMERGER_CUT;
975  gOPT.zmax_thr = ZMAX_THR;
976 
977  gOPT.dump_wave = DUMP_WAVE;
978  gOPT.dump_cluster = DUMP_CLUSTER;
979 
980  if(options.CompareTo("")!=0) {
981  cout << options << endl;
982  if(!options.Contains("--")) { // parameters are used only by cwb_inet
983 
984  TObjArray* token = TString(options).Tokenize(TString(' '));
985  for(int j=0;j<token->GetEntries();j++){
986 
987  TObjString* tok = (TObjString*)token->At(j);
988  TString stok = tok->GetString();
989 
990  if(stok.Contains("ebbh_chi2_thr=")) {
991  TString ebbh_chi2_thr=stok;
992  ebbh_chi2_thr.Remove(0,ebbh_chi2_thr.Last('=')+1);
993  if(ebbh_chi2_thr.IsFloat()) gOPT.chi2_thr=ebbh_chi2_thr.Atof();
994  }
995 
996  if(stok.Contains("ebbh_tmerger_cut=")) {
997  TString ebbh_tmerger_cut=stok;
998  ebbh_tmerger_cut.Remove(0,ebbh_tmerger_cut.Last('=')+1);
999  if(ebbh_tmerger_cut.IsFloat()) gOPT.tmerger_cut=ebbh_tmerger_cut.Atof();
1000  }
1001 
1002  if(stok.Contains("ebbh_zmax_thr=")) {
1003  TString ebbh_zmax_thr=stok;
1004  ebbh_zmax_thr.Remove(0,ebbh_zmax_thr.Last('=')+1);
1005  if(ebbh_zmax_thr.IsFloat()) gOPT.zmax_thr=ebbh_zmax_thr.Atof();
1006  }
1007 
1008  if(stok.Contains("ebbh_dump_wave=")) {
1009  TString opt=stok;
1010  opt.Remove(0,opt.Last('=')+1);
1011  if(opt=="false") gOPT.dump_wave=false;
1012  if(opt=="true") gOPT.dump_wave=true;
1013  }
1014 
1015  if(stok.Contains("ebbh_dump_cluster=")) {
1016  TString opt=stok;
1017  opt.Remove(0,opt.Last('=')+1);
1018  if(opt=="false") gOPT.dump_cluster=false;
1019  if(opt=="true") gOPT.dump_cluster=true;
1020  }
1021  }
1022  }
1023  }
1024 }
1025 
1027 
1028  cout << "-----------------------------------------" << endl;
1029  cout << "eBBH Plugin config options" << endl;
1030  cout << "-----------------------------------------" << endl << endl;
1031  cout << endl;
1032  cout << "chi2_thr : " << gOPT.chi2_thr << endl;
1033  cout << "tmerger_cut : " << gOPT.tmerger_cut << endl;
1034  cout << "zmax_thr : " << gOPT.zmax_thr << endl;
1035  cout << "dump_wave : " << gOPT.dump_wave << endl;
1036  cout << "dump_cluster : " << gOPT.dump_cluster << endl;
1037  cout << endl;
1038 
1039 }
1040 
double GetTimeBoundaries(wavearray< double > x, double P, double &bT, double &eT)
monster wdmMRA
list of pixel pointers for MRA
Definition: network.hh:651
wavearray< double > gwREC[NIFO_MAX]
std::vector< char * > ifoName
Definition: network.hh:609
detector * getifo(size_t n)
param: detector index
Definition: network.hh:436
CWB::config * cfg
size_t rateANA
Definition: cwb.hh:279
static float _sse_abs_ps(__m128 *_a)
Definition: watsse.hh:137
void PlotWaveform(TString ofname, TString title, CWB::config *cfg, wavearray< double > *wf1, wavearray< double > *wf2, wavearray< double > *wf3, wavearray< double > *wref, bool fft=false, TString pdir="", double P=0.99)
double iwindow
Definition: config.hh:200
wavearray< double > GetAlignedWaveform(wavearray< double > *wf1, wavearray< double > *wf2)
#define NIFO
Definition: wat.hh:74
size_t nLag
Definition: network.hh:573
Float_t * rho
biased null statistics
Definition: netevent.hh:112
char xtitle[1024]
std::vector< netpixel > gCLUSTER
double M
Definition: DrawEBHH.C:13
void CWB_Plugin(TFile *jfile, CWB::config *cfg, network *NET, WSeries< double > *x, TString ifo, int type)
std::vector< wavearray< double > > GetRecWaveform(network *NET, netevent *EVT, int id)
#define ZMAX_THR
std::vector< wavearray< double > > wREC[MAX_TRIALS]
void setSLags(float *slag)
Definition: netevent.cc:426
std::vector< double > * getmdcTime()
Definition: network.hh:422
wavearray< double > GetDifWaveform(wavearray< double > *wf1, wavearray< double > *wf2)
virtual void rate(double r)
Definition: wavearray.hh:141
float factor
TString opt
wavearray< double > a(hp.size())
int n
Definition: cwb_net.C:28
static void _sse_zero_ps(__m128 *_p)
Definition: watsse.hh:44
bool invert
Definition: Toolbox.hh:88
wavearray< float > a_90
buffer for cluster sky 00 amplitude
Definition: network.hh:653
TTree * setTree()
Definition: netevent.cc:434
TString("c")
char ifo[32]
Definition: Toolbox.hh:84
int ID
Definition: TestMDC.C:70
Int_t * size
cluster volume
Definition: netevent.hh:80
float gCRATE
float mchirp
Definition: netcluster.hh:84
ofstream out
Definition: cwb_merge.C:214
double gINRATE
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<< "\ layers : "<< nLAYERS<< "\ dF(hz) : "<< dF<< "\ 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:2207
std::vector< wavearray< double > * > RWFP
Definition: detector.hh:382
bool c4
Definition: Toolbox.hh:89
int _sse_mra_ps(float *, float *, float, int)
Definition: network.hh:1276
double fLow
Definition: config.hh:140
double gSEGGPS
#define CHI2_THR
int pattern
Definition: config.hh:241
netpixel pix(nifo)
netcluster * pwc
Definition: cwb_job_obj.C:38
std::vector< TGraph * > graph
Definition: watplot.hh:194
CWB_STAGE istage
Definition: cwb.hh:247
bool cedDump
Definition: config.hh:297
void ResetPCA(network *NET, CWB::config *cfg, netcluster *pwc, std::vector< netpixel > *vPIX, int ID)
dqfile DQF[DQF_MAX]
Definition: config.hh:349
std::vector< wavearray< double > > vINJ
wavearray< double > gsREC[NIFO_MAX]
waveform wf
Long_t size
WSeries< double > waveBand
Definition: detector.hh:356
int m
Definition: cwb_net.C:28
virtual void start(double s)
Definition: wavearray.hh:137
int j
Definition: cwb_net.C:28
bool getMRAwave(size_t ID, size_t lag, char atype='S', int mode=0, bool tof=false)
Definition: network.cc:3666
i drho i
double rate
Definition: netcluster.hh:378
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:150
std::vector< detector * > ifoList
Definition: network.hh:608
#define N
bool outPlugin
Definition: config.hh:369
bool core
Definition: netpixel.hh:120
char ifo[NIFO_MAX][8]
void dopen(const char *fname, char *mode, bool header=true)
Definition: netevent.hh:312
Definition: cwb2G.hh:33
TTree * gTREE
size_t ifoListSize()
Definition: network.hh:431
std::vector< wavearray< double > > GetInjWaveform(network *NET, netevent *EVT, int id, double factor)
int Psave
Definition: config.hh:193
wavearray< double > AddWaveforms(wavearray< double > *wf1, wavearray< double > *wf2)
#define nIFO
CWB_CAT cat
Definition: Toolbox.hh:86
wavearray< float > pNRG
buffers for cluster residual energy
Definition: network.hh:655
virtual size_t size() const
Definition: wavearray.hh:145
void PrintUserOptions()
double GetFrequencyBoundaries(wavearray< double > x, double P, double &bF, double &eF)
TCanvas * canvas
Definition: watplot.hh:192
wavearray< double > gwINJ[NIFO_MAX]
char file[1024]
Definition: Toolbox.hh:85
cwb gCWB
jfile
Definition: cwb_job_obj.C:43
void output2G(TTree *, network *, size_t, int, double)
Definition: netevent.cc:711
wavearray< double > gwDAT[NIFO_MAX]
TString gOUTPUT
std::vector< netpixel > DoPCA(network *NET, CWB::config *cfg, int lag, int id)
void ClearWaveforms(detector *ifo)
Int_t Psave
number of detectors
Definition: netevent.hh:66
int simulation
Definition: config.hh:199
int nDQF
Definition: config.hh:348
#define IMPORT(TYPE, VAR)
Definition: cwb.hh:69
xtalk getXTalk(int nLay1, size_t indx1, int nLay2, size_t indx2)
param: numbers of layers identifying the resolution of the first pixel param: TF map index of the fir...
Definition: monster.cc:351
Double_t * gps
average center_of_gravity time
Definition: netevent.hh:95
i() int(T_cor *100))
void DumpOutputFile(network *NET, netevent *&EVT, CWB::config *cfg, int ID, int k, int factor)
network NET
Definition: cwb_dump_inj.C:30
const int NIFO_MAX
Definition: wat.hh:22
int BATCH
Definition: config.hh:245
char tmp_dir[1024]
Definition: config.hh:325
char parPlugin[1024]
Definition: config.hh:363
wavearray< float > a_00
wdm multi-resolution analysis
Definition: network.hh:652
void SetEventWindow(CWB::config *cfg, double gps)
double start
Definition: netcluster.hh:379
printf("total live time: non-zero lags = %10.1f \, liveTot)
TIter next(twave->GetListOfBranches())
char fname[1024]
wavearray< double > gsINJ[NIFO_MAX]
double shift
Definition: Toolbox.hh:87
std::vector< int > RWFID
Definition: detector.hh:381
size_t mdc__IDSize()
Definition: network.hh:414
int gRATEANA
int k
uoptions gOPT
Float_t * netcc
effective correlated SNR
Definition: netevent.hh:113
_EBBH gEBBH
std::vector< wavearray< double > * > IWFP
Definition: detector.hh:379
double mchirp(int ID, double=2.5, double=1.e20, double=0)
Definition: netcluster.cc:1424
TObjArray * token
float chi2chirp
Definition: netcluster.hh:86
double acor
Definition: network.hh:585
double F
std::vector< int > IWFID
Definition: detector.hh:378
wavearray< double > gHOT[NIFO_MAX]
double GetCentralFrequency(wavearray< double > x)
#define TMERGER_CUT
std::vector< netpixel > GetCluster(network *NET, CWB::config *cfg, int lag, int id)
wavearray< double > sINJ[NIFO_MAX]
double GetCentralTime(wavearray< double > x)
TFile * froot
netpixel * getPixel(size_t n, size_t i)
Definition: netcluster.hh:413
virtual void FFTW(int=1)
Definition: wavearray.cc:896
#define NCHIRP_MAX
void ReadUserOptions(CWB::config *cfg)
char options[256]
#define DUMP_CLUSTER
float mchirperr
Definition: netcluster.hh:85
void GetEccentricityIndex(network *NET, int lag, int id)
char title[256]
Definition: SSeriesExample.C:1
Float_t * phi
sqrt(h+*h+ + hx*hx)
Definition: netevent.hh:87
double gps
netevent EVT(itree, nifo)
Float_t * theta
[0]-reconstructed, [1]-injected phi angle, [2]-RA
Definition: netevent.hh:88
Double_t * time
beam pattern coefficients for hx
Definition: netevent.hh:94
double T
Definition: testWDM_4.C:11
std::vector< int > sCuts
Definition: netcluster.hh:392
WSeries< double > waveForm
Definition: detector.hh:355
void dclose()
Definition: netevent.hh:321
#define DUMP_WAVE
double fHigh
Definition: config.hh:141
char Name[16]
Definition: detector.hh:327
char ifo[NIFO_MAX][8]
Definition: config.hh:124
double fabs(const Complex &x)
Definition: numpy.cc:55
wavearray< double > GetWaveform(int ifoId, network *NET, int lag, int id, char type, bool shift=true)
int gIFACTOR
strcpy(RunLabel, RUN_LABEL)
Float_t likelihood
Definition: netevent.hh:124
netcluster * getwc(size_t n)
param: delay index
Definition: network.hh:439
sprintf(tfres,"(1/%g)x(%g) (sec)x(Hz)", 2 *df, df)
std::vector< clusterdata > cData
Definition: netcluster.hh:391
int nIFO
Definition: config.hh:120
Definition: cwb.hh:136
DataType_t * data
Definition: wavearray.hh:319
wavearray< float > rNRG
buffer for cluster sky 90 amplitudes
Definition: network.hh:654
Long_t id
double factors[FACTORS_MAX]
Definition: config.hh:202
wavearray< double > wINJ[NIFO_MAX]
void SetOutputFile(network *NET, netevent *&EVT, CWB::config *cfg, bool dump_ebbh)
char data_label[1024]
Definition: config.hh:332
Bool_t fill_in(network *, int, bool=true)
Definition: injection.cc:349
WSeries< double > waveNull
Definition: detector.hh:357
size_t getmdc__ID(size_t n)
Definition: network.hh:426
double shift[NIFO_MAX]
wavearray< double > ** pwf
[x1,y1,z1,x2,y2,z2] components of spin vector
Definition: injection.hh:101
int check
size_t inRate
Definition: config.hh:132
detector ** pD
exit(0)
bool setdata(double a, char type='R', size_t n=0)
Definition: netpixel.hh:58
double avr
bool dump
Definition: config.hh:295