Fork me on GitHub

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

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

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

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