Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
CWB_Plugin_WavePeaks.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 SAVE_WHT_PLOT // enable event WHITE plots
21 
22 #define nPEAK 10 // number of peaks to be extracted from the reconstructed waveform
23 
24 void GetGlitchParams(wavearray<double>* wf, int ifoID, float PF[][nPEAK],
25  float PB[][nPEAK], float PA[][nPEAK], float PE[][nPEAK]);
27  CWB::config* cfg, bool fft=false, bool strain=false);
28 
29 void
31 //!MISCELLANEA
32 // Extract whitened reconstructed waveforms, and compute the first nPEAK parameters
33 // Save peak parameters to the output wave root file
34 
35  cout << endl;
36  cout << "-----> CWB_Plugin_WavePeaks.C" << endl;
37  cout << "ifo " << ifo.Data() << endl;
38  cout << "type " << type << endl;
39  cout << endl;
40 
41  int nP = nPEAK; // number of peaks
42  float PF[NIFO_MAX][nPEAK]; // peak frequency
43  float PB[NIFO_MAX][nPEAK]; // peak bandwidth
44  float PA[NIFO_MAX][nPEAK]; // peak amplitude/(max signal amplitude)
45  float PE[NIFO_MAX][nPEAK]; // peak energy/(signal energy)
46 
47  if(type==CWB_PLUGIN_CONFIG) {
48  cfg->outPlugin=true; // disable built-in output root file
49  }
50 
51  if(type==CWB_PLUGIN_ILIKELIHOOD) {
52  NET->wfsave=true; // enable waveform save
53  }
54 
55  if(type==CWB_PLUGIN_OLIKELIHOOD) {
56 
57  if(TString(cfg->analysis)!="2G") {
58  cout << "CWB_Plugin_WavePeaks.C -> "
59  << "CWB_PLUGIN_OLIKELIHOOD implemented only for 2G" << endl;
60  gSystem->Exit(1);
61  }
62 
63  // import ifactor
64  int gIFACTOR=-1; IMPORT(int,gIFACTOR)
65  cout << "-----> CWB_Plugin_WavePeaks.C -> "
66  << " gIFACTOR : " << gIFACTOR << endl;
67 
68  // import slagShift
69  float* gSLAGSHIFT=NULL; IMPORT(float*,gSLAGSHIFT)
70 
71  int nIFO = NET->ifoListSize(); // number of detectors
72  int K = NET->nLag; // number of time lag
73  netevent* EVT;
75  double factor = cfg->simulation==3||cfg->simulation==4 ? -gIFACTOR : cfg->factors[gIFACTOR];
76  int rate = 0; // select all resolutions
77 
78  // search output root file in the system list
79  TFile* froot = NULL;
80  TList *files = (TList*)gROOT->GetListOfFiles();
81  TString outDump="";
82  if (files) {
83  TIter next(files);
84  TSystemFile *file;
85  TString fname;
86  bool check=false;
87  while ((file=(TSystemFile*)next())) {
88  fname = file->GetName();
89  // set output root file as the current file
90  if(fname.Contains("wave_")) {
91  froot=(TFile*)file;froot->cd();
92  outDump=fname;
93  outDump.ReplaceAll(".root.tmp",".txt");
94  //cout << "output file name : " << fname << endl;
95  }
96  }
97  if(!froot) {
98  cout << "CWB_Plugin_WavePeaks.C : Error - output root file not found" << endl;
99  gSystem->Exit(1);
100  }
101  } else {
102  cout << "CWB_Plugin_WavePeaks.C : Error - output root file not found" << endl;
103  gSystem->Exit(1);
104  }
105 
106  TTree* net_tree = (TTree *) froot->Get("waveburst");
107  if(net_tree!=NULL) {
108  EVT = new netevent(net_tree,nIFO);
109  net_tree->SetBranchAddress("nP",&nP);
110  net_tree->SetBranchAddress("PF",PF);
111  net_tree->SetBranchAddress("PB",PB);
112  net_tree->SetBranchAddress("PA",PA);
113  net_tree->SetBranchAddress("PE",PE);
114  } else {
115  EVT = new netevent(nIFO);
116  net_tree = EVT->setTree();
117  net_tree->Branch("nP",&nP,"nP/I");
118  net_tree->Branch("PF",PF,TString::Format("PF[%i][%i]/F",cfg->nIFO,nPEAK));
119  net_tree->Branch("PB",PB,TString::Format("PB[%i][%i]/F",cfg->nIFO,nPEAK));
120  net_tree->Branch("PA",PA,TString::Format("PA[%i][%i]/F",cfg->nIFO,nPEAK));
121  net_tree->Branch("PE",PE,TString::Format("PE[%i][%i]/F",cfg->nIFO,nPEAK));
122  }
123  EVT->setSLags(gSLAGSHIFT); // set slags into netevent
124 
125  for(int k=0; k<K; k++) { // loop over the lags
126 
127  id = NET->getwc(k)->get(const_cast<char*>("ID"), 0, 'L', rate);
128 
129  for(int j=0; j<(int)id.size(); j++) { // loop over cluster index
130 
131  int ID = size_t(id.data[j]+0.5);
132 
133  if(NET->getwc(k)->sCuts[ID-1]!=-1) continue; // skip rejected/processed clusters
134 
135  double ofactor=0;
136  if(cfg->simulation==4) ofactor=-factor;
137  else if(cfg->simulation==3) ofactor=-gIFACTOR;
138  else ofactor=factor;
139 
140  EVT->output2G(NULL,NET,ID,k,ofactor); // get reconstructed parameters
141 
142  wavearray<double>** pwfREC = new wavearray<double>*[nIFO];
143  detector* pd = NET->getifo(0);
144  int idSize = pd->RWFID.size();
145 
146  int wfIndex=-1;
147  for (int mm=0; mm<idSize; mm++) if (pd->RWFID[mm]==ID) wfIndex=mm;
148  if(wfIndex==-1) continue;
149 
150  // extract whitened reconstructed waveforms
151  for(int n=0; n<nIFO; n++) {
152 
153  pd = NET->getifo(n);
154 
155  pwfREC[n] = pd->RWFP[wfIndex];
156  wavearray<double>* wfREC = pwfREC[n]; // array of reconstructed waveforms
157 
158 #ifdef SAVE_WHT_PLOT
159  //PlotWaveform(NET->ifoName[n], wfREC, cfg, false, false);
160  PlotWaveform(NET->ifoName[n], wfREC, cfg, true, false);
161 #endif
162 
163  for(int i=0;i<nPEAK;i++) {
164  GetGlitchParams(wfREC, n, PF, PB, PA, PE);
165 // printf("%d -> \t %5.1f \t %5.1f \t %3.2f \t %3.2f\n",
166 // i,PF[n][i],PB[n][i],PA[n][i],PE[n][i]);
167  }
168  }
169  delete [] pwfREC;
170 
171  std::vector<int> sCuts = NET->getwc(k)->sCuts; // save cCuts
172  // set sCuts=1 to the events which must be not copied with cps to pwc
173  for(int i=0; i<(int)sCuts.size(); i++) if(i!=ID-1) NET->getwc(k)->sCuts[i]=1;
174 
175  // ID can not be used to get the event, to get event use ID=1 (ex: for watplot)
176  NET->getwc(k)->sCuts = sCuts; // restore cCuts
177 
178  if(cfg->dump) EVT->dopen(outDump.Data(),const_cast<char*>("a"),false);
179  EVT->output2G(net_tree,NET,ID,k,ofactor); // get reconstructed parameters
180  if(cfg->dump) EVT->dclose();
181  if(!cfg->cedDump) NET->getwc(k)->sCuts[ID-1]=1; // mark as processed
182  }
183  }
184 
185  jfile->cd();
186  if(EVT) delete EVT;
187  }
188  return;
189 }
190 
191 void
193  float PB[][nPEAK], float PA[][nPEAK], float PE[][nPEAK]) {
194 
195  int size = wf->size()/2;
196  double Fs=((double)wf->rate()/(double)(2*size));
197 
198  wavearray<float> ener(size);
199 
200  wf->FFTW(1);
201  for(int j=0;j<size;j++) ener[j]=pow(wf->data[2*j],2)+pow(wf->data[2*j+1],2);
202  wf->FFTW(-1);
203  wf->resetFFTW();
204 
205  // compute total signal energy
206  double ET=0;
207  for(int j=0;j<size;j++) ET+=ener[j];
208  // compute max signal amplitude
209  double MA=0;
210  for(int j=0;j<size;j++) if(ener[j]>MA) MA=ener[j];
211  MA=sqrt(MA);
212 
213  for(int k=0;k<nPEAK;k++) {
214 
215  // find frequency @ max energy
216  float ME=0; // peak max energy
217  int jME=0; // peak max energy index
218  for(int j=0;j<size;j++) if(ener[j]>ME) {ME=ener[j];jME=j;}
219 
220  // compute peak energy
221  float fmean=ME*Fs*jME;
222  float frms=ME*pow(Fs*jME,2);
223  double EP=ME;
224  double ge=ME;
225  ener[jME]=0;
226  float finf=0;
227  for(int j=jME-1;j>=0;j--) {
228  float E = ener[j];
229  float F = Fs*j;
230  if(E<ge) {EP+=E;finf=F;fmean+=E*F;frms+=E*F*F;ener[j]=0;} else break;
231  ge=E;
232  }
233  ge = ME;
234  float fsup=size*Fs;
235  for(int j=jME+1;j<size;j++) {
236  float E = ener[j];
237  float F = Fs*j;
238  if(E<ge) {EP+=E;fsup=F;fmean+=E*F;frms+=E*F*F;ener[j]=0;} else break;
239  ge=E;
240  }
241 
242  fmean = fmean/EP;
243  frms = frms/EP-fmean*fmean;
244  frms = frms>0 ? sqrt(frms) : 0;
245 
246  PF[ifoID][k] = fmean;
247  PA[ifoID][k] = sqrt(ME)/MA;
248  PB[ifoID][k] = frms;
249  PE[ifoID][k] = ET>0 ? EP/ET : 0;
250  }
251 }
252 
253 void
255  CWB::config* cfg, bool fft, bool strain) {
256 
257  watplot PTS(const_cast<char*>("ptswrc"),200,20,800,500);
258 
259  //cout << "Print " << fname << endl;
260  double tmin = wfREC->start();
261  wfREC->start(wfREC->start()-tmin);
262  if(fft) {
263  PTS.plot(wfREC, const_cast<char*>("ALP"), 1, 0, 0, true, cfg->fLow, cfg->fHigh);
264  } else {
265  PTS.plot(wfREC, const_cast<char*>("ALP"), 1, 0, 0);
266  }
267  PTS.graph[0]->SetLineWidth(1);
268  wfREC->start(wfREC->start()+tmin);
269 
270  char label[64]="";
271  if(fft) sprintf(label,"%s","fft");
272  else sprintf(label,"%s","time");
273  if(strain) sprintf(label,"%s_%s",label,"strain");
274  else sprintf(label,"%s_%s",label,"white");
275 
276  char fname[1024];
277  sprintf(fname, "%s_wf_%s_rec_gps_%d.root",ifo.Data(),label,int(tmin));
278  PTS.canvas->Print(fname);
279  cout << "write : " << fname << endl;
280  //PTS.canvas->Write(REPLACE(fname,dirCED,gtype));
281 }
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
void setSLags(float *slag)
Definition: netevent.cc:404
#define nPEAK
virtual void rate(double r)
Definition: wavearray.hh:123
void PlotWaveform(TString ifo, wavearray< double > *wfREC, CWB::config *cfg, bool fft=false, bool strain=false)
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
std::vector< wavearray< double > * > RWFP
Definition: detector.hh:364
double fLow
Definition: config.hh:122
std::vector< TGraph * > graph
Definition: watplot.hh:176
bool cedDump
Definition: config.hh:279
waveform wf
Long_t size
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
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
#define nIFO
TCanvas * canvas
Definition: watplot.hh:174
double factor
jfile
Definition: cwb_job_obj.C:25
void GetGlitchParams(wavearray< double > *wf, int ifoID, float PF[][nPEAK], float PB[][nPEAK], float PA[][nPEAK], float PE[][nPEAK])
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
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
virtual void resetFFTW()
Definition: wavearray.cc:959
double F
TFile * froot
virtual void FFTW(int=1)
Definition: wavearray.cc:878
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
string file
Definition: cwb_online.py:385
netcluster * getwc(size_t n)
param: delay index
Definition: network.hh:421
sprintf(tfres,"(1/%g)x(%g) (sec)x(Hz)", 2 *df, df)
int nIFO
Definition: config.hh:102
DataType_t * data
Definition: wavearray.hh:301
Long_t id
double factors[FACTORS_MAX]
Definition: config.hh:184
list files
Definition: cwb_online.py:529
int check
void CWB_Plugin(TFile *jfile, CWB::config *cfg, network *NET, WSeries< double > *x, TString ifo, int type)
COHERENCE.
bool dump
Definition: config.hh:277