source: trunk/test/ExRootHEPEVTConverter.cpp@ 11

Last change on this file since 11 was 2, checked in by Pavel Demin, 16 years ago

first commit

File size: 7.6 KB
Line 
1
2#include <iostream>
3#include <utility>
4#include <deque>
5
6#include "TROOT.h"
7#include "TApplication.h"
8
9#include "TFile.h"
10#include "TTree.h"
11#include "TBranch.h"
12#include "TLeaf.h"
13#include "TString.h"
14#include "TLorentzVector.h"
15
16#include "ExRootAnalysis/ExRootClasses.h"
17
18#include "ExRootAnalysis/ExRootTreeReader.h"
19#include "ExRootAnalysis/ExRootTreeWriter.h"
20#include "ExRootAnalysis/ExRootTreeBranch.h"
21#include "ExRootAnalysis/ExRootUtilities.h"
22
23using namespace std;
24
25//------------------------------------------------------------------------------
26
27struct HEPEvent
28{
29 Int_t Nevhep;
30 Int_t Nhep;
31 Int_t *Idhep; //[Nhep]
32 Int_t *Jsmhep; //[Nhep]
33 Int_t *Jsdhep; //[Nhep]
34 Float_t *Phep; //[Nhep][5]
35 Float_t *Vhep; //[Nhep][4]
36 Int_t Irun;
37 Int_t Ievt;
38 Float_t Weight;
39 Float_t Xsecn;
40 Int_t Ifilter;
41 Int_t Nparam;
42 Float_t *Param; //[Nparam]
43};
44
45//------------------------------------------------------------------------------
46
47class HEPTreeReader
48{
49public:
50
51 HEPTreeReader(TTree *tree, HEPEvent *event);
52 ~HEPTreeReader();
53
54 Long64_t GetEntries() const { return fChain ? static_cast<Long64_t>(fChain->GetEntries()) : 0; }
55
56 Bool_t ReadEntry(Long64_t element);
57
58private:
59
60 void Notify();
61
62 TTree *fChain; // pointer to the analyzed TTree or TChain
63 Int_t fCurrentTree; // current Tree number in a TChain
64
65 HEPEvent *fEvent;
66
67 deque< pair<TString, TBranch*> > fBranches;
68
69};
70
71//------------------------------------------------------------------------------
72
73HEPTreeReader::HEPTreeReader(TTree *tree, HEPEvent *event) : fChain(tree), fCurrentTree(-1), fEvent(event)
74{
75 if(!fChain) return;
76
77 TBranch *branch;
78 TLeaf *leaf;
79 TString name;
80 Int_t i;
81
82 name = "Nhep";
83 branch = fChain->GetBranch(name);
84 branch->SetAddress(&(event->Nhep));
85 fBranches.push_back(make_pair(name, branch));
86
87 TString intNames[3] = {"Idhep", "Jsmhep", "Jsdhep"};
88 Int_t **intData[3] = {&event->Idhep, &event->Jsmhep, &event->Jsdhep};
89
90 for(i = 0; i < 3; ++i)
91 {
92 name = intNames[i];
93 branch = fChain->GetBranch(name);
94 leaf = branch->GetLeaf(name);
95 *intData[i] = new Int_t[leaf->GetNdata()];
96 branch->SetAddress(*intData[i]);
97 fBranches.push_back(make_pair(name, branch));
98 }
99
100 TString floatNames[2] = {"Phep", "Vhep"};
101 Float_t **floatData[2] = {&event->Phep, &event->Vhep};
102
103 for(i = 0; i < 2; ++i)
104 {
105 name = floatNames[i];
106 branch = fChain->GetBranch(name);
107 leaf = branch->GetLeaf(name);
108 *floatData[i] = new Float_t[leaf->GetNdata()];
109 branch->SetAddress(*floatData[i]);
110 fBranches.push_back(make_pair(name, branch));
111 }
112}
113
114//------------------------------------------------------------------------------
115
116HEPTreeReader::~HEPTreeReader()
117{
118 if(!fChain) return;
119
120 delete[] fEvent->Idhep;
121 delete[] fEvent->Jsmhep;
122 delete[] fEvent->Jsdhep;
123 delete[] fEvent->Phep;
124 delete[] fEvent->Vhep;
125}
126
127//------------------------------------------------------------------------------
128
129Bool_t HEPTreeReader::ReadEntry(Long64_t element)
130{
131 if(!fChain) return kFALSE;
132
133 Int_t treeEntry = fChain->LoadTree(element);
134 if(treeEntry < 0) return kFALSE;
135
136 if(fChain->IsA() == TChain::Class())
137 {
138 TChain *chain = static_cast<TChain*>(fChain);
139 if(chain->GetTreeNumber() != fCurrentTree)
140 {
141 fCurrentTree = chain->GetTreeNumber();
142 Notify();
143 }
144 }
145
146 deque< pair<TString, TBranch*> >::iterator it_deque;
147 TBranch *branch;
148
149 for(it_deque = fBranches.begin(); it_deque != fBranches.end(); ++it_deque)
150 {
151 branch = it_deque->second;
152 if(branch)
153 {
154 branch->GetEntry(treeEntry);
155 }
156 }
157
158 return kTRUE;
159}
160
161//------------------------------------------------------------------------------
162
163void HEPTreeReader::Notify()
164{
165 // Called when loading a new file.
166 // Get branch pointers.
167 if(!fChain) return;
168
169 deque< pair<TString, TBranch*> >::iterator it_deque;
170
171 for(it_deque = fBranches.begin(); it_deque != fBranches.end(); ++it_deque)
172 {
173 it_deque->second = fChain->GetBranch(it_deque->first);
174 if(!it_deque->second)
175 {
176 cout << "** WARNING: cannot get branch '" << it_deque->first << "'" << endl;
177 }
178 }
179}
180
181
182//------------------------------------------------------------------------------
183
184int main(int argc, char *argv[])
185{
186 char *appName = "ExRootHEPEVTConverter";
187
188 if(argc != 3)
189 {
190 cout << " Usage: " << appName << " input_file" << " output_file" << endl;
191 cout << " input_file - list of HEPEVT files in ROOT format ('h101' tree)," << endl;
192 cout << " output_file - output file in ROOT format." << endl;
193 return 1;
194 }
195
196 gROOT->SetBatch();
197
198 int appargc = 1;
199 char *appargv[] = {appName};
200 TApplication app(appName, &appargc, appargv);
201
202 TString inputFileList(argv[1]);
203 TString outputFileName(argv[2]);
204 string buffer;
205
206 TChain *chain = new TChain("h101");
207 if(FillChain(chain, inputFileList))
208 {
209
210 HEPEvent event;
211 HEPTreeReader *treeReader = new HEPTreeReader(chain, &event);
212
213 TFile *outputFile = TFile::Open(outputFileName, "RECREATE");
214 ExRootTreeWriter *treeWriter = new ExRootTreeWriter(outputFile, "Analysis");
215
216 ExRootTreeBranch *branchGen = treeWriter->NewBranch("Gen", ExRootGenParticle::Class());
217
218 Long64_t entry, allEntries = treeReader->GetEntries();
219
220 cout << "** Chain contains " << allEntries << " events" << endl;
221
222 Int_t address, particle;
223 Double_t signPz;
224
225 TLorentzVector momentum;
226
227 ExRootGenParticle *element;
228
229 // Loop over all events
230 for(entry = 0; entry < allEntries; ++entry)
231 {
232 // Load selected branches with data from specified event
233 treeReader->ReadEntry(entry);
234 treeWriter->Clear();
235
236 if((entry % 100) == 0 && entry > 0 )
237 {
238 cout << "** Processing element # " << entry << endl;
239 }
240
241 for(particle = 0; particle < event.Nhep; ++particle)
242 {
243 element = (ExRootGenParticle*) branchGen->NewEntry();
244
245 element->PID = event.Idhep[particle];
246 element->Status = event.Jsmhep[particle]/16000000
247 + event.Jsdhep[particle]/16000000*100;
248 element->M1 = (event.Jsmhep[particle]%16000000)/4000 - 1;
249 element->M2 = event.Jsmhep[particle]%4000 - 1;
250 element->D1 = (event.Jsdhep[particle]%16000000)/4000 - 1;
251 element->D2 = event.Jsdhep[particle]%4000 - 1;
252
253 address = particle*5;
254 element->E = event.Phep[address + 3];
255 element->Px = event.Phep[address + 0];
256 element->Py = event.Phep[address + 1];
257 element->Pz = event.Phep[address + 2];
258
259 momentum.SetPxPyPzE(element->Px, element->Py, element->Pz, element->E);
260 element->PT = momentum.Perp();
261 signPz = (element->Pz >= 0.0) ? 1.0 : -1.0;
262 element->Eta = element->PT == 0.0 ? signPz*999.9 : momentum.Eta();
263 element->Phi = momentum.Phi();
264
265 element->Rapidity = element->PT == 0.0 ? signPz*999.9 : momentum.Rapidity();
266
267 address = particle*4;
268 element->T = event.Vhep[address + 3];
269 element->X = event.Vhep[address + 0];
270 element->Y = event.Vhep[address + 1];
271 element->Z = event.Vhep[address + 2];
272 }
273 treeWriter->Fill();
274 }
275
276 treeWriter->Write();
277
278 cout << "** Exiting..." << endl;
279
280 delete treeWriter;
281 delete outputFile;
282 delete treeReader;
283 }
284 delete chain;
285}
286
287
288
Note: See TracBrowser for help on using the repository browser.