Fork me on GitHub

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

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

minor change: cleaning of the 'cout' and 'cerr'

File size: 8.8 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
41#include "TROOT.h"
42#include "TTree.h"
43#include "TBranch.h"
44#include "TLeaf.h"
45#include "TString.h"
46#include "TLorentzVector.h"
47#include "TChain.h"
48
49#include "BlockClasses.h"
50#include "ExRootTreeWriter.h"
51#include "ExRootTreeBranch.h"
52
53#include "HEPEVTConverter.h"
54
55using namespace std;
56
57//------------------------------------------------------------------------------
58
59HEPTreeReader::HEPTreeReader(TTree *tree, HEPEvent *event) : fChain(tree), fCurrentTree(-1), fEvent(event)
60{
61 if(!fChain) return;
62
63 TBranch *branch;
64 TLeaf *leaf;
65 TString name;
66 Int_t i;
67
68 name = "Nhep";
69 branch = fChain->GetBranch(name);
70 branch->SetAddress(&(event->Nhep));
71 fBranches.push_back(make_pair(name, branch));
72
73 TString intNames[3] = {"Idhep", "Jsmhep", "Jsdhep"};
74 Int_t **intData[3] = {&event->Idhep, &event->Jsmhep, &event->Jsdhep};
75
76 for(i = 0; i < 3; ++i)
77 {
78 name = intNames[i];
79 branch = fChain->GetBranch(name);
80 leaf = branch->GetLeaf(name);
81 *intData[i] = new Int_t[leaf->GetNdata()];
82 branch->SetAddress(*intData[i]);
83 fBranches.push_back(make_pair(name, branch));
84 }
85
86 TString floatNames[2] = {"Phep", "Vhep"};
87 Float_t **floatData[2] = {&event->Phep, &event->Vhep};
88
89 for(i = 0; i < 2; ++i)
90 {
91 name = floatNames[i];
92 branch = fChain->GetBranch(name);
93 leaf = branch->GetLeaf(name);
94 *floatData[i] = new Float_t[leaf->GetNdata()];
95 branch->SetAddress(*floatData[i]);
96 fBranches.push_back(make_pair(name, branch));
97 }
98}
99
100//------------------------------------------------------------------------------
101
102HEPTreeReader::~HEPTreeReader()
103{
104 if(!fChain) return;
105
106 delete[] fEvent->Idhep;
107 delete[] fEvent->Jsmhep;
108 delete[] fEvent->Jsdhep;
109 delete[] fEvent->Phep;
110 delete[] fEvent->Vhep;
111}
112
113//------------------------------------------------------------------------------
114
115Bool_t HEPTreeReader::ReadEntry(Long64_t element)
116{
117 if(!fChain) return kFALSE;
118
119 Int_t treeEntry = fChain->LoadTree(element);
120 if(treeEntry < 0) return kFALSE;
121
122 if(fChain->IsA() == TChain::Class())
123 {
124 TChain *chain = static_cast<TChain*>(fChain);
125 if(chain->GetTreeNumber() != fCurrentTree)
126 {
127 fCurrentTree = chain->GetTreeNumber();
128 Notify();
129 }
130 }
131
132 deque< pair<TString, TBranch*> >::iterator it_deque;
133 TBranch *branch;
134
135 for(it_deque = fBranches.begin(); it_deque != fBranches.end(); ++it_deque)
136 {
137 branch = it_deque->second;
138 if(branch)
139 {
140 branch->GetEntry(treeEntry);
141 }
142 }
143
144 return kTRUE;
145}
146
147//------------------------------------------------------------------------------
148
149void HEPTreeReader::Notify()
150{
151 // Called when loading a new file.
152 // Get branch pointers.
153 if(!fChain) return;
154
155 deque< pair<TString, TBranch*> >::iterator it_deque;
156
157 for(it_deque = fBranches.begin(); it_deque != fBranches.end(); ++it_deque)
158 {
159 it_deque->second = fChain->GetBranch(it_deque->first);
160 if(!it_deque->second)
161 {
162 cerr << left << setw(40) <<"** WARNING: cannot get branch: "<<""
163 << left << setw(20) << it_deque->first <<""
164 << right << setw(9) <<" **"<<""<<endl;
165 }
166 }
167}
168
169
170
171HEPEVTConverter::HEPEVTConverter(const string& inputFileList, const string& outputFileName, const PdgTable& pdg, const int& Nevents) : DataConverter(pdg,Nevents)
172{
173 string buffer;
174 int nevt_already_processed=0;
175
176 gROOT->SetBatch();
177
178 ExRootTreeWriter *treeWriter = new ExRootTreeWriter(outputFileName, "GEN");
179 ExRootTreeBranch *branchGen = treeWriter->NewBranch("Particle", TRootC::GenParticle::Class());
180
181 ifstream infile(inputFileList.c_str());
182 if(!infile.is_open())
183 {
184 cerr << left << setw(30) <<"** ERROR: Can't open "<<""
185 << left << setw(20) << inputFileList <<""
186 << right << setw(19) <<"for input **"<<""<<endl;
187 exit(1);
188 }
189 while(1)
190 {
191 if (Nevt>0 && nevt_already_processed >=Nevt) break; // enough events already processed
192 TChain *chain = new TChain("h101");
193 infile >> buffer;
194 ifstream checking_the_file(buffer.c_str());
195 if(!checking_the_file.good())
196 {
197 cerr << left << setw(30) <<"** ERROR: Can't find file "<<""
198 << left << setw(36) << buffer <<""
199 << right << setw(3) <<" **"<<""<<endl;
200 continue;
201 }
202 else checking_the_file.close();
203
204 if(!infile.good()) break; // end of the list
205 chain->Add(buffer.c_str());
206
207 HEPEvent event;
208 HEPTreeReader *treeReader = new HEPTreeReader(chain, &event);
209
210 Long64_t entry, allEntries = treeReader->GetEntries();
211 Int_t address, particle;
212 Float_t signEta;
213
214 TRootC::GenParticle *element;
215 allEntries = (Nevt<0)?allEntries:min((int)allEntries,Nevt);
216 // Loop over all events
217 for(entry = 0; entry < allEntries; ++entry)
218 {
219 if (Nevt>0 && nevt_already_processed >=Nevt) break; // enough events already processed
220 // Load selected branches with data from specified event
221 treeReader->ReadEntry(entry);
222 treeWriter->Clear();
223
224 for(particle = 0; particle < event.Nhep; ++particle)
225 {
226 element = (TRootC::GenParticle*) branchGen->NewEntry();
227
228 element->PID = event.Idhep[particle];
229 element->Status = event.Jsmhep[particle]/16000000
230 + event.Jsdhep[particle]/16000000*100;
231 element->M1 = (event.Jsmhep[particle]%16000000)/4000 -1;
232 element->M2 = event.Jsmhep[particle]%4000 -1;
233 element->D1 = (event.Jsdhep[particle]%16000000)/4000 -1;
234 element->D2 = event.Jsdhep[particle]%4000 -1;
235 address = particle*5;
236 element->E = event.Phep[address + 3];
237 element->Px = event.Phep[address + 0];
238 element->Py = event.Phep[address + 1];
239 element->Pz = event.Phep[address + 2];
240 element->M = event.Phep[address + 4];
241
242 TLorentzVector vector(element->Px, element->Py, element->Pz, element->E);
243 element->PT = vector.Perp();
244 signEta = (element->Pz >= 0.0) ? 1.0 : -1.0;
245 element->Eta = vector.CosTheta()*vector.CosTheta() == 1.0 ? signEta*999.9 : vector.Eta();
246 element->Phi = vector.Phi();
247
248 address = particle*4;
249 element->T = event.Vhep[address + 3];
250 element->X = event.Vhep[address + 0];
251 element->Y = event.Vhep[address + 1];
252 element->Z = event.Vhep[address + 2];
253 }
254 treeWriter->Fill();
255 ++nevt_already_processed;
256 }
257 delete chain;
258 delete treeReader;
259 }
260 treeWriter->Write();
261
262 delete treeWriter;
263
264}
265
266
267
268
Note: See TracBrowser for help on using the repository browser.