Fork me on GitHub

source: svn/trunk/src/HEPEVTConverter.cc@ 982

Last change on this file since 982 was 537, checked in by Xavier Rouby, 15 years ago

adding #include<cstdlib> for compatibility with newer versions of g++

File size: 8.9 KB
Line 
1/***********************************************************************
2** **
3** /----------------------------------------------\ **
4** | Delphes, a framework for the fast simulation | **
5** | of a generic collider experiment | **
6** \------------- arXiv:0903.2225v1 ------------/ **
7** **
8** **
9** This package uses: **
10** ------------------ **
11** ROOT: Nucl. Inst. & Meth. in Phys. Res. A389 (1997) 81-86 **
12** FastJet algorithm: Phys. Lett. B641 (2006) [hep-ph/0512210] **
13** Hector: JINST 2:P09005 (2007) [physics.acc-ph:0707.1198v2] **
14** FROG: [hep-ex/0901.2718v1] **
15** HepMC: Comput. Phys. Commun.134 (2001) 41 **
16** **
17** ------------------------------------------------------------------ **
18** **
19** Main authors: **
20** ------------- **
21** **
22** Severine Ovyn Xavier Rouby **
23** severine.ovyn@uclouvain.be xavier.rouby@cern **
24** **
25** Center for Particle Physics and Phenomenology (CP3) **
26** Universite catholique de Louvain (UCL) **
27** Louvain-la-Neuve, Belgium **
28** **
29** Copyright (C) 2008-2009, **
30** All rights reserved. **
31** **
32***********************************************************************/
33
34#include <iostream>
35#include <utility>
36#include <deque>
37#include <sstream>
38#include <fstream>
39#include <cmath>
40#include <cstdlib> // for exit()
41
42#include "TROOT.h"
43#include "TTree.h"
44#include "TBranch.h"
45#include "TLeaf.h"
46#include "TString.h"
47#include "TLorentzVector.h"
48#include "TChain.h"
49
50#include "BlockClasses.h"
51#include "ExRootTreeWriter.h"
52#include "ExRootTreeBranch.h"
53
54#include "HEPEVTConverter.h"
55
56using namespace std;
57
58//------------------------------------------------------------------------------
59
60HEPTreeReader::HEPTreeReader(TTree *tree, HEPEvent *event) : fChain(tree), fCurrentTree(-1), fEvent(event)
61{
62 if(!fChain) return;
63
64 TBranch *branch;
65 TLeaf *leaf;
66 TString name;
67 Int_t i;
68
69 name = "Nhep";
70 branch = fChain->GetBranch(name);
71 branch->SetAddress(&(event->Nhep));
72 fBranches.push_back(make_pair(name, branch));
73
74 TString intNames[3] = {"Idhep", "Jsmhep", "Jsdhep"};
75 Int_t **intData[3] = {&event->Idhep, &event->Jsmhep, &event->Jsdhep};
76
77 for(i = 0; i < 3; ++i)
78 {
79 name = intNames[i];
80 branch = fChain->GetBranch(name);
81 leaf = branch->GetLeaf(name);
82 *intData[i] = new Int_t[leaf->GetNdata()];
83 branch->SetAddress(*intData[i]);
84 fBranches.push_back(make_pair(name, branch));
85 }
86
87 TString floatNames[2] = {"Phep", "Vhep"};
88 Float_t **floatData[2] = {&event->Phep, &event->Vhep};
89
90 for(i = 0; i < 2; ++i)
91 {
92 name = floatNames[i];
93 branch = fChain->GetBranch(name);
94 leaf = branch->GetLeaf(name);
95 *floatData[i] = new Float_t[leaf->GetNdata()];
96 branch->SetAddress(*floatData[i]);
97 fBranches.push_back(make_pair(name, branch));
98 }
99}
100
101//------------------------------------------------------------------------------
102
103HEPTreeReader::~HEPTreeReader()
104{
105 if(!fChain) return;
106
107 delete[] fEvent->Idhep;
108 delete[] fEvent->Jsmhep;
109 delete[] fEvent->Jsdhep;
110 delete[] fEvent->Phep;
111 delete[] fEvent->Vhep;
112}
113
114//------------------------------------------------------------------------------
115
116Bool_t HEPTreeReader::ReadEntry(Long64_t element)
117{
118 if(!fChain) return kFALSE;
119
120 Int_t treeEntry = fChain->LoadTree(element);
121 if(treeEntry < 0) return kFALSE;
122
123 if(fChain->IsA() == TChain::Class())
124 {
125 TChain *chain = static_cast<TChain*>(fChain);
126 if(chain->GetTreeNumber() != fCurrentTree)
127 {
128 fCurrentTree = chain->GetTreeNumber();
129 Notify();
130 }
131 }
132
133 deque< pair<TString, TBranch*> >::iterator it_deque;
134 TBranch *branch;
135
136 for(it_deque = fBranches.begin(); it_deque != fBranches.end(); ++it_deque)
137 {
138 branch = it_deque->second;
139 if(branch)
140 {
141 branch->GetEntry(treeEntry);
142 }
143 }
144
145 return kTRUE;
146}
147
148//------------------------------------------------------------------------------
149
150void HEPTreeReader::Notify()
151{
152 // Called when loading a new file.
153 // Get branch pointers.
154 if(!fChain) return;
155
156 deque< pair<TString, TBranch*> >::iterator it_deque;
157
158 for(it_deque = fBranches.begin(); it_deque != fBranches.end(); ++it_deque)
159 {
160 it_deque->second = fChain->GetBranch(it_deque->first);
161 if(!it_deque->second)
162 {
163 cerr << left << setw(40) <<"** WARNING: cannot get branch: "<<""
164 << left << setw(20) << it_deque->first <<""
165 << right << setw(9) <<" **"<<""<<endl;
166 }
167 }
168}
169
170
171
172HEPEVTConverter::HEPEVTConverter(const string& inputFileList, const string& outputFileName, const PdgTable& pdg, const int& Nevents) : DataConverter(pdg,Nevents)
173{
174 string buffer;
175 int nevt_already_processed=0;
176
177 gROOT->SetBatch();
178
179 ExRootTreeWriter *treeWriter = new ExRootTreeWriter(outputFileName, "GEN");
180 ExRootTreeBranch *branchGen = treeWriter->NewBranch("Particle", TRootC::GenParticle::Class());
181
182 ifstream infile(inputFileList.c_str());
183 if(!infile.is_open())
184 {
185 cerr << left << setw(30) <<"** ERROR: Can't open "<<""
186 << left << setw(20) << inputFileList <<""
187 << right << setw(19) <<"for input **"<<""<<endl;
188 exit(1);
189 }
190 while(1)
191 {
192 if (Nevt>0 && nevt_already_processed >=Nevt) break; // enough events already processed
193 TChain *chain = new TChain("h101");
194 infile >> buffer;
195 ifstream checking_the_file(buffer.c_str());
196 if(!checking_the_file.good())
197 {
198 cerr << left << setw(30) <<"** ERROR: Can't find file "<<""
199 << left << setw(36) << buffer <<""
200 << right << setw(3) <<" **"<<""<<endl;
201 continue;
202 }
203 else checking_the_file.close();
204
205 if(!infile.good()) break; // end of the list
206 chain->Add(buffer.c_str());
207
208 HEPEvent event;
209 HEPTreeReader *treeReader = new HEPTreeReader(chain, &event);
210
211 Long64_t entry, allEntries = treeReader->GetEntries();
212 Int_t address, particle;
213 Float_t signEta;
214
215 TRootC::GenParticle *element;
216 allEntries = (Nevt<0)?allEntries:min((int)allEntries,Nevt);
217 // Loop over all events
218 for(entry = 0; entry < allEntries; ++entry)
219 {
220 if (Nevt>0 && nevt_already_processed >=Nevt) break; // enough events already processed
221 // Load selected branches with data from specified event
222 treeReader->ReadEntry(entry);
223 treeWriter->Clear();
224
225 for(particle = 0; particle < event.Nhep; ++particle)
226 {
227 element = (TRootC::GenParticle*) branchGen->NewEntry();
228
229 element->PID = event.Idhep[particle];
230 element->Status = event.Jsmhep[particle]/16000000
231 + event.Jsdhep[particle]/16000000*100;
232 element->M1 = (event.Jsmhep[particle]%16000000)/4000 -1;
233 element->M2 = event.Jsmhep[particle]%4000 -1;
234 element->D1 = (event.Jsdhep[particle]%16000000)/4000 -1;
235 element->D2 = event.Jsdhep[particle]%4000 -1;
236 address = particle*5;
237 element->E = event.Phep[address + 3];
238 element->Px = event.Phep[address + 0];
239 element->Py = event.Phep[address + 1];
240 element->Pz = event.Phep[address + 2];
241 element->M = event.Phep[address + 4];
242
243 TLorentzVector vector(element->Px, element->Py, element->Pz, element->E);
244 element->PT = vector.Perp();
245 signEta = (element->Pz >= 0.0) ? 1.0 : -1.0;
246 element->Eta = vector.CosTheta()*vector.CosTheta() == 1.0 ? signEta*999.9 : vector.Eta();
247 element->Phi = vector.Phi();
248
249 address = particle*4;
250 element->T = event.Vhep[address + 3];
251 element->X = event.Vhep[address + 0];
252 element->Y = event.Vhep[address + 1];
253 element->Z = event.Vhep[address + 2];
254 }
255 treeWriter->Fill();
256 ++nevt_already_processed;
257 }
258 delete chain;
259 delete treeReader;
260 }
261 treeWriter->Write();
262
263 delete treeWriter;
264
265}
266
267
268
269
Note: See TracBrowser for help on using the repository browser.