Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
ReadFileMDC.C
Go to the documentation of this file.
1 //
2 // This example show how to read event parameters from root mdc file
3 // The full parameter list is defined in trunk/wat/injection.hh
4 // Author : Gabriele Vedovato
5 
6 
7 void ReadFileMDC(TString fName, int nIFO=0) {
8 
9  injection event(fName.Data(),nIFO); // open wave root file
10  int nevt = event.GetEntries(); // get number of detected events
11  if(nevt==0) {
12  cout << "events are not present in the file : " << fName.Data() << endl;
13  gSystem->Exit(1);
14  } else cout << "nevt : " << nevt << endl;
15 
16  // if detector objects are present then print detector infos
17  TList* ifoList = event.fChain->GetUserInfo();
18  detectorParams dParams[NIFO_MAX];
19  for (int n=0;n<ifoList->GetSize();n++) {
20  detector* pDetector = (detector*)ifoList->At(n);
21  dParams[n] = pDetector->getDetectorParams();
22  cout << dParams[n].name << endl;
23  pDetector->print();
24  }
25 
26  // if nIFO is not provided (nIFO=0) then it is retrieved from the input root file
27  nIFO = nIFO==0 ? ifoList->GetSize() : nIFO;
28  if(nIFO==0) {
29  cout << "nIFO is not contained in the root file, must be declared as the input macro parameter" << endl;
30  gSystem->Exit(1);
31  }
32 
33  for(int n=0;n<nevt;n++) { // loop over the detected events
34 
35  event.GetEntry(n); // load event #n
36 
37  cout.precision(14);
38  cout << "-----------------------------------" << endl;
39  cout << " Event Injection Parameters " << n << endl;
40  cout << "-----------------------------------" << endl;
41  cout << endl;
42  cout << "mdc index : " << event.type << endl;
43  cout << "phi : " << event.phi[0] << endl;
44  cout << "theta : " << event.theta[0] << endl;
45  cout << endl;
46  for(int i=0;i<nIFO;i++) {
47  cout << "-----------------------------------" << endl;
48  if(ifoList->GetSize()==nIFO)
49  cout << "ifo : " << dParams[i].name << endl;
50  else
51  cout << "ifo : " << i << endl;
52  cout << "-----------------------------------" << endl;
53  cout << endl;
54  cout << "time : " << event.time[i] << endl;
55  cout << "hrss : " << event.hrss[i] << endl;
56  cout << endl;
57  }
58  }
59 
60  exit(0);
61 }
detectorParams getDetectorParams()
Definition: detector.cc:201
char name[32]
Definition: detector.hh:32
par[0] name
int n
Definition: cwb_net.C:10
TString("c")
int nevt
i drho i
#define nIFO
const int NIFO_MAX
Definition: wat.hh:4
void print()
Definition: detector.cc:1768
void ReadFileMDC(TString fName, int nIFO=0)
Definition: ReadFileMDC.C:7
char fName[256]
exit(0)