Fork me on GitHub

source: svn/trunk/src/HepMCConverter.cc@ 376

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

bug fix: the gen-level charge was wrong in GEN tree

File size: 10.9 KB
Line 
1/***********************************************************************
2** **
3** /----------------------------------------------\ **
4** | Delphes, a framework for the fast simulation | **
5** | of a generic collider experiment | **
6** \----------------------------------------------/ **
7** **
8** **
9** This package uses: **
10** ------------------ **
11** FastJet algorithm: Phys. Lett. B641 (2006) [hep-ph/0512210] **
12** Hector: JINST 2:P09005 (2007) [physics.acc-ph:0707.1198v2] **
13** FROG: [hep-ex/0901.2718v1] **
14** **
15** ------------------------------------------------------------------ **
16** **
17** Main authors: **
18** ------------- **
19** **
20** Severine Ovyn Xavier Rouby **
21** severine.ovyn@uclouvain.be xavier.rouby@cern **
22** **
23** Center for Particle Physics and Phenomenology (CP3) **
24** Universite catholique de Louvain (UCL) **
25** Louvain-la-Neuve, Belgium **
26** **
27** Copyright (C) 2008-2009, **
28** All rights reserved. **
29** **
30***********************************************************************/
31
32
33#include <iostream>
34#include <fstream>
35#include "TLorentzVector.h"
36#include "BlockClasses.h"
37
38#include "SmearUtil.h"
39#include "ExRootTreeWriter.h"
40#include "ExRootTreeBranch.h"
41#include "HepMCConverter.h"
42
43#include "GenParticle.h"
44#include "GenVertex.h"
45#include "IO_Ascii.h"
46#include "IO_GenEvent.h"
47
48//-------------------------------------------------------------------------
49int HepMCConverter::find_in_map( const std::map<HepMC::GenParticle*,int>& m, HepMC::GenParticle *p) const
50{
51 std::map<HepMC::GenParticle*,int>::const_iterator iter = m.find(p);
52 return (iter == m.end()) ? 0 : iter->second;
53}
54
55//--------------------------------------------------------------------------
56void HepMCConverter::ReadStats()
57{
58
59 unsigned int particle_counter=0;
60 index_to_particle.reserve(evt->particles_size());
61 index_to_particle[0] = 0;
62 HepMC::GenEvent::vertex_const_iterator v;
63 for (v = evt->vertices_begin(); v != evt->vertices_end(); ++v )
64 {
65 // making a list of incoming particles of the vertices
66 // so that the mother indices in HEPEVT can be filled properly
67 HepMC::GenVertex::particles_out_const_iterator p1;
68 for (p1 = (*v)->particles_in_const_begin();p1 != (*v)->particles_in_const_end(); ++p1 )
69 {
70
71 ++particle_counter;
72 //particle_counter can be very large for heavy ions
73 if(particle_counter >= index_to_particle.size() )
74 {
75 //make it large enough to hold up to this index
76 index_to_particle.resize(particle_counter+1);
77 }
78 index_to_particle[particle_counter] = *p1;
79 particle_to_index[*p1] = particle_counter;
80 }
81 // daughters are entered only if they aren't a mother of
82 // another vertex
83 HepMC::GenVertex::particles_out_const_iterator p2;
84 for (p2 = (*v)->particles_out_const_begin();p2 != (*v)->particles_out_const_end(); ++p2)
85 {
86 if (!(*p2)->end_vertex())
87 {
88 ++particle_counter;
89 //particle_counter can be very large for heavy ions
90 if(particle_counter >= index_to_particle.size() )
91 {
92 //make it large enough to hold up to this index
93 index_to_particle.resize(particle_counter+1);
94 }
95 index_to_particle[particle_counter] = *p2;
96 particle_to_index[*p2] = particle_counter;
97 }
98 }
99
100 }
101}
102
103
104//-------------------------------------------------------------------------
105void HepMCConverter::getStatsFromTuple(int &mo1, int &mo2, int &da1, int &da2, int &status, int &pid, int j) const
106{
107 if (!evt)
108 {
109 cout << "HepMCFileReader: Got no event :-( Game over already ?" <<endl;
110 }
111 else
112 {
113 status = index_to_particle[j]->status();
114 pid = index_to_particle[j]->pdg_id();
115 if ( index_to_particle[j]->production_vertex() )
116 {
117 int num_mothers = index_to_particle[j]->production_vertex()->particles_in_size();
118 if (num_mothers ==0) {
119 mo1 = 0;
120 mo2 = 0;
121 }
122 else {
123 int first_mother = find_in_map( particle_to_index,*(index_to_particle[j]->production_vertex()->particles_in_const_begin()));
124 int last_mother = first_mother + num_mothers - 1;
125 if ( first_mother == 0 ) last_mother = 0;
126 mo1=first_mother;
127 mo2=last_mother;
128 } // if num_mothers !=0
129 }
130 else // no data on production_vertex
131 {
132 mo1 =0;
133 mo2 =0;
134 }
135 if (index_to_particle[j]->end_vertex())
136 {
137 //find # of 1. daughter
138 int first_daughter = find_in_map( particle_to_index,*(index_to_particle[j]->end_vertex()->particles_begin(HepMC::children)));
139 //cout <<"first_daughter "<< first_daughter << "num_daughters " << num_daughters << endl;
140 HepMC::GenVertex::particle_iterator ic;
141 int last_daughter=0;
142 //find # of last daughter
143 for (ic = index_to_particle[j]->end_vertex()->particles_begin(HepMC::children);ic != index_to_particle[j]->end_vertex()->particles_end(HepMC::children); ++ic)
144 last_daughter= find_in_map( particle_to_index,*ic);
145
146 if (first_daughter== 0) last_daughter = 0;
147 da1=first_daughter;
148 da2=last_daughter;
149 }
150 else
151 {
152 da1=0;
153 da2=0;
154 }
155 }
156}
157
158
159
160//---------------------------------------------------------------------------
161
162void HepMCConverter::AnalyseEvent(ExRootTreeBranch *branch,HepMC::GenEvent& evt,const Long64_t eventNumber)
163{
164 TRootLHEFEvent *element;
165
166 element = static_cast<TRootLHEFEvent*>(branch->NewEntry());
167 element->Number = eventNumber;
168 element->Nparticles = evt.particles_size();
169 element->ProcessID = evt.signal_process_id();
170 // element->Weight = hepeup.XWGTUP;
171 element->ScalePDF = evt.event_scale();
172 element->CouplingQED = evt.alphaQED();
173 element->CouplingQCD = evt.alphaQCD();
174
175}
176
177//---------------------------------------------------------------------------
178
179void HepMCConverter::AnalyseParticles(ExRootTreeBranch *branch, const HepMC::GenEvent& evt)
180{
181 TRootC::GenParticle *element;
182
183 TLorentzVector momentum;
184 Double_t signPz;
185
186 ReadStats();
187 for(int n=1; n<=evt.particles_size(); n++)
188 {
189 getStatsFromTuple( mo1,mo2,da1,da2,status,pid,n);
190
191 element = static_cast<TRootC::GenParticle*>(branch->NewEntry());
192
193 element->PID = pid;
194 element->Status = status;
195 element->M1 = mo1;
196 element->M2 = mo2;
197 element->D1 = da1;
198 element->D2 = da2;
199 element->Charge = ChargeVal(pid);
200
201 element->E = index_to_particle[n]->momentum().e();
202 element->Px = index_to_particle[n]->momentum().px();
203 element->Py = index_to_particle[n]->momentum().py();
204 element->Pz = index_to_particle[n]->momentum().pz();
205 element->M = index_to_particle[n]->momentum().m();
206
207 float PT = sqrt(pow(element->Px,2)+pow(element->Py,2));
208 element->PT = PT;
209
210 momentum.SetPxPyPzE(element->Px, element->Py, element->Pz, element->E);
211 signPz = (element->Pz >= 0.0) ? 1.0 : -1.0;
212 //element->Eta = element->PT == 0.0 ? signPz*999.9 : momentum.Eta(); to avoid a warning from ROOT, replace the "==0" by "< 1e-6"
213 element->Eta = element->PT < 1e-6 ? signPz*999.9 : momentum.Eta();
214 element->Phi = index_to_particle[n]->momentum().phi();
215
216 //In particle at vertex
217 HepMC::GenVertex* vrtI = (index_to_particle[n])->production_vertex();
218 HepMC::GenVertex::particles_in_const_iterator partI;
219
220 if(vrtI)
221 {
222 element->T = vrtI->position().t();
223 element->X = vrtI->position().x();
224 element->Y = vrtI->position().y();
225 element->Z = vrtI->position().z();
226 }
227 else
228 {
229 element->T = 0.;
230 element->X = 0.;
231 element->Y = 0.;
232 element->Z = 0.;
233 }
234 }
235}
236
237//
238
239//------------------------------------------------------------------------------
240
241HepMCConverter::~HepMCConverter()
242{
243 delete evt;
244
245 /*cout << "delete index to particle" << endl;
246 std::vector<HepMC::GenParticle*>::iterator i;
247 / for (i=index_to_particle.begin();i != index_to_particle.end(); i++) delete (*i);
248 cout << "done" << endl;
249
250 std::map<HepMC::GenParticle*,int>::iterator j;
251 for (j=particle_to_index.begin(); j != particle_to_index.end(); j++) delete j->first;
252 */
253
254}
255
256//------------------------------------------------------------------------------
257
258HepMCConverter::HepMCConverter(const string& inputFileList, const string& outputFileName, const int& Nevents) : DataConverter(Nevents)
259{
260
261 ExRootTreeWriter *treeWriter = new ExRootTreeWriter(outputFileName, "GEN");
262 // information about generated event
263 ExRootTreeBranch *branchGenEvent = treeWriter->NewBranch("Event", TRootLHEFEvent::Class());
264 // generated particles from HEPEVT
265 ExRootTreeBranch *branchGenParticle = treeWriter->NewBranch("Particle", TRootC::GenParticle::Class());
266
267 // Open a stream connected to an event file:
268 ifstream infile(inputFileList.c_str());
269 string filename;
270 if(!infile.is_open()) {
271 cerr << left << setw(30) <<"** ERROR: Can't open "<<""
272 << left << setw(20) << inputFileList <<""
273 << right << setw(19) <<"for input **"<<""<<endl;
274 exit(1);
275 }
276
277 Long64_t entry = 0;
278 while(1)
279 {
280 infile >> filename;
281 if(!infile.good()) break;
282 ifstream checking_the_file(filename.c_str());
283 if(!checking_the_file.good())
284 {
285 cerr << left << setw(30) <<"** ERROR: Can't find file "<<""
286 << left << setw(20) << filename <<""
287 << right << setw(19) <<"for input **"<<""<<endl;
288 continue;
289 }
290 else checking_the_file.close();
291
292 //Open the file
293 HepMC::IO_GenEvent ascii_in(filename,std::ios::in);
294 // get the first event
295 evt = ascii_in.read_next_event();
296
297 while ( evt )
298 {
299 treeWriter->Clear();
300 AnalyseEvent(branchGenEvent, *evt,entry+1);
301 AnalyseParticles(branchGenParticle, *evt);
302 delete evt;
303 // read the next event
304 ascii_in >> evt;
305 treeWriter->Fill();
306 ++entry;
307 }
308 }
309
310 treeWriter->Write();
311 delete treeWriter;
312
313}
314
315
316
317
Note: See TracBrowser for help on using the repository browser.