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