Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
CWB_Plugin_Qveto.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);
29  CWB::config* cfg, bool fft=false, bool strain=false);
30 
31 void
33 //!MISCELLANEA
34 // Extract whitened reconstructed waveforms, and compute the Qveto & Lveto parameters
35 
36  cout << endl;
37  cout << "-----> CWB_Plugin_Qveto.C" << endl;
38  cout << "ifo " << ifo.Data() << endl;
39  cout << "type " << type << endl;
40  cout << endl;
41 
42  float Qveto[NIFO_MAX]; // Qveto
43  float Lveto[3]; // Lveto
44 
45  if(type==CWB_PLUGIN_CONFIG) {
46  cfg->outPlugin=true; // disable built-in output root file
47  }
48 
49  if(type==CWB_PLUGIN_ILIKELIHOOD) {
50  NET->wfsave=true; // enable waveform save
51 
52  // search output root file in the system list
53  TFile* froot = NULL;
54  TList *files = (TList*)gROOT->GetListOfFiles();
55  TString outDump="";
56  netevent* EVT;
57  int nIFO = NET->ifoListSize(); // number of detectors
58  if (files) {
59  TIter next(files);
60  TSystemFile *file;
61  TString fname;
62  bool check=false;
63  while ((file=(TSystemFile*)next())) {
64  fname = file->GetName();
65  // set output root file as the current file
66  if(fname.Contains("wave_")) {
67  froot=(TFile*)file;froot->cd();
68  outDump=fname;
69  outDump.ReplaceAll(".root.tmp",".txt");
70  //cout << "output file name : " << fname << endl;
71  }
72  }
73  if(!froot) {
74  cout << "CWB_Plugin_Qveto.C : Error - output root file not found" << endl;
75  gSystem->Exit(1);
76  }
77  } else {
78  cout << "CWB_Plugin_Qveto.C : Error - output root file not found" << endl;
79  gSystem->Exit(1);
80  }
81 
82  TTree* net_tree = (TTree *) froot->Get("waveburst");
83  if(net_tree==NULL) {
84  EVT = new netevent(nIFO);
85  net_tree = EVT->setTree();
86  net_tree->Branch("Qveto",Qveto,TString::Format("Qveto[%i]/F",cfg->nIFO));
87  net_tree->Branch("Lveto",Lveto,TString::Format("Lveto[%i]/F",3));
88  }
89  }
90 
91  if(type==CWB_PLUGIN_OLIKELIHOOD) {
92 
93  if(TString(cfg->analysis)!="2G") {
94  cout << "CWB_Plugin_Qveto.C -> "
95  << "CWB_PLUGIN_OLIKELIHOOD implemented only for 2G" << endl;
96  gSystem->Exit(1);
97  }
98 
99  // import ifactor
100  int gIFACTOR=-1; IMPORT(int,gIFACTOR)
101  cout << "-----> CWB_Plugin_Qveto.C -> "
102  << " gIFACTOR : " << gIFACTOR << endl;
103 
104  // import slagShift
105  float* gSLAGSHIFT=NULL; IMPORT(float*,gSLAGSHIFT)
106 
107  int nIFO = NET->ifoListSize(); // number of detectors
108  int K = NET->nLag; // number of time lag
109  netevent* EVT;
111  //double factor = cfg->simulation==3||cfg->simulation==4 ? -gIFACTOR : cfg->factors[gIFACTOR];
112  double factor = cfg->factors[gIFACTOR];
113  int rate = 0; // select all resolutions
114 
115  // search output root file in the system list
116  TFile* froot = NULL;
117  TList *files = (TList*)gROOT->GetListOfFiles();
118  TString outDump="";
119  if (files) {
120  TIter next(files);
121  TSystemFile *file;
122  TString fname;
123  bool check=false;
124  while ((file=(TSystemFile*)next())) {
125  fname = file->GetName();
126  // set output root file as the current file
127  if(fname.Contains("wave_")) {
128  froot=(TFile*)file;froot->cd();
129  outDump=fname;
130  outDump.ReplaceAll(".root.tmp",".txt");
131  //cout << "output file name : " << fname << endl;
132  }
133  }
134  if(!froot) {
135  cout << "CWB_Plugin_Qveto.C : Error - output root file not found" << endl;
136  gSystem->Exit(1);
137  }
138  } else {
139  cout << "CWB_Plugin_Qveto.C : Error - output root file not found" << endl;
140  gSystem->Exit(1);
141  }
142 
143  TTree* net_tree = (TTree *) froot->Get("waveburst");
144  if(net_tree!=NULL) {
145  EVT = new netevent(net_tree,nIFO);
146  net_tree->SetBranchAddress("Qveto",Qveto);
147  net_tree->SetBranchAddress("Lveto",Lveto);
148  } else {
149  EVT = new netevent(nIFO);
150  net_tree = EVT->setTree();
151  net_tree->Branch("Qveto",Qveto,TString::Format("Qveto[%i]/F",cfg->nIFO));
152  net_tree->Branch("Lveto",Lveto,TString::Format("Lveto[%i]/F",3));
153  }
154  EVT->setSLags(gSLAGSHIFT); // set slags into netevent
155 
156  for(int k=0; k<K; k++) { // loop over the lags
157 
158  id = NET->getwc(k)->get(const_cast<char*>("ID"), 0, 'L', rate);
159 
160  for(int j=0; j<(int)id.size(); j++) { // loop over cluster index
161 
162  int ID = size_t(id.data[j]+0.5);
163 
164  if(NET->getwc(k)->sCuts[ID-1]!=-1) continue; // skip rejected/processed clusters
165 
166  double ofactor=0;
167  if(cfg->simulation==4) ofactor=-factor;
168  else if(cfg->simulation==3) ofactor=-gIFACTOR;
169  else ofactor=factor;
170 
171  EVT->output2G(NULL,NET,ID,k,ofactor); // get reconstructed parameters
172 
173  wavearray<double>** pwfREC = new wavearray<double>*[nIFO];
174  detector* pd = NET->getifo(0);
175  int idSize = pd->RWFID.size();
176 
177  int wfIndex=-1;
178  for (int mm=0; mm<idSize; mm++) if (pd->RWFID[mm]==ID) wfIndex=mm;
179  if(wfIndex==-1) continue;
180 
181  netcluster* pwc = NET->getwc(k);
182  cout << endl << "----------------------------------------------------------------" << endl;
183  GetLveto(pwc, ID, nIFO, Lveto);
184  cout << "Lveto : " << "fmean : " << Lveto[0] << " frms : " << Lveto[1]
185  << " Energy Ratio : " << Lveto[2] << endl << endl;
186 
187  // extract whitened reconstructed waveforms
188  for(int n=0; n<nIFO; n++) {
189 
190  pd = NET->getifo(n);
191 
192  pwfREC[n] = pd->RWFP[wfIndex];
193  wavearray<double>* wfREC = pwfREC[n]; // array of reconstructed waveforms
194 
195 #ifdef PLOT_WHITENED_WAVEFORMS
196  //PlotWaveform(NET->ifoName[n], wfREC, cfg, false, false);
197  PlotWaveform(NET->ifoName[n], wfREC, cfg, true, false);
198 #endif
199 
200  Qveto[n] = GetQveto(wfREC);
201  cout << "Qveto : " << pd->Name << " Qveto = " << Qveto[n] << endl;
202  }
203  cout << "----------------------------------------------------------------" << endl;
204  delete [] pwfREC;
205 
206  std::vector<int> sCuts = NET->getwc(k)->sCuts; // save cCuts
207  // set sCuts=1 to the events which must be not copied with cps to pwc
208  for(int i=0; i<(int)sCuts.size(); i++) if(i!=ID-1) NET->getwc(k)->sCuts[i]=1;
209 
210  // ID can not be used to get the event, to get event use ID=1 (ex: for watplot)
211  NET->getwc(k)->sCuts = sCuts; // restore cCuts
212 
213  if(cfg->dump) EVT->dopen(outDump.Data(),const_cast<char*>("a"),false);
214  EVT->output2G(net_tree,NET,ID,k,ofactor); // get reconstructed parameters
215  if(cfg->dump) {
216  // add Qveto to dump file
217  fprintf(EVT->fP,"Qveto: ");
218  for(int i=0; i<nIFO; i++) fprintf(EVT->fP,"%f ",Qveto[i]);
219  fprintf(EVT->fP,"\n");
220  // add Lveto to dump file
221  fprintf(EVT->fP,"Lveto: ");
222  for(int i=0; i<3; i++) fprintf(EVT->fP,"%f ",Lveto[i]);
223  fprintf(EVT->fP,"\n");
224  }
225  if(cfg->dump) EVT->dclose();
226  if(!cfg->cedDump) NET->getwc(k)->sCuts[ID-1]=1; // mark as processed
227  }
228  }
229 
230  jfile->cd();
231  if(EVT) delete EVT;
232  }
233  return;
234 }
235 
236 float
238 
239  wavearray<double> x = *wf;;
240 
241  // resample data by a factor 4
242  int xsize=x.size();
243  x.FFTW(1);
244  x.resize(4*x.size());
245  x.rate(4*x.rate());
246  for(int j=xsize;j<x.size();j++) x[j]=0;
247  x.FFTW(-1);
248 
249  // extract max/min values and save the absolute values in the array 'a'
250  wavearray<double> a(x.size());
251  int size=0;
252  double dt = 1./x.rate();
253  double prev=x[0];
254  double xmax=0;
255  for (int i=1;i<x.size();i++) {
256  if(fabs(x[i])>xmax) xmax=fabs(x[i]);
257  if(prev*x[i]<0) {
258  a[size]=xmax;
259  size++;
260  xmax=0;
261  }
262  prev=x[i];
263  }
264 
265  // find max value/index ans save on amax/imax
266  int imax=-1;
267  double amax=0;
268  for (int i=1;i<size;i++) {
269  if(a[i]>amax) {amax=a[i];imax=i;}
270  }
271 
272 /*
273  cout << endl;
274  cout << "a[imax-2] " << a[imax-2] << endl;
275  cout << "a[imax-1] " << a[imax-1] << endl;
276  cout << "a[imax] " << a[imax] << endl;
277  cout << "a[imax+1] " << a[imax+1] << endl;
278  cout << "a[imax+2] " << a[imax+2] << endl;
279  cout << endl;
280 */
281 
282  // compute Qveto
283  double ein=0; // energy of max values inside NTHR
284  double eout=0; // energy of max values outside NTHR
285  for (int i=0;i<size;i++) {
286  if(abs(imax-i)<=NTHR) {
287  ein+=a[i]*a[i];
288  //cout << i << " ein " << a[i] << " " << amax << endl;
289  } else {
290  if(a[i]>amax/ATHR) eout+=a[i]*a[i];
291  //if(a[i]>amax/ATHR) cout << i << " eout " << a[i] << " " << amax << endl;
292  }
293  }
294  float Qveto = ein>0 ? eout/ein : 0.;
295  //cout << "Qveto : " << Qveto << " ein : " << ein << " eout : " << eout << endl;
296 
297  return Qveto;
298 }
299 
300 void
301 GetLveto(netcluster* pwc, int cid, int nifo, float* Lveto) {
302 //
303 // input
304 // pwc : pointer to netcluster object
305 // cid : cluster id
306 // nifo : number of detectors
307 // output
308 // Lveto[0] : line frequency
309 // Lveto[1] : line bandwitdh
310 // Lveto[2] : line enery ratio (line_energy / total_energy)
311 //
312 
313  Lveto[0] = Lveto[1] = Lveto[2] = 0;
314 
315  std::vector<int>* vint = &(pwc->cList[cid-1]); // pixel list
316  int V = vint->size(); // cluster size
317  if(!V) return;
318 
319  // ------------------------------------------------------------------
320  // Find max pixel parameters
321  // ------------------------------------------------------------------
322 
323  double likeMax=0; // maximum pixel's energy
324  double likeTot=0; // total cluster energy
325  double freqMax; // frequency of the pixel with max energy
326  double dfreqMax; // df of the pixel with max energy
327  for(int n=0; n<V; n++) {
328  netpixel* pix = pwc->getPixel(cid,n);
329  if(pix->layers%2==0) {
330  cout << "CWB_Plugin_Qveto.C - Error : is enabled only for WDM (2G)" << endl;
331  exit(1);
332  }
333  if(!pix->core) continue; // select only the principal components pixels
334 
335  double likePix=0;
336  for(int m=0; m<nifo; m++) {
337  likePix += pow(pix->getdata('S',m),2); // like whitened reconstructed signal 00
338  likePix += pow(pix->getdata('P',m),2); // like whitened reconstructed signal 90
339  }
340 
341  double freq = pix->frequency*pix->rate/2.;
342  double df = pix->rate/2.;
343 
344  likeTot+=likePix;
345  if(likePix>likeMax) {likeMax=likePix;freqMax=freq;dfreqMax=df;}
346  }
347  //cout << "likeMax : " << likeMax << " likeTot : " << likeTot
348  // << " freqMax : " << freqMax << " dfreqMax : " << dfreqMax << endl;
349 
350  // ------------------------------------------------------------------
351  // Compute Lveto parameters
352  // ------------------------------------------------------------------
353 
354  double fmean=0; // line mean frequency
355  double frms=0; // line bandwidth
356  double likeLin=0; // line energy
357  for(int n=0; n<V; n++) {
358  netpixel* pix = pwc->getPixel(cid,n);
359  if(!pix->core) continue; // select only the principal components pixels
360 
361  double likePix=0;
362  for(int m=0; m<nifo; m++) {
363  likePix += pow(pix->getdata('S',m),2); // like whitened reconstructed signal 00
364  likePix += pow(pix->getdata('P',m),2); // like whitened reconstructed signal 90
365  }
366 
367  // the estimation is done for all pixels
368  // with freq in the range [freqMax-dfreqMax, freqMax+dfreqMax]
369  double freq = pix->frequency*pix->rate/2.;
370  if(fabs(freq-freqMax)<=dfreqMax) {
371  likeLin += likePix;
372  fmean += likePix*freq;
373  frms += likePix*freq*freq;
374  }
375  }
376 
377  fmean = fmean/likeLin;
378  frms = frms/likeLin-fmean*fmean;
379  frms = frms>0 ? sqrt(frms) : 0.;
380 
381  if(frms<dfreqMax/2.) frms=dfreqMax/2.;
382 
383  // ------------------------------------------------------------------
384  // Save Lveto parameters
385  // ------------------------------------------------------------------
386 
387  Lveto[0] = fmean; // line mean frequency
388  Lveto[1] = frms; // line bandwidth
389  Lveto[2] = likeTot>0. ? likeLin/likeTot : 0.; // energy ratio energy inside_line/total
390 
391  // ------------------------------------------------------------------
392  // plot time-frequency energy
393  // ------------------------------------------------------------------
394 
395 #if defined PLOT_LIKELIHOOD
396  watplot WTS(const_cast<char*>("wts"));
397  WTS.plot(pwc, cid, nifo, 'L', 0, const_cast<char*>("COLZ"));
398  WTS.canvas->Print("l_tfmap_scalogram.png");
399 #endif
400 
401  return;
402 }
403 
404 void
406  CWB::config* cfg, bool fft, bool strain) {
407 
408  watplot PTS(const_cast<char*>("ptswrc"),200,20,800,500);
409 
410  //cout << "Print " << fname << endl;
411  double tmin = wfREC->start();
412  wfREC->start(wfREC->start()-tmin);
413  if(fft) {
414  PTS.plot(wfREC, const_cast<char*>("ALP"), 1, 0, 0, true, cfg->fLow, cfg->fHigh);
415  } else {
416  PTS.plot(wfREC, const_cast<char*>("ALP"), 1, 0, 0);
417  }
418  PTS.graph[0]->SetLineWidth(1);
419  wfREC->start(wfREC->start()+tmin);
420 
421  char label[64]="";
422  if(fft) sprintf(label,"%s","fft");
423  else sprintf(label,"%s","time");
424  if(strain) sprintf(label,"%s_%s",label,"strain");
425  else sprintf(label,"%s_%s",label,"white");
426 
427  char fname[1024];
428  //sprintf(fname, "%s_wf_%s_rec_gps_%d.root",ifo.Data(),label,int(tmin));
429  sprintf(fname, "%s_wf_%s_rec_gps_%d.png",ifo.Data(),label,int(tmin));
430  PTS.canvas->Print(fname);
431  cout << "write : " << fname << endl;
432  //PTS.canvas->Write(REPLACE(fname,dirCED,gtype));
433 }
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
#define ATHR
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
wavearray< double > a(hp.size())
size_t frequency
Definition: netpixel.hh:93
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")
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
void PlotWaveform(TString ifo, wavearray< double > *wfREC, CWB::config *cfg, bool fft=false, bool strain=false)
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
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
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
void CWB_Plugin(TFile *jfile, CWB::config *cfg, network *NET, WSeries< double > *x, TString ifo, int type)
COHERENCE.
#define nIFO
TCanvas * canvas
Definition: watplot.hh:174
double factor
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
#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
TIter next(twave->GetListOfBranches())
char fname[1024]
std::vector< int > RWFID
Definition: detector.hh:363
int k
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
void dclose()
Definition: netevent.hh:300
double fHigh
Definition: config.hh:123
char Name[16]
Definition: detector.hh:309
double fabs(const Complex &x)
Definition: numpy.cc:37
#define NTHR
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)
float GetQveto(wavearray< double > *wf)
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
void GetLveto(netcluster *pwc, int cid, int nifo, float *Lveto)
virtual void resize(unsigned int)
Definition: wavearray.cc:445
int check
exit(0)
bool dump
Definition: config.hh:277