Logo coherent WaveBurst  
Library Reference Guide
Logo
CWB_Plugin_QLveto.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 
44 double GetQfactor(wavearray<double>* wf, double frequency, bool fixAmax);
46 void GetGaussFitPars(wavearray<double>* wf, double& mean, double& sigma, bool doenv=true);
47 void GetGaussFitPars2(wavearray<double>* wf, double& mean, double& sigma, bool fixAmax);
48 
49 void GetQveto(wavearray<double>* wf, float &Qveto, float &Qfactor);
50 void GetLveto(netcluster* pwc, int cid, int nifo, float* Lveto);
52  CWB::config* cfg, bool fft=false, bool strain=false);
54 
55 std::vector<netpixel> DoPCA(network* NET, CWB::config* cfg, int lag, int id);
56 void ResetPCA(network* NET, CWB::config* cfg, netcluster* pwc, std::vector<netpixel>* vPIX, int ID);
57 
58 
59 void
61 //!MISCELLANEA
62 // Extract whitened reconstructed waveforms, and compute the Qveto, Lveto parameters
63 
64 // cout << endl;
65 // cout << "-----> CWB_Plugin_QLveto.C" << endl;
66 // cout << "ifo " << ifo.Data() << endl;
67 // cout << "type " << type << endl;
68 // cout << endl;
69 
70  float Qveto[4]; // Qveto
71  float Lveto[3]; // Lveto
72 
73  if(type==CWB_PLUGIN_CONFIG) {
74  cout << endl;
75  cout << "-----> CWB_Plugin_QLveto.C" << endl;
76  cout << endl;
77  cfg->outPlugin=true; // disable built-in output root file
78  }
79 
80  if(type==CWB_PLUGIN_ILIKELIHOOD) {
81  NET->wfsave=true; // enable waveform save
82 
83  // search output root file in the system list
84  TFile* froot = NULL;
85  TList *files = (TList*)gROOT->GetListOfFiles();
86  TString outDump="";
87  netevent* EVT;
88  int nIFO = NET->ifoListSize(); // number of detectors
89  if (files) {
90  TIter next(files);
91  TSystemFile *file;
92  TString fname;
93  bool check=false;
94  while ((file=(TSystemFile*)next())) {
95  fname = file->GetName();
96  // set output root file as the current file
97  if(fname.Contains("wave_")) {
98  froot=(TFile*)file;froot->cd();
99  outDump=fname;
100  outDump.ReplaceAll(".root.tmp",".txt");
101  //cout << "output file name : " << fname << endl;
102  }
103  }
104  if(!froot) {
105  cout << "CWB_Plugin_QLveto.C : Error - output root file not found" << endl;
106  gSystem->Exit(1);
107  }
108  } else {
109  cout << "CWB_Plugin_QLveto.C : Error - output root file not found" << endl;
110  gSystem->Exit(1);
111  }
112 
113  TTree* net_tree = (TTree *) froot->Get("waveburst");
114  if(net_tree==NULL) {
115  EVT = new netevent(nIFO);
116  net_tree = EVT->setTree();
117  net_tree->Branch("Qveto",Qveto,TString::Format("Qveto[%i]/F",4));
118  net_tree->Branch("Lveto",Lveto,TString::Format("Lveto[%i]/F",3));
119  } else {
120  TBranch* branch;
121  bool qveto_exists=false;
122  bool lveto_exists=false;
123  TIter next(net_tree->GetListOfBranches());
124  while ((branch=(TBranch*)next())) {
125  if(TString("Qveto").CompareTo(branch->GetName())==0) qveto_exists=true;
126  if(TString("Lveto").CompareTo(branch->GetName())==0) lveto_exists=true;
127  }
128  next.Reset();
129  if(!qveto_exists) net_tree->Branch("Qveto",Qveto,TString::Format("Qveto[%i]/F",4));
130  if(!lveto_exists) net_tree->Branch("Lveto",Lveto,TString::Format("Lveto[%i]/F",3));
131  }
132  }
133 
134  if(type==CWB_PLUGIN_OLIKELIHOOD) {
135 
136  if(TString(cfg->analysis)!="2G") {
137  cout << "CWB_Plugin_QLveto.C -> "
138  << "CWB_PLUGIN_OLIKELIHOOD implemented only for 2G" << endl;
139  gSystem->Exit(1);
140  }
141 
142  // import ifactor
143  int gIFACTOR=-1; IMPORT(int,gIFACTOR)
144  //cout << "-----> CWB_Plugin_QLveto.C -> "
145  // << " gIFACTOR : " << gIFACTOR << endl;
146 
147  // import slagShift
148  float* gSLAGSHIFT=NULL; IMPORT(float*,gSLAGSHIFT)
149 
150  int nIFO = NET->ifoListSize(); // number of detectors
151  int K = NET->nLag; // number of time lag
152  netevent* EVT;
154  //double factor = cfg->simulation==3||cfg->simulation==4 ? -gIFACTOR : cfg->factors[gIFACTOR];
155  double factor = cfg->factors[gIFACTOR];
156  int rate = 0; // select all resolutions
157 
158  // search output root file in the system list
159  TFile* froot = NULL;
160  TList *files = (TList*)gROOT->GetListOfFiles();
161  TString outDump="";
162  if (files) {
163  TIter next(files);
164  TSystemFile *file;
165  TString fname;
166  bool check=false;
167  while ((file=(TSystemFile*)next())) {
168  fname = file->GetName();
169  // set output root file as the current file
170  if(fname.Contains("wave_")) {
171  froot=(TFile*)file;froot->cd();
172  outDump=fname;
173  outDump.ReplaceAll(".root.tmp",".txt");
174  //cout << "output file name : " << fname << endl;
175  }
176  }
177  if(!froot) {
178  cout << "CWB_Plugin_QLveto.C : Error - output root file not found" << endl;
179  gSystem->Exit(1);
180  }
181  } else {
182  cout << "CWB_Plugin_QLveto.C : Error - output root file not found" << endl;
183  gSystem->Exit(1);
184  }
185 
186  TTree* net_tree = (TTree *) froot->Get("waveburst");
187  if(net_tree!=NULL) {
188  EVT = new netevent(net_tree,nIFO);
189  net_tree->SetBranchAddress("Qveto",Qveto);
190  net_tree->SetBranchAddress("Lveto",Lveto);
191  } else {
192  EVT = new netevent(nIFO);
193  net_tree = EVT->setTree();
194  net_tree->Branch("Qveto",Qveto,TString::Format("Qveto[%i]/F",4));
195  net_tree->Branch("Lveto",Lveto,TString::Format("Lveto[%i]/F",3));
196  }
197  EVT->setSLags(gSLAGSHIFT); // set slags into netevent
198 
199  for(int k=0; k<K; k++) { // loop over the lags
200 
201  id = NET->getwc(k)->get(const_cast<char*>("ID"), 0, 'L', rate);
202 
203  for(int j=0; j<(int)id.size(); j++) { // loop over cluster index
204 
205  int ID = size_t(id.data[j]+0.5);
206 
207  if(NET->getwc(k)->sCuts[ID-1]!=-1) continue; // skip rejected/processed clusters
208 
209  double ofactor=0;
210  if(cfg->simulation==4) ofactor=-factor;
211  else if(cfg->simulation==3) ofactor=-gIFACTOR;
212  else ofactor=factor;
213 
214  EVT->output2G(NULL,NET,ID,k,ofactor); // get reconstructed parameters
215 
216  wavearray<double>** pwfREC = new wavearray<double>*[nIFO];
217  detector* pd = NET->getifo(0);
218  int idSize = pd->RWFID.size();
219 
220  int wfIndex=-1;
221  for (int mm=0; mm<idSize; mm++) if (pd->RWFID[mm]==ID) wfIndex=mm;
222  if(wfIndex==-1) continue;
223 
224  netcluster* pwc = NET->getwc(k);
225  cout << endl << "----------------------------------------------------------------" << endl;
226 
227  // extract whitened reconstructed waveforms
228  Qveto[0]=Qveto[1]=1.e20;
229  double Qfactor1[NIFO_MAX];
230  double Qfactor2[NIFO_MAX];
231  for(int n=0; n<nIFO; n++) {
232 
233  pd = NET->getifo(n);
234 
235  pwfREC[n] = pd->RWFP[wfIndex];
236  wavearray<double>* wfREC = pwfREC[n]; // array of reconstructed waveforms
237 
238 #ifdef PLOT_WHITENED_WAVEFORMS
239  //PlotWaveform(NET->ifoName[n], wfREC, cfg, false, false);
240  PlotWaveform(NET->ifoName[n], wfREC, cfg, true, false);
241 #endif
242  // store in Qveto[0/1] the minimum velue of Qveto/Qfactor between the reconstructed,whitened waveforms in all ifos
243  float qveto,qfactor;
244 
245  // reconstructed whitened waveform
246  NET->getMRAwave(ID,k,'S',0,true);
247  GetQveto(&(pd->waveForm), qveto, qfactor);
248  if(qveto<Qveto[0]) Qveto[0]=qveto;
249  if(qfactor<Qveto[1]) Qveto[1]=qfactor;
250 
251  Qfactor1[n] = qfactor;
252  double fpeak = GetPeakFrequency(&(pd->waveForm));
253  Qfactor2[n] = GetQfactor(&(pd->waveForm), fpeak, true); // get qfactor with amax fixed;
254  cout << endl;
255  cout << "Qfactor1/2 : " << pd->Name << " frequency[0] = " << EVT->frequency[0] << " fpeak = " << fpeak
256  << " -> " << Qfactor1[n] << " / " << Qfactor2[n] << endl;
257  cout << endl;
258 
259  // whitened waveform
260  NET->getMRAwave(ID,k,'W',0,true);
261  GetQveto(&(pd->waveBand), qveto, qfactor);
262  if(qveto<Qveto[0]) Qveto[0]=qveto;
263  if(qfactor<Qveto[1]) Qveto[1]=qfactor;
264 
265  cout << "Qveto : " << pd->Name << " Qveto[0] = " << Qveto[0]
266  << " Qveto[1] = " << Qveto[1] << endl;
267 
268  if(!cfg->simulation) ClearWaveforms(pd); // release waveform memory
269  }
270  delete [] pwfREC;
271 
272  // compute Qveto[2] using the weighted of Qfactor1 values
273  Qveto[2]=0.;
274  for(int n=0; n<nIFO; n++) Qveto[2]+=EVT->sSNR[n]*Qfactor1[n];
275  Qveto[2]/=EVT->likelihood;
276  // compute Qveto[3] using the weighted of Qfactor2 values
277  Qveto[3]=0.;
278  for(int n=0; n<nIFO; n++) Qveto[3]+=EVT->sSNR[n]*Qfactor2[n];
279  Qveto[3]/=EVT->likelihood;
280  cout << "Qveto :" << " Qveto[2] = " << Qveto[2]
281  << " Qveto[3] = " << Qveto[3] << endl;
282 
283  std::vector<netpixel> vPIX;
284  if(cfg->pattern>0) vPIX = DoPCA(NET, cfg, k, ID); // do PCA analysis
285  GetLveto(pwc, ID, nIFO, Lveto);
286  if(cfg->pattern>0) ResetPCA(NET, cfg, pwc, &vPIX, ID); // restore WP pwc
287 
288  cout << endl;
289  cout << "Lveto : " << "fmean : " << Lveto[0] << " frms : " << Lveto[1]
290  << " Energy Ratio : " << Lveto[2] << endl << endl;
291  cout << "----------------------------------------------------------------" << endl;
292 
293  std::vector<int> sCuts = NET->getwc(k)->sCuts; // save cCuts
294  // set sCuts=1 to the events which must be not copied with cps to pwc
295  for(int i=0; i<(int)sCuts.size(); i++) if(i!=ID-1) NET->getwc(k)->sCuts[i]=1;
296 
297  // ID can not be used to get the event, to get event use ID=1 (ex: for watplot)
298  NET->getwc(k)->sCuts = sCuts; // restore cCuts
299 
300  if(cfg->dump) EVT->dopen(outDump.Data(),const_cast<char*>("a"),false);
301  EVT->output2G(net_tree,NET,ID,k,ofactor); // get reconstructed parameters
302  if(cfg->dump) {
303  // add Qveto to dump file
304  fprintf(EVT->fP,"Qveto: ");
305  for(int i=0; i<2*nIFO; i++) fprintf(EVT->fP,"%f ",Qveto[i]);
306  fprintf(EVT->fP,"\n");
307  // add Lveto to dump file
308  fprintf(EVT->fP,"Lveto: ");
309  for(int i=0; i<3; i++) fprintf(EVT->fP,"%f ",Lveto[i]);
310  fprintf(EVT->fP,"\n");
311  }
312  if(cfg->dump) EVT->dclose();
313  if(!cfg->cedDump) NET->getwc(k)->sCuts[ID-1]=1; // mark as processed
314  }
315  }
316 
317  jfile->cd();
318  if(EVT) delete EVT;
319  }
320  return;
321 }
322 
323 void
324 GetQveto(wavearray<double>* wf, float &Qveto, float &Qfactor) {
325 
326  wavearray<double> x = *wf;;
327 
328  // resample data by a factor 4
329  int xsize=x.size();
330  x.FFTW(1);
331  x.resize(4*x.size());
332  x.rate(4*x.rate());
333  for(int j=xsize;j<x.size();j++) x[j]=0;
334  x.FFTW(-1);
335 
336  // extract max/min values and save the absolute values in the array 'a'
337  wavearray<double> a(x.size());
338  int size=0;
339  double dt = 1./x.rate();
340  double prev=x[0];
341  double xmax=0;
342  for (int i=1;i<x.size();i++) {
343  if(fabs(x[i])>xmax) xmax=fabs(x[i]);
344  if(prev*x[i]<0) {
345  a[size]=xmax;
346  size++;
347  xmax=0;
348  }
349  prev=x[i];
350  }
351 
352  // find max value/index ans save on amax/imax
353  int imax=-1;
354  double amax=0;
355  for (int i=1;i<size;i++) {
356  if(a[i]>amax) {amax=a[i];imax=i;}
357  }
358 
359 /*
360  cout << endl;
361  cout << "a[imax-2] " << a[imax-2] << endl;
362  cout << "a[imax-1] " << a[imax-1] << endl;
363  cout << "a[imax] " << a[imax] << endl;
364  cout << "a[imax+1] " << a[imax+1] << endl;
365  cout << "a[imax+2] " << a[imax+2] << endl;
366  cout << endl;
367 */
368 
369  // compute Qveto
370  double ein=0; // energy of max values inside NTHR
371  double eout=0; // energy of max values outside NTHR
372  for (int i=0;i<size;i++) {
373  if(abs(imax-i)<=NTHR) {
374  ein+=a[i]*a[i];
375  //cout << i << " ein " << a[i] << " " << amax << endl;
376  } else {
377  if(a[i]>amax/ATHR) eout+=a[i]*a[i];
378  //if(a[i]>amax/ATHR) cout << i << " eout " << a[i] << " " << amax << endl;
379  }
380  }
381  Qveto = ein>0 ? eout/ein : 0.;
382  //cout << "Qveto : " << Qveto << " ein : " << ein << " eout : " << eout << endl;
383 
384  // compute Qfactor
385  float R = (a[imax-1]+a[imax+1])/a[imax]/2.;
386  Qfactor = sqrt(-pow(TMath::Pi(),2)/log(R)/2.);
387  //cout << "Qfactor : " << Qfactor << endl;
388 
389  return;
390 }
391 
392 void
393 GetLveto(netcluster* pwc, int cid, int nifo, float* Lveto) {
394 //
395 // input
396 // pwc : pointer to netcluster object
397 // cid : cluster id
398 // nifo : number of detectors
399 // output
400 // Lveto[0] : line frequency
401 // Lveto[1] : line bandwitdh
402 // Lveto[2] : line enery ratio (line_energy / total_energy)
403 //
404 
405  Lveto[0] = Lveto[1] = Lveto[2] = 0;
406 
407  std::vector<int>* vint = &(pwc->cList[cid-1]); // pixel list
408  int V = vint->size(); // cluster size
409  if(!V) return;
410 
411  // ------------------------------------------------------------------
412  // Find max pixel parameters
413  // ------------------------------------------------------------------
414 
415  double likeMax=0; // maximum pixel's energy
416  double likeTot=0; // total cluster energy
417  double freqMax; // frequency of the pixel with max energy
418  double dfreqMax; // df of the pixel with max energy
419  for(int n=0; n<V; n++) {
420  netpixel* pix = pwc->getPixel(cid,n);
421  if(pix->layers%2==0) {
422  cout << "CWB_Plugin_QLveto.C - Error : is enabled only for WDM (2G)" << endl;
423  exit(1);
424  }
425  if(!pix->core) continue; // select only the principal components pixels
426 
427  double likePix=0;
428  for(int m=0; m<nifo; m++) {
429  likePix += pow(pix->getdata('S',m),2); // like whitened reconstructed signal 00
430  likePix += pow(pix->getdata('P',m),2); // like whitened reconstructed signal 90
431  }
432 
433  double freq = pix->frequency*pix->rate/2.;
434  double df = pix->rate/2.;
435 
436  likeTot+=likePix;
437  if(likePix>likeMax) {likeMax=likePix;freqMax=freq;dfreqMax=df;}
438  }
439  //cout << "likeMax : " << likeMax << " likeTot : " << likeTot
440  // << " freqMax : " << freqMax << " dfreqMax : " << dfreqMax << endl;
441 
442  // ------------------------------------------------------------------
443  // Compute Lveto parameters
444  // ------------------------------------------------------------------
445 
446  double fmean=0; // line mean frequency
447  double frms=0; // line bandwidth
448  double likeLin=0; // line energy
449  for(int n=0; n<V; n++) {
450  netpixel* pix = pwc->getPixel(cid,n);
451  if(!pix->core) continue; // select only the principal components pixels
452 
453  double likePix=0;
454  for(int m=0; m<nifo; m++) {
455  likePix += pow(pix->getdata('S',m),2); // like whitened reconstructed signal 00
456  likePix += pow(pix->getdata('P',m),2); // like whitened reconstructed signal 90
457  }
458 
459  // the estimation is done for all pixels
460  // with freq in the range [freqMax-dfreqMax, freqMax+dfreqMax]
461  double freq = pix->frequency*pix->rate/2.;
462  if(fabs(freq-freqMax)<=dfreqMax) {
463  likeLin += likePix;
464  fmean += likePix*freq;
465  frms += likePix*freq*freq;
466  }
467  }
468 
469  fmean = fmean/likeLin;
470  frms = frms/likeLin-fmean*fmean;
471  frms = frms>0 ? sqrt(frms) : 0.;
472 
473  if(frms<dfreqMax/2.) frms=dfreqMax/2.;
474 
475  // ------------------------------------------------------------------
476  // Save Lveto parameters
477  // ------------------------------------------------------------------
478 
479  Lveto[0] = fmean; // line mean frequency
480  Lveto[1] = frms; // line bandwidth
481  Lveto[2] = likeTot>0. ? likeLin/likeTot : 0.; // energy ratio energy inside_line/total
482 
483  // ------------------------------------------------------------------
484  // plot time-frequency energy
485  // ------------------------------------------------------------------
486 
487 #if defined PLOT_LIKELIHOOD
488  watplot WTS(const_cast<char*>("wts"));
489  WTS.plot(pwc, cid, nifo, 'L', 0, const_cast<char*>("COLZ"));
490  WTS.canvas->Print("l_tfmap_scalogram.png");
491 #endif
492 
493  return;
494 }
495 
496 void
498  CWB::config* cfg, bool fft, bool strain) {
499 
500  watplot PTS(const_cast<char*>("ptswrc"),200,20,800,500);
501 
502  //cout << "Print " << fname << endl;
503  double tmin = wfREC->start();
504  wfREC->start(wfREC->start()-tmin);
505  if(fft) {
506  PTS.plot(wfREC, const_cast<char*>("ALP"), 1, 0, 0, true, cfg->fLow, cfg->fHigh);
507  } else {
508  PTS.plot(wfREC, const_cast<char*>("ALP"), 1, 0, 0);
509  }
510  PTS.graph[0]->SetLineWidth(1);
511  wfREC->start(wfREC->start()+tmin);
512 
513  char label[64]="";
514  if(fft) sprintf(label,"%s","fft");
515  else sprintf(label,"%s","time");
516  if(strain) sprintf(label,"%s_%s",label,"strain");
517  else sprintf(label,"%s_%s",label,"white");
518 
519  char fname[1024];
520  //sprintf(fname, "%s_wf_%s_rec_gps_%d.root",ifo.Data(),label,int(tmin));
521  sprintf(fname, "%s_wf_%s_rec_gps_%d.png",ifo.Data(),label,int(tmin));
522  PTS.canvas->Print(fname);
523  cout << "write : " << fname << endl;
524  //PTS.canvas->Write(REPLACE(fname,dirCED,gtype));
525 }
526 
527 void
529 
530  int n;
531 
532  n = ifo->IWFP.size();
533  for (int i=0;i<n;i++) {
535  delete wf;
536  }
537  ifo->IWFP.clear();
538  ifo->IWFID.clear();
539 
540  n = ifo->RWFP.size();
541  for (int i=0;i<n;i++) {
543  delete wf;
544  }
545  ifo->RWFP.clear();
546  ifo->RWFID.clear();
547 }
548 
549 std::vector<netpixel> DoPCA(network* NET, CWB::config* cfg, int lag, int id) {
550 
551  double ee;
552 
553  size_t nIFO = NET->ifoList.size();
554 
555  float En = 2*NET->acor*NET->acor*nIFO; // network energy threshold in the sky loop
556 
557  int size = NET->a_00.size();
558  int f_ = NIFO/4;
559 
560  netpixel* pix;
561  netcluster* pwc = NET->getwc(lag);
562  std::vector<netpixel> vPIX;
563 
564  std::vector<int> pI = NET->wdmMRA.getXTalk(pwc, id, false); // buffer for pixel IDs
565  int V = pI.size();
566  if(V>cfg->BATCH) return vPIX; // attach TD amp to pixels < V
567 
568  wavearray<float> xi(size); xi=0; // PC 00 array
569  wavearray<float> XI(size); XI=0; // PC 90 array
570 
571  __m128* _xi = (__m128*) xi.data; // set pointer to PC 00 array
572  __m128* _XI = (__m128*) XI.data; // set pointer to PC 90 array
573 
574  __m128* _aa = (__m128*) NET->a_00.data; // set pointer to 00 array
575  __m128* _AA = (__m128*) NET->a_90.data; // set pointer to 90 array
576 
577  int nPC = 0;
578  for(int j=0; j<V; j++) {
579  int jf = j*f_; // source sse pointer increment
580  _sse_zero_ps(_xi+jf); // zero MRA amplitudes
581  _sse_zero_ps(_XI+jf); // zero MRA amplitudes
582  ee = _sse_abs_ps(_aa+jf,_AA+jf); // total pixel energy / quadrature
583  if(ee>En) nPC++; else ee=0.; // count core pixels
584  NET->rNRG.data[j] = ee; // init residual energy array
585  NET->pNRG.data[j] = NET->rNRG.data[j]; // init residual energy array
586  }
587 
588  nPC = NET->_sse_mra_ps(xi.data,XI.data,En,nPC); // get principal components
589 
590  for(int j=0; j<V; j++) { // loop over principal components
591  pix = pwc->getPixel(id,pI[j]);
592  vPIX.push_back(*pix); // save original pixels
593  pix->core = false;
594  ee = NET->pNRG.data[j]; // total pixel energy
595  if(ee<En) continue;
596  pix->core = true;
597  for(int i=0; i<nIFO; i++) {
598  pix->setdata(double(xi.data[j*NIFO+i]),'S',i); // store 00 whitened response PC
599  pix->setdata(double(XI.data[j*NIFO+i]),'P',i); // store 90 whitened response PC
600  }
601  }
602 
603  return vPIX;
604 }
605 
606 void ResetPCA(network* NET, CWB::config* cfg, netcluster* pwc, std::vector<netpixel>* vPIX, int ID) {
607 
608  std::vector<int> pI = NET->wdmMRA.getXTalk(pwc, ID, false); // buffer for pixel IDs
609  int V = pI.size();
610  if(V>cfg->BATCH) return; // attach TD amp to pixels < V
611  for(int j=0; j<V; j++) {
612  netpixel* pix = pwc->getPixel(ID,pI[j]);
613  *pix = (*vPIX)[j];
614  }
615 
616  while(!vPIX->empty()) vPIX->pop_back();
617  vPIX->clear(); std::vector<netpixel>().swap(*vPIX);
618 }
619 
620 double
621 GetQfactor(wavearray<double>* wf, double frequency, bool fixAmax) {
622 //
623 // wf -> input waveform
624 // frequency -> input peak frequency
625 // fixAmax -> true/false : fix/not-fix amplitude of the gaussian fit
626 //
627 // -> the input waveform wf is resampled 4x
628 // -> get mean and sigma from a gaussian fit of the waveform envelope
629 // -> compute the range [xmin,xmax] = [mean-3*sigma,mean-3*sigma]]
630 // -> compute area of the envelope in the rage [xmin,xmax]
631 // -> compute qfactor using area, max_amplitude and frequency
632 // -> return qfactor
633 //
634 
635  wavearray<double> x = *wf;;
636 
637  // resample data by a factor 4
638  int xsize=x.size();
639 
640  x.FFTW(1);
641  x.resize(4*x.size());
642  x.rate(4*x.rate());
643  for(int j=xsize;j<x.size();j++) x[j]=0;
644  x.FFTW(-1);
645 
646  // compute range [xmin,xmax] used for final fit using the square of x
647  double gmean, gsigma;
648  GetGaussFitPars2(&x, gmean, gsigma, fixAmax);
649 // cout << "GetQfactor -> gmean " << gmean << " gsigma " << gsigma << endl;
650 
651  double xmin = gmean-3.0*gsigma;
652  double xmax = gmean+3.0*gsigma;
653 // cout << "GetQfactor -> xmin " << xmin << " xmax " << xmax << endl;
654 
655  x = CWB::mdc::GetEnvelope(&x); // get envelope
656 
657  // extract max/min values and save the absolute values in the array 'a'
658  double dt = 1./x.rate();
659  // find max value/index and save on amax/imax
660  int imax=-1;
661  double amax=0;
662  for (int i=1;i<x.size();i++) {
663  if(x[i]>amax) {amax=x[i];imax=i;}
664  }
665 
666  double area=0.;
667  for (int i=1;i<x.size();i++) {
668  double t=i*dt;
669  if(t>xmin && t<xmax) area+=x[i]*dt;
670  }
671 
672  double sigma = area/amax/sqrt(2*TMath::Pi());
673 // cout << "sigma integral is " << sigma << endl;
674  double qfactor = sigma*(TMath::TwoPi()*frequency); // compute qfactor
675 
676  return qfactor;
677 }
678 
679 double
681 //
682 // wf -> input waveform
683 //
684 // -> multiply the input waveform wf by its envelope
685 // -> perform FFT transform
686 // -> compute max frequency -> peak
687 // -> return frequency peak
688 //
689 
691  wavearray<double> ex = CWB::mdc::GetEnvelope(&x); // get envelope
692  int N = x.size();
693  for(int i=0;i<N;i++) x[i]*=ex[i];
694 
695  x.FFTW(1);
696 
697  // fill time/amp graph
698  double df=(x.rate()/(double)N)/2.;
699  double amax=0.;
700  double fmax=0.;
701  for(int i=0;i<N;i+=2) {
702  double freq=i*df;
703  double amp=x[i]*x[i]+x[i+1]*x[i+1];
704  if(amp>amax) {fmax=freq;amax=amp;}
705  }
706 
707  return fmax;
708 }
709 
710 void
711 GetGaussFitPars(wavearray<double>* wf, double& mean, double& sigma, bool doenv) {
712 //
713 // wf -> input waveform
714 // doenv -> true/false do/not-do waveform envelope
715 //
716 // mean, sigma -> output gaussian fit parameters
717 //
718 // -> get evelope from the input waveform wf and perform fit with gaussian function
719 // -> return mean and gigma of gaussian fit
720 //
721 
722  wavearray<double> x = doenv ? CWB::mdc::GetEnvelope(wf) : *wf; // get envelope
723 
724  // fill time/amp graph
725  int N = x.size();
726  double* time = new double[N];
727  double* amp = new double[N];
728  double dt = 1./x.rate();
729  for(int i=0;i<N;i++) {
730  time[i]=i*dt;
731  amp[i]=x[i];
732  }
733  TGraph gr(N,time,amp);
734 
735  // fit signal envelope with gaussian function
736  gr.Fit("gaus","q0");
737  TF1* gaus = (TF1*)gr.GetFunction("gaus");
738 
739  mean = gaus->GetParameter("Mean"); // extract sigma from gauss fit
740  sigma = gaus->GetParameter("Sigma"); // extract sigma from gauss fit
741 
742  delete [] time;
743  delete [] amp;
744 
745  return;
746 }
747 
748 void
749 GetGaussFitPars2(wavearray<double>* wf, double& mean, double& sigma, bool fixAmax) {
750 //
751 // wf -> input waveform
752 // fixAmax -> true/false : fix/not-fix amplitude of the gaussian fit
753 //
754 // mean, sigma -> output gaussian fit parameters
755 //
756 //
757 // -> get the square of the evelope from the input waveform wf
758 // -> perform fit with gaussian function and get mean and gigma of gaussian fit
759 // -> compute the range [xmin,xmax] = [mean-3*sigma,mean-3*sigma]]
760 //
761 // -> get the evelope from the input waveform wf
762 // -> perform fit with gaussian function in the range [xmin,xmax] fixing the amplitude to max envelope (according to fixAmax value)
763 // -> return mean and gigma of gaussian fit
764 //
765 
766  wavearray<double> x = CWB::mdc::GetEnvelope(wf); // get envelope
767 
768  // compute range [xmin,xmax] used for final fit using the square of x
769  // x*x highlight the main peak
770  wavearray<double> x2 = x; // getx*x
771  for(int i=0;i<x.size();i++) x2[i]=x[i]*x[i];
772 
773  double gmean, gsigma;
774  GetGaussFitPars(&x2, gmean, gsigma, false);
775 // cout << "GetGaussFitPars2 -> gmean x2 " << gmean << " gsigma x2 " << gsigma << endl;
776 
777  double xmin = gmean-3.0*gsigma;
778  double xmax = gmean+3.0*gsigma;
779 // cout << "GetGaussFitPars2 -> xmin " << xmin << " xmax " << xmax << endl;
780 
781  // fill time/amp graph
782  int N = x.size();
783  double* time = new double[N];
784  double* amp = new double[N];
785  double dt = 1./x.rate();
786  double amax=0;
787  for(int i=0;i<N;i++) {
788  time[i]=i*dt;
789  amp[i]=x[i];
790  if(fabs(x[i])>amax) amax=fabs(x[i]);
791  }
792  TGraph gr(N,time,amp);
793 // cout << "GetGaussFitPars2 -> amax " << amax << endl;
794 
795  // fit signal envelope with gaussian function
796  TF1 *xgaus = new TF1("xgaus","gaus",xmin,xmax);
797  if(fixAmax) xgaus->SetParLimits(0,amax*0.999,amax*1.001); // fix gauss constant par to amax
798  gr.Fit("xgaus","RQ");
799 
800  mean = xgaus->GetParameter("Mean"); // extract sigma from gauss fit
801  sigma = xgaus->GetParameter("Sigma"); // extract sigma from gauss fit
802 
803  delete xgaus;
804  delete [] time;
805  delete [] amp;
806 
807  return;
808 }
809 
wavearray< double > t(hp.size())
monster wdmMRA
list of pixel pointers for MRA
Definition: network.hh:651
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
static float _sse_abs_ps(__m128 *_a)
Definition: watsse.hh:137
TString outDump
#define NIFO
Definition: wat.hh:74
size_t nLag
Definition: network.hh:573
void GetLveto(netcluster *pwc, int cid, int nifo, float *Lveto)
#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
void ResetPCA(network *NET, CWB::config *cfg, netcluster *pwc, std::vector< netpixel > *vPIX, int ID)
int n
Definition: cwb_net.C:28
static void _sse_zero_ps(__m128 *_p)
Definition: watsse.hh:44
wavearray< float > a_90
buffer for cluster sky 00 amplitude
Definition: network.hh:653
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
double frequency
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
int _sse_mra_ps(float *, float *, float, int)
Definition: network.hh:1276
double fLow
Definition: config.hh:140
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
bool cedDump
Definition: config.hh:297
waveform wf
Long_t size
WSeries< double > waveBand
Definition: detector.hh:356
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
bool getMRAwave(size_t ID, size_t lag, char atype='S', int mode=0, bool tof=false)
Definition: network.cc:3666
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
std::vector< detector * > ifoList
Definition: network.hh:608
#define N
bool outPlugin
Definition: config.hh:369
std::vector< netpixel > DoPCA(network *NET, CWB::config *cfg, int lag, int id)
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())
#define nIFO
double GetPeakFrequency(wavearray< double > *wf)
#define NTHR
wavearray< float > pNRG
buffers for cluster residual energy
Definition: network.hh:655
virtual size_t size() const
Definition: wavearray.hh:145
void CWB_Plugin(TFile *jfile, CWB::config *cfg, network *NET, WSeries< double > *x, TString ifo, int type)
COHERENCE.
TCanvas * canvas
Definition: watplot.hh:192
Float_t * frequency
GPS stop time of the cluster.
Definition: netevent.hh:102
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
TGraph * gr
int simulation
Definition: config.hh:199
#define IMPORT(TYPE, VAR)
Definition: cwb.hh:69
static wavearray< double > GetEnvelope(wavearray< double > *x)
Definition: mdc.cc:7828
watplot * WTS
Definition: ChirpMass.C:119
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
i() int(T_cor *100))
network NET
Definition: cwb_dump_inj.C:30
TString label
Definition: MergeTrees.C:21
double Pi
void ClearWaveforms(detector *ifo)
const int NIFO_MAX
Definition: wat.hh:22
bool log
Definition: WaveMDC.C:41
int BATCH
Definition: config.hh:245
double GetQfactor(wavearray< double > *wf, double frequency, bool fixAmax)
wavearray< float > a_00
wdm multi-resolution analysis
Definition: network.hh:652
TIter next(twave->GetListOfBranches())
char fname[1024]
std::vector< int > RWFID
Definition: detector.hh:381
int k
std::vector< wavearray< double > * > IWFP
Definition: detector.hh:379
double acor
Definition: network.hh:585
double sigma
std::vector< int > IWFID
Definition: detector.hh:378
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
void PlotWaveform(TString ifo, wavearray< double > *wfREC, CWB::config *cfg, bool fft=false, bool strain=false)
bool wfsave
Definition: network.hh:600
netevent EVT(itree, nifo)
std::vector< int > sCuts
Definition: netcluster.hh:392
WSeries< double > waveForm
Definition: detector.hh:355
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
Float_t * sSNR
data-signal correlation Xk*Sk
Definition: netevent.hh:136
TBranch * branch
int gIFACTOR
Float_t likelihood
Definition: netevent.hh:124
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)
void GetQveto(wavearray< double > *wf, float &Qveto, float &Qfactor)
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
void GetGaussFitPars2(wavearray< double > *wf, double &mean, double &sigma, bool fixAmax)
float rate
Definition: netpixel.hh:113
void GetGaussFitPars(wavearray< double > *wf, double &mean, double &sigma, bool doenv=true)
double getdata(char type='R', size_t n=0)
Definition: netpixel.hh:74
virtual void resize(unsigned int)
Definition: wavearray.cc:463
int check
exit(0)
bool setdata(double a, char type='R', size_t n=0)
Definition: netpixel.hh:58
bool dump
Definition: config.hh:295