Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
injection.cc
Go to the documentation of this file.
1 #include <fstream>
2 #include "injection.hh"
3 #include "TH2.h"
4 #include "TStyle.h"
5 #include "TCanvas.h"
6 
7 #define MDC_TREE_NAME "mdc"
8 
9 ClassImp(injection) // used by THtml doc
10 
11 TTree* injection::Init(TString fName, int n)
12 {
13  iFile = TFile::Open(fName);
14  if((iFile==NULL) || (iFile!=NULL && !iFile->IsOpen())) {
15  cout << "injection::Init : Error opening root file " << fName.Data() << endl;
16  exit(1);
17  }
18 
19  TTree* tree = (TTree *) iFile->Get(MDC_TREE_NAME);
20  if(tree) {
21  ndim = tree->GetUserInfo()->GetSize(); // get number of detectors
22  if(ndim>0) {
23  if(n>0 && ndim!=n) {
24  cout << "injection::Init : number of detectors declared in the constructor (" << n
25  << ") are not equals to the one ("<<ndim<<") declared in the root file : "
26  << fName.Data() << endl;
27  exit(1);
28  }
29  } else ndim=n;
30  } else {
31  cout << "injection::Init : object tree " << MDC_TREE_NAME
32  << " not present in the root file " << fName.Data() << endl;
33  exit(1);
34  }
35 
36  if(ndim==0) {
37  cout << "injection::Init : number of detectors is not declared in the constructor or"
38  << " not present in the root file : " << endl << fName.Data() << endl;
39  exit(1);
40  }
41 
42  return tree;
43 }
44 
45 // Set branch addresses
46 void injection::Init(TTree *tree)
47 {
48  if (tree == 0) return;
49  fChain = tree;
50  fCurrent = -1;
51  fChain->SetMakeClass(1);
52 
53  fChain->SetBranchAddress("run",&run);
54  fChain->SetBranchAddress("ndim",&ndim);
55  fChain->SetBranchAddress("nevent",&nevent);
56  fChain->SetBranchAddress("eventID",&eventID);
57  fChain->SetBranchAddress("type",&type);
58  fChain->SetBranchAddress("name",&name);
59 
60  fChain->SetBranchAddress("factor",&factor);
61  fChain->SetBranchAddress("distance",&distance);
62  fChain->SetBranchAddress("mchirp",&mchirp);
63  fChain->SetBranchAddress("rp0",&rp0);
64  fChain->SetBranchAddress("e0",&e0);
65  fChain->SetBranchAddress("redshift",&redshift);
66  fChain->SetBranchAddress("gps",&gps);
67  fChain->SetBranchAddress("strain",&strain);
68  fChain->SetBranchAddress("psi",psi);
69  fChain->SetBranchAddress("iota",iota);
70  fChain->SetBranchAddress("phi",phi);
71  fChain->SetBranchAddress("theta",theta);
72  fChain->SetBranchAddress("bp",bp);
73  fChain->SetBranchAddress("bx",bx);
74 
75  fChain->SetBranchAddress("time",time);
76  fChain->SetBranchAddress("duration",duration);
77 
78  fChain->SetBranchAddress("frequency",frequency);
79  fChain->SetBranchAddress("bandwidth",bandwidth);
80  fChain->SetBranchAddress("hrss",hrss);
81  fChain->SetBranchAddress("snr",snr);
82  fChain->SetBranchAddress("Deff",Deff);
83  fChain->SetBranchAddress("mass",mass);
84  fChain->SetBranchAddress("spin",spin);
85 
86  Notify();
87 }
88 
89 // allocate memory
91 {
92  if (!name) name= new string();
93  else {delete name;name= new string();}
94  if (!psi) psi= (Float_t*)malloc(2*sizeof(Float_t));
95  else psi= (Float_t*)realloc(psi,2*sizeof(Float_t));
96  if (!iota) iota= (Float_t*)malloc(2*sizeof(Float_t));
97  else iota= (Float_t*)realloc(iota,2*sizeof(Float_t));
98  if (!phi) phi= (Float_t*)malloc(2*sizeof(Float_t));
99  else phi= (Float_t*)realloc(phi,2*sizeof(Float_t));
100  if (!theta) theta= (Float_t*)malloc(2*sizeof(Float_t));
101  else theta= (Float_t*)realloc(theta,2*sizeof(Float_t));
102  if (!bp) bp= (Float_t*)malloc(ndim*sizeof(Float_t));
103  else bp= (Float_t*)realloc(bp,ndim*sizeof(Float_t));
104  if (!bx) bx= (Float_t*)malloc(ndim*sizeof(Float_t));
105  else bx= (Float_t*)realloc(bx,ndim*sizeof(Float_t));
106 
107  if (!time) time= (Double_t*)malloc(ndim*sizeof(Double_t));
108  else time= (Double_t*)realloc(time,ndim*sizeof(Double_t));
109  if (!duration) duration= (Float_t*)malloc(ndim*sizeof(Float_t));
110  else duration= (Float_t*)realloc(duration,ndim*sizeof(Float_t));
111 
112  if (!frequency) frequency=(Float_t*)malloc(ndim*sizeof(Float_t));
113  else frequency=(Float_t*)realloc(frequency,ndim*sizeof(Float_t));
114  if (!bandwidth) bandwidth=(Float_t*)malloc(ndim*sizeof(Float_t));
115  else bandwidth=(Float_t*)realloc(bandwidth,ndim*sizeof(Float_t));
116  if (!hrss) hrss= (Double_t*)malloc(ndim*sizeof(Double_t));
117  else hrss= (Double_t*)realloc(hrss,ndim*sizeof(Double_t));
118  if (!snr) snr= (Float_t*)malloc(ndim*sizeof(Float_t));
119  else snr= (Float_t*)realloc(snr,ndim*sizeof(Float_t));
120  if (!Deff) Deff= (Float_t*)malloc(ndim*sizeof(Float_t));
121  else Deff= (Float_t*)realloc(Deff,ndim*sizeof(Float_t));
122  if (!mass) mass= (Float_t*)malloc(2*sizeof(Float_t));
123  else mass= (Float_t*)realloc(mass,2*sizeof(Float_t));
124  if (!spin) spin= (Float_t*)malloc(6*sizeof(Float_t));
125  else spin= (Float_t*)realloc(spin,6*sizeof(Float_t));
126  if (!pwf) pwf= (wavearray<double>**)malloc(ndim*sizeof(wavearray<double>*));
127  else pwf= (wavearray<double>**)realloc(pwf,ndim*sizeof(wavearray<double>*));
128 
129  for(int n=0; n<ndim; n++) pwf[n] = NULL;
130 
131  return;
132 }
133 
134 // init array
136 {
137 // for(int i=0;i<NIFO_MAX;i++) pD[i]=NULL;
138  return;
139 }
140 
142 {
143  // Called when loading a new file.
144  // Get branch pointers.
145  b_ndim = fChain->GetBranch("ndim");
146  b_run = fChain->GetBranch("run");
147  b_nevent = fChain->GetBranch("nevent");
148  b_eventID = fChain->GetBranch("eventID");
149  b_type = fChain->GetBranch("type");
150  b_name = fChain->GetBranch("name");
151  b_factor = fChain->GetBranch("factor");
152  b_distance = fChain->GetBranch("distance");
153  b_mchirp = fChain->GetBranch("mchirp");
154  b_rp0 = fChain->GetBranch("rp0");
155  b_e0 = fChain->GetBranch("e0");
156  b_redshift = fChain->GetBranch("redshift");
157  b_gps = fChain->GetBranch("gps");
158  b_strain = fChain->GetBranch("strain");
159  b_psi = fChain->GetBranch("psi");
160  b_iota = fChain->GetBranch("iota");
161  b_phi = fChain->GetBranch("phi");
162  b_theta = fChain->GetBranch("theta");
163  b_bp= fChain->GetBranch("bx");
164  b_bx = fChain->GetBranch("bp");
165  b_time = fChain->GetBranch("time");
166  b_duration = fChain->GetBranch("duration");
167  b_frequency = fChain->GetBranch("frequency");
168  b_bandwidth = fChain->GetBranch("bandwidth");
169  b_hrss = fChain->GetBranch("hrss");
170  b_snr = fChain->GetBranch("snr");
171  b_Deff = fChain->GetBranch("Deff");
172  b_mass = fChain->GetBranch("mass");
173  b_spin = fChain->GetBranch("spin");
174 
175  return kTRUE;
176 }
177 
179 {
180  if (!fChain) return 0;
181  return fChain->GetEntry(entry);
182 };
184 {
185  if (!fChain) return 0;
186  return fChain->GetEntries();
187 };
188 
190 {
191 // Print contents of entry.
192 // If entry is not specified, print current entry
193  if (!fChain) return;
194  fChain->Show(entry);
195 }
196 
197 
198 //++++++++++++++++++++++++++++++++++++++++++++++
199 // set single event tree
200 //++++++++++++++++++++++++++++++++++++++++++++++
202 {
203  TTree* waveTree = new TTree(MDC_TREE_NAME,MDC_TREE_NAME);
204 
205  char cpsi[16];
206  char ciota[16];
207  char cphi[16];
208  char ctheta[16];
209  char cbp[16];
210  char cbx[16];
211  char ctime[16];
212  char cduration[16];
213  char cfrequency[16];
214  char cbandwidth[16];
215  char chrss[16];
216  char csnr[16];
217  char cDeff[16];
218  char cmass[16];
219  char cspin[16];
220 
221  sprintf(cpsi, "psi[2]/F");
222  sprintf(ciota, "iota[2]/F");
223  sprintf(cphi, "phi[2]/F");
224  sprintf(ctheta, "theta[2]/F");
225  sprintf(cbp, "bp[%1d]/F",ndim);
226  sprintf(cbx, "bx[%1d]/F",ndim);
227  sprintf(ctime, "time[%1d]/D",ndim);
228  sprintf(cduration, "duration[%1d]/F",ndim);
229  sprintf(cfrequency, "frequency[%1d]/F",ndim);
230  sprintf(cbandwidth, "bandwidth[%1d]/F",ndim);
231  sprintf(chrss, "hrss[%1d]/D",ndim);
232  sprintf(csnr, "snr[%1d]/F",ndim);
233  sprintf(cDeff, "Deff[%1d]/F",ndim);
234  sprintf(cmass, "mass[2]/F");
235  sprintf(cspin, "spin[6]/F");
236 
237  //==================================
238  // Define trigger tree
239  //==================================
240 
241  waveTree->Branch("ndim", &ndim, "ndim/I");
242  waveTree->Branch("run", &run, "run/I");
243  waveTree->Branch("nevent", &nevent, "nevent/I");
244  waveTree->Branch("eventID", &eventID, "eventID/I");
245  waveTree->Branch("type", &type, "type/I");
246  waveTree->Branch("name", name);
247  waveTree->Branch("factor", &factor, "factor/F");
248  waveTree->Branch("distance", &distance, "distance/F");
249  waveTree->Branch("mchirp", &mchirp, "mchirp/F");
250  waveTree->Branch("rp0", &rp0, "rp0/F");
251  waveTree->Branch("e0", &e0, "e0/F");
252  waveTree->Branch("redshift", &redshift, "redshift/F");
253  waveTree->Branch("gps", &gps, "gps/D");
254  waveTree->Branch("strain", &strain, "strain/D");
255  waveTree->Branch("psi", psi, cpsi);
256  waveTree->Branch("iota", iota, ciota);
257  waveTree->Branch("phi", phi, cphi);
258  waveTree->Branch("theta", theta, ctheta);
259  waveTree->Branch("bp", bp, cbp);
260  waveTree->Branch("bx", bx, cbx);
261  waveTree->Branch("time", time, ctime);
262  waveTree->Branch("duration", duration, cduration);
263  waveTree->Branch("frequency", frequency, cfrequency);
264  waveTree->Branch("bandwidth", bandwidth, cbandwidth);
265  waveTree->Branch("hrss", hrss, chrss);
266  waveTree->Branch("snr", snr, csnr);
267  waveTree->Branch("Deff", Deff, cDeff);
268  waveTree->Branch("mass", mass, cmass);
269  waveTree->Branch("spin", spin, cspin);
270 
271  return waveTree;
272 }
273 
275 {
276  int i;
277 
278  if(this->ndim != a.ndim) { this->ndim=a.ndim; this->allocate(); }
279 
280  this->run= a.run;
281  this->nevent= a.nevent;
282  this->eventID= a.eventID;
283  this->factor= a.factor;
284  this->distance= a.distance;
285  this->mchirp= a.mchirp;
286  this->rp0= a.rp0;
287  this->e0= a.e0;
288  this->redshift= a.redshift;
289  this->gps= a.gps;
290  this->type= a.type;
291  *this->name= *a.name;
292  this->strain= a.strain;
293  this->psi[0]= a.psi[0];
294  this->iota[0]= a.iota[0];
295  this->phi[0]= a.phi[0];
296  this->theta[0]= a.theta[0];
297  this->psi[1]= a.psi[1];
298  this->phi[1]= a.phi[1];
299  this->theta[1]= a.theta[1];
300  this->mass[0]= a.mass[1];
301  this->mass[1]= a.mass[1];
302 
303  for(i=0; i<6; i++) this->spin[i] = a.spin[i];
304 
305  for(i=0; i<ndim; i++){
306 
307  this->bp[i]= a.bp[i];
308  this->bx[i]= a.bx[i];
309  this->time[i]= a.time[i];
310  this->duration[i]= a.duration[i];
311  this->frequency[i]= a.frequency[i];
312  this->bandwidth[i]= a.bandwidth[i];
313  this->hrss[i]= a.hrss[i];
314  this->snr[i]= a.snr[i];
315  this->Deff[i]= a.Deff[i];
316  }
317  return *this;
318 }
319 
320 
321 
322 //++++++++++++++++++++++++++++++++++++++++++++++++++++
323 // fill this with data
324 //++++++++++++++++++++++++++++++++++++++++++++++++++++
325 Bool_t injection::fill_in(network* net, int id, bool checkEdges)
326 {
327  bool save = true;
328 
329  size_t N = net->mdc__IDSize();
330  size_t I = net->mdcTypeSize();
331  size_t M = net->ifoListSize();
332 
333  int ID = id<0 ? abs(id+1) : abs(id);
334 
335  if(!N || !M || !I || ID<0) return false;
336 
337  size_t i,m,nst;
338  string itag;
339  char* p;
340  char ch[1024];
341  double hphp, hxhx, hphx, T0, To, Eo;
342  double Pi = 3.14159265358979312;
343  detector* pd;
344 
345  this->gps = net->getifo(0)->getTFmap()->start();
346  size_t K = net->getifo(0)->getTFmap()->size();
347  double R = net->getifo(0)->getTFmap()->wavearray<double>::rate();
348  double sTARt = this->gps + net->Edge + 1.;
349  double sTOp = this->gps + K/R - net->Edge - 1.;
350 
351  if((net->getmdcTime(ID)<sTARt || net->getmdcTime(ID)>sTOp) && checkEdges) return false;
352 
353  this->eventID = ID;
354 
355  string str(net->getmdcList(ID));
356  sprintf(ch,"%s",str.c_str());
357 
358  if((p = strtok(ch," \t")) == NULL) return false;
359 
360  p = strtok(NULL," \t"); // get MDC strain
361  this->strain = atof(p);
362 
363  p = strtok(NULL," \t"); // skip
364  p = strtok(NULL," \t"); // skip
365  p = strtok(NULL," \t"); // get internal phi
366  this->iota[0] = atof(p); // get iota
367  this->phi[1] = atof(p);
368  p = strtok(NULL," \t"); // get internal psi
369  this->psi[1] = atof(p);
370  p = strtok(NULL," \t"); // get external theta
371  if(fabs(atof(p))>1) {
372  cout<<"injection:fill_in error: external theta not valid, must be [-1,1]\n"<<endl;
373  exit(1);
374  }
375  this->theta[0] = acos(atof(p));
376  this->theta[0]*= 180/Pi;
377  p = strtok(NULL," \t"); // get external phi
378  this->phi[0] = atof(p) > 0 ? atof(p) : 2*Pi+atof(p);
379  this->phi[0]*= 180/Pi;
380  p = strtok(NULL," \t"); // get external psi
381  this->psi[0] = atof(p);
382  this->psi[0]*= 180/Pi;
383 
384  p = strtok(NULL," \t"); // skip
385  p = strtok(NULL," \t"); // injection time
386 
387  // printf("%12.3f %12.3f %12.3f \n",this->time[0],sTARt,sTOp);
388 
389  p = strtok(NULL," \t"); // injection name
390 
391  this->name->assign(p, strlen(p));
392 
393  for(i=0; i<I; i++) {
394  itag = net->getmdcType(i);
395  if(itag.find(p) == string::npos) continue;
396  this->type = i+1; break;
397  }
398 
399  p = strtok(NULL," \t"); // h+h+
400  hphp = atof(p);
401  p = strtok(NULL," \t"); // hxhx
402  hxhx = atof(p);
403  p = strtok(NULL," \t"); // h+hx
404  hphx = atof(p);
405 
406  save = true;
407  for(m=0; m<M; m++) { // loop over detectors
408  nst = str.find(net->getifo(m)->Name);
409  if(nst >= str.length()) {
410  cout<<"injection:fill_in error: no injections for detector "
411  << net->getifo(m)->Name <<" was found in the injection list !!!\n\n";
412  save = false;
413  exit(1);
414  }
415 
416  itag = str.substr(nst);
417  sprintf(ch,"%s",itag.c_str());
418  if((p = strtok(ch," \t")) == NULL) continue; // detector name
419 
420  p = strtok(NULL," \t"); // time
421  this->time[m] = atof(p);
422 // if(this->time[m]<sTARt || this->time[m]>sTOp) save = false;
423 
424 // printf("%13.3f %13.3f %13.3f\n",sTARt,time[m],sTOp);
425 
426  p = strtok(NULL," \t"); // F+
427  this->bp[m] = atof(p);
428  p = strtok(NULL," \t"); // Fx
429  this->bx[m] = atof(p);
430 
431  nst = str.find("insp") < str.find("ebbh") ? str.find("insp") : str.find("ebbh");
432  this->Deff[m] = 0.;
433  if(nst < str.length()) { // inspiral MDC log
434  p = strtok(NULL," \t"); // effective distance
435  this->Deff[m] = atof(p);
436  }
437  else { // standard burst MDC log
438  this->hrss[m] = hphp*bp[m]*bp[m] + hxhx*bx[m]*bx[m] + 2*hphx*bp[m]*bx[m];
439  this->hrss[m] = this->hrss[m]>0. ? sqrt(this->hrss[m]) : 0.;
440  }
441  }
442 
443  if(!save) return save;
444 
445  nst = str.find("distance");
446  this->distance = 0.;
447  if(nst < str.length()) {
448  itag = str.substr(nst);
449  sprintf(ch,"%s",itag.c_str());
450  if((p = strtok(ch," \t")) != NULL)
451  this->distance = atof(strtok(NULL," \t"));
452  }
453 
454  nst = str.find("mass1");
455  this->mass[0] = 0.;
456  if(nst < str.length()) {
457  itag = str.substr(nst);
458  sprintf(ch,"%s",itag.c_str());
459  if((p = strtok(ch," \t")) != NULL)
460  this->mass[0] = atof(strtok(NULL," \t"));
461  }
462 
463  nst = str.find("mass2");
464  this->mass[1] = 0.;
465  if(nst < str.length()) {
466  itag = str.substr(nst);
467  sprintf(ch,"%s",itag.c_str());
468  if((p = strtok(ch," \t")) != NULL)
469  this->mass[1] = atof(strtok(NULL," \t"));
470  }
471 
472  nst = str.find("mchirp");
473  this->mchirp = 0.;
474  if(nst < str.length()) {
475  itag = str.substr(nst);
476  sprintf(ch,"%s",itag.c_str());
477  if((p = strtok(ch," \t")) != NULL)
478  this->mchirp = atof(strtok(NULL," \t"));
479  }
480 
481  nst = str.find("rp0");
482  this->rp0 = 0.;
483  if(nst < str.length()) {
484  itag = str.substr(nst);
485  sprintf(ch,"%s",itag.c_str());
486  if((p = strtok(ch," \t")) != NULL)
487  this->rp0 = atof(strtok(NULL," \t"));
488  }
489 
490  nst = str.find("e0");
491  this->e0 = 0.;
492  if(nst < str.length()) {
493  itag = str.substr(nst);
494  sprintf(ch,"%s",itag.c_str());
495  if((p = strtok(ch," \t")) != NULL)
496  this->e0 = atof(strtok(NULL," \t"));
497  }
498 
499  nst = str.find("redshift");
500  this->redshift = 0.;
501  if(nst < str.length()) {
502  itag = str.substr(nst);
503  sprintf(ch,"%s",itag.c_str());
504  if((p = strtok(ch," \t")) != NULL)
505  this->redshift = atof(strtok(NULL," \t"));
506  }
507 
508  nst = str.find("spin1");
509  this->spin[0] = 0.;
510  this->spin[1] = 0.;
511  this->spin[2] = 0.;
512  if(nst < str.length()) {
513  itag = str.substr(nst);
514  sprintf(ch,"%s",itag.c_str());
515  if((p = strtok(ch," \t")) != NULL) {
516  this->spin[0] = atof(strtok(NULL," \t"));
517  this->spin[1] = atof(strtok(NULL," \t"));
518  this->spin[2] = atof(strtok(NULL," \t"));
519  }
520  }
521 
522  nst = str.find("spin2");
523  this->spin[3] = 0.;
524  this->spin[4] = 0.;
525  this->spin[5] = 0.;
526  if(nst < str.length()) {
527  itag = str.substr(nst);
528  sprintf(ch,"%s",itag.c_str());
529  if((p = strtok(ch," \t")) != NULL) {
530  this->spin[3] = atof(strtok(NULL," \t"));
531  this->spin[4] = atof(strtok(NULL," \t"));
532  this->spin[5] = atof(strtok(NULL," \t"));
533  }
534  }
535 
536 // read simulation parameters from the detector objects
537 
538  To = Eo = 0.;
539  for(m=0; m<M; m++) { // loop over detectors
540  pd = net->getifo(m);
541  if(!pd->HRSS.size()) continue;
542  To += pd->TIME.data[ID]*pd->ISNR.data[ID]*pd->ISNR.data[ID];
543  Eo += pd->ISNR.data[ID]*pd->ISNR.data[ID];
544  if(pd->TIME.data[ID] < 1.) save = false;
545  }
546 
547  if(!save || Eo<=0.) return save;
548 
549  To /= Eo; // central injection time
550  T0 = this->time[0];
551  for(m=0; m<M; m++) { // loop over detectors
552  pd = net->getifo(m);
553  if(!pd->HRSS.size()) continue;
554  this->frequency[m] = pd->FREQ.data[ID];
555  this->bandwidth[m] = pd->BAND.data[ID];
556  this->time[m] += To - T0;
557  this->duration[m] = pd->TDUR.data[ID];
558  this->hrss[m] = pd->HRSS.data[ID];
559  this->snr[m] = sqrt(pd->ISNR.data[ID]);
560 
561  int idSize = pd->IWFID.size();
562  int wfIndex=-1;
563  for (int mm=0; mm<idSize; mm++) if (pd->IWFID[mm]==id) wfIndex=mm;
564  this->pwf[m] = wfIndex>=0 ? pd->IWFP[wfIndex] : NULL;
565  }
566 
567  return save;
568 }
569 
570 //++++++++++++++++++++++++++++++++++++++++++++++++++++
571 // output injection metadata
572 //++++++++++++++++++++++++++++++++++++++++++++++++++++
573 void injection::output(TTree* waveTree, network* net, double factor, bool checkEdges)
574 {
575  size_t N = net->mdc__IDSize();
576  size_t I = net->mdcTypeSize();
577  size_t M = net->ifoListSize();
578  if(!N || !M || !I) return;
579  double FACTOR = fabs(factor); // used to tag events in root file
580  factor = factor<=0 ? 1 : fabs(factor); // used to rescale injected events
581 
582 
583  size_t n,m;
584 
585  // add detectors to tree user info
586  if(waveTree!=NULL) {
587  int dsize = waveTree->GetUserInfo()->GetSize();
588  if(dsize!=0 && dsize!=M) {
589  cout<<"injection::output(): wrong user detector list in header tree"<<endl; exit(1);
590  }
591  if(dsize==0) {
592  for(int n=0;n<M;n++) {
593  // must be done a copy because detector object
594  // is destroyed when waveTree is destroyed
595  detector* pD = new detector(*net->getifo(n));
596  waveTree->GetUserInfo()->Add(pD);
597  }
598  }
599  }
600 
601  this->run = net->nRun;
602  this->nevent = 0;
603 
604  for(n=0; n<N; n++) {
605  this->eventID = net->getmdc__ID(n);
606  if(this->fill_in(net,this->eventID,checkEdges)) {
607  this->nevent++;
608  this->factor = FACTOR;
609  this->strain *= factor;
610  this->distance /= factor;
611  for(m=0; m<M; m++) {
612  this->hrss[m] *= factor;
613  this->snr[m] *= factor;
614  this->Deff[m] /= factor;
615  }
616  waveTree->Fill();
617  }
618  }
619 }
620 
621 
622 
623 /*
624 Int_t injection::LoadTree(Int_t entry)
625 {
626 // Set the environment to read one entry
627  if (!fChain) return -5;
628  Int_t centry = fChain->LoadTree(entry);
629  if (centry < 0) return centry;
630  if (fChain->IsA() != TChain::Class()) return centry;
631  TChain *chain = (TChain*)fChain;
632  if (chain->GetTreeNumber() != fCurrent) {
633  fCurrent = chain->GetTreeNumber();
634  Notify();
635  }
636  return centry;
637 }
638 
639 Int_t injection::Cut(Int_t entry)
640 {
641 // This function may be called from Loop.
642 // returns 1 if entry is accepted.
643 // returns -1 otherwise.
644  return 1;
645 }
646 
647 void injection::Loop()
648 {
649 // In a ROOT session, you can do:
650 // Root > .L injection.C
651 // Root > injection t
652 // Root > t.GetEntry(12); // Fill t data members with entry number 12
653 // Root > t.Show(); // Show values of entry 12
654 // Root > t.Show(16); // Read and show values of entry 16
655 // Root > t.Loop(); // Loop on all entries
656 //
657 
658 // This is the loop skeleton where:
659 // jentry is the global entry number in the chain
660 // ientry is the entry number in the current Tree
661 // Note that the argument to GetEntry must be:
662 // jentry for TChain::GetEntry
663 // ientry for TTree::GetEntry and TBranch::GetEntry
664 //
665 // To read only selected branches, Insert statements like:
666 // METHOD1:
667 // fChain->SetBranchStatus("*",0); // disable all branches
668 // fChain->SetBranchStatus("branchname",1); // activate branchname
669 // METHOD2: replace line
670 // fChain->GetEntry(jentry); //read all branches
671 //by b_branchname->GetEntry(ientry); //read only this branch
672  if (fChain == 0) return;
673 
674  Int_t nentries = Int_t(fChain->GetEntriesFast());
675 
676  Int_t nbytes = 0, nb = 0;
677  for (Int_t jentry=0; jentry<nentries;jentry++) {
678  Int_t ientry = LoadTree(jentry);
679  if (ientry < 0) break;
680  nb = fChain->GetEntry(jentry); nbytes += nb;
681  // if (Cut(ientry) < 0) continue;
682  }
683 }
684 */
685 
686 
687 
688 
689 
690 
691 
692 
Int_t fCurrent
pointer to the analyzed TTree or TChain
Definition: injection.hh:42
TTree * tree
Definition: TimeSortTree.C:20
size_t mdcTypeSize()
Definition: network.hh:392
detector * getifo(size_t n)
param: detector index
Definition: network.hh:418
TBranch * b_nevent
Definition: injection.hh:88
TBranch * b_hrss
Definition: injection.hh:114
TTree * Init(TString fName, int n)
Definition: injection.cc:11
double sTARt
virtual size_t size() const
Definition: wavearray.hh:127
Double_t * time
beam pattern coefficients for hx
Definition: injection.hh:70
Float_t rp0
Definition: injection.hh:57
string * name
Definition: injection.hh:51
TBranch * b_psi
Definition: injection.hh:101
Float_t e0
Definition: injection.hh:58
void Init()
Definition: ChirpMass.C:280
TBranch * b_bandwidth
Definition: injection.hh:112
Int_t GetEntry(Int_t)
Definition: injection.cc:178
std::vector< double > * getmdcTime()
Definition: network.hh:404
TBranch * b_distance
Definition: injection.hh:94
wavearray< double > a(hp.size())
int n
Definition: cwb_net.C:10
Bool_t Notify()
Definition: injection.cc:141
TString("c")
void output(TTree *, network *, double, bool=true)
Definition: injection.cc:573
size_t nRun
Definition: network.hh:554
int ID
Definition: TestMDC.C:70
wavearray< double > HRSS
Definition: detector.hh:353
Float_t * snr
injected hrss in the detectors
Definition: injection.hh:77
TBranch * b_Deff
Definition: injection.hh:116
TBranch * b_name
Definition: injection.hh:91
string getmdcType(size_t n)
Definition: network.hh:402
injection() allocate()
Float_t * iota
source psi angle
Definition: injection.hh:64
string getmdcList(size_t n)
Definition: network.hh:400
TBranch * b_duration
Definition: injection.hh:109
Float_t * duration
injection gps time
Definition: injection.hh:71
Float_t * bandwidth
average center_of_hrss frequency
Definition: injection.hh:74
TTree * setTree()
Definition: injection.cc:201
#define M
Definition: UniqSLagsList.C:3
int m
Definition: cwb_net.C:10
virtual void start(double s)
Definition: wavearray.hh:119
Double_t strain
Definition: injection.hh:62
i drho i
TBranch * b_snr
Definition: injection.hh:115
TBranch * b_theta
Definition: injection.hh:104
#define N
TBranch * b_bp
Definition: injection.hh:105
Float_t * theta
source phi angle
Definition: injection.hh:66
Double_t gps
Definition: injection.hh:61
Int_t ndim
current Tree number in a TChain
Definition: injection.hh:43
TBranch * b_redshift
Definition: injection.hh:98
double Edge
Definition: network.hh:560
network ** net
NOISE_MDC_SIMULATION.
size_t ifoListSize()
Definition: network.hh:413
TBranch * b_e0
Definition: injection.hh:97
TBranch * b_iota
Definition: injection.hh:102
TBranch * b_spin
Definition: injection.hh:118
double factor
Float_t redshift
Definition: injection.hh:59
char str[1024]
Float_t * psi
Definition: injection.hh:63
void Show(Int_t entry=-1)
Definition: injection.cc:189
Float_t factor
injection name
Definition: injection.hh:53
TBranch * b_type
Definition: injection.hh:90
Double_t * hrss
estimated bandwidth
Definition: injection.hh:76
virtual injection & operator=(const injection &)
Definition: injection.cc:274
wavearray< double > TIME
Definition: detector.hh:357
double Pi
Float_t mchirp
Definition: injection.hh:55
TBranch * b_frequency
Definition: injection.hh:111
TBranch * b_mass
Definition: injection.hh:117
Int_t GetEntries()
Definition: injection.cc:183
Float_t * mass
detector specific effective distance
Definition: injection.hh:79
size_t mdc__IDSize()
Definition: network.hh:396
TBranch * b_phi
Definition: injection.hh:103
TBranch * b_rp0
Definition: injection.hh:96
TBranch * b_factor
Definition: injection.hh:93
Int_t run
number of detectors
Definition: injection.hh:47
TBranch * b_eventID
Definition: injection.hh:89
std::vector< wavearray< double > * > IWFP
Definition: detector.hh:361
wavearray< double > TDUR
Definition: detector.hh:358
double * entry
Definition: cwb_setcuts.C:206
std::vector< int > IWFID
Definition: detector.hh:360
Float_t * phi
source iota angle
Definition: injection.hh:65
Int_t type
Definition: injection.hh:50
wavearray< double > BAND
Definition: detector.hh:356
TBranch * b_run
Definition: injection.hh:87
Float_t * bx
beam pattern coefficients for hp
Definition: injection.hh:68
#define MDC_TREE_NAME
Definition: injection.cc:7
WSeries< double > * getTFmap()
param: no parameters
Definition: detector.hh:161
Float_t distance
Definition: injection.hh:54
double sTOp
char Name[16]
Definition: detector.hh:309
double fabs(const Complex &x)
Definition: numpy.cc:37
TBranch * b_strain
Definition: injection.hh:100
sprintf(tfres,"(1/%g)x(%g) (sec)x(Hz)", 2 *df, df)
TBranch * b_gps
Definition: injection.hh:99
Float_t * bp
source theta angle
Definition: injection.hh:67
DataType_t * data
Definition: wavearray.hh:301
TBranch * b_bx
Definition: injection.hh:106
wavearray< double > FREQ
Definition: detector.hh:355
Int_t eventID
Definition: injection.hh:49
bool save
double ctime
Float_t * frequency
estimated duration
Definition: injection.hh:73
Bool_t fill_in(network *, int, bool=true)
Definition: injection.cc:325
size_t getmdc__ID(size_t n)
Definition: network.hh:408
char fName[256]
TBranch * b_mchirp
Definition: injection.hh:95
wavearray< double > ** pwf
[x1,y1,z1,x2,y2,z2] components of spin vector
Definition: injection.hh:82
TBranch * b_ndim
pointer to the reconstructed waveform
Definition: injection.hh:86
Float_t * spin
[m1,m2], binary mass parameters
Definition: injection.hh:80
TTree * fChain
root input file cointainig the mdc TTree
Definition: injection.hh:41
wavearray< double > ISNR
Definition: detector.hh:354
detector ** pD
Int_t nevent
Definition: injection.hh:48
Float_t * Deff
injected snr in the detectors
Definition: injection.hh:78
exit(0)
TBranch * b_time
Definition: injection.hh:108