Fork me on GitHub

source: svn/trunk/src/LHEFConverter.cc@ 441

Last change on this file since 441 was 419, checked in by Xavier Rouby, 16 years ago

Nevent implemented in DetectorCard

File size: 6.8 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
36#include "TLorentzVector.h"
37
38#include "ExRootTreeWriter.h"
39#include "BlockClasses.h"
40#include "LHEFConverter.h"
41#include "SmearUtil.h"
42#include "PdgParticle.h"
43
44using namespace std;
45
46
47//------------------------------------------------------------------------------
48
49void LHEFConverter::AnalyseEvent(LHEF::Reader *reader, ExRootTreeBranch *branch, const Long64_t eventNumber)
50{
51 const LHEF::HEPEUP &hepeup = reader->hepeup;
52
53 TRootLHEFEvent *element;
54
55 element = (TRootLHEFEvent*) branch->NewEntry();
56
57 element->Number = eventNumber;
58 element->Nparticles = hepeup.NUP;
59 element->ProcessID = hepeup.IDPRUP;
60 element->Weight = hepeup.XWGTUP;
61 element->ScalePDF = hepeup.SCALUP;
62 element->CouplingQED = hepeup.AQEDUP;
63 element->CouplingQCD = hepeup.AQCDUP;
64}
65
66//---------------------------------------------------------------------------
67
68void LHEFConverter::AnalyseParticles(LHEF::Reader *reader, ExRootTreeBranch *branch)
69{
70 const LHEF::HEPEUP &hepeup = reader->hepeup;
71 Double_t signPz;
72 TLorentzVector momentum;
73 TRootLHEFParticle *element;
74
75 for(Int_t particle = 0; particle < hepeup.NUP; ++particle)
76 {
77 element = (TRootLHEFParticle*) branch->NewEntry();
78
79 element->PID = hepeup.IDUP[particle];
80 element->Status = hepeup.ISTUP[particle];
81 element->Mother1 = hepeup.MOTHUP[particle].first;
82 element->Mother2 = hepeup.MOTHUP[particle].second;
83 element->ColorLine1 = hepeup.ICOLUP[particle].first;
84 element->ColorLine2 = hepeup.ICOLUP[particle].second;
85 element->Px = hepeup.PUP[particle][0];
86 element->Py = hepeup.PUP[particle][1];
87 element->Pz = hepeup.PUP[particle][2];
88 element->E = hepeup.PUP[particle][3];
89 //element->M = hepeup.PUP[particle][4];
90 //element->Charge = ChargeVal(element->PID);
91
92 PdgParticle pdg_part = pdg_table[element->PID];
93 element->Charge = pdg_part.charge();
94 element->M = pdg_part.mass();
95
96 momentum.SetPxPyPzE(element->Px, element->Py, element->Pz, element->E);
97 element->PT = momentum.Perp();
98 signPz = (element->Pz >= 0.0) ? 1.0 : -1.0;
99 //element->Eta = element->PT == 0.0 ? signPz*999.9 : momentum.Eta(); to avoid a warning from ROOT, replace the "==0" by "< 1e-6"
100 element->Eta = element->PT < 1E-6 ? signPz*999.9 : momentum.Eta();
101 element->Phi = momentum.Phi();
102 element->Rapidity = element->PT < 1E-6 ? signPz*999.9 : momentum.Rapidity();
103 element->LifeTime = hepeup.VTIMUP[particle];
104 element->Spin = hepeup.SPINUP[particle];
105 }
106}
107
108LHEFConverter::~LHEFConverter()
109{
110}
111
112//------------------------------------------------------------------------------
113LHEFConverter::LHEFConverter(const string& inputFileList, const string& outputFileName, const PdgTable& pdg, const int& Nevents) :
114 DataConverter(pdg,Nevents) {
115
116 ExRootTreeWriter *treeWriter = new ExRootTreeWriter(outputFileName, "GEN");
117
118 // generated event from LHEF
119 ExRootTreeBranch *branchEvent = treeWriter->NewBranch("Event", TRootLHEFEvent::Class());
120
121 // generated partons from LHEF
122 ExRootTreeBranch *branchParticle = treeWriter->NewBranch("Particle", TRootLHEFParticle::Class());
123
124 // Open a stream connected to an event file:
125 ifstream infile(inputFileList.c_str());
126 string filename;
127 if(!infile.is_open()) {
128 cerr << left << setw(30) <<"** ERROR: Can't open "<<""
129 << left << setw(20) << inputFileList <<""
130 << right << setw(19) <<"for input **"<<""<<endl;
131 exit(1);
132 }
133
134 while(1) { // parsing the list of files
135
136 infile >> filename;
137 if(!infile.good()) break;
138 ifstream checking_the_file(filename.c_str());
139 if(!checking_the_file.good())
140 {
141 cerr << left << setw(30) <<"** ERROR: Can't find file "<<""
142 << left << setw(20) << filename <<""
143 << right << setw(19) <<"for input **"<<""<<endl;
144 continue;
145 }
146 else checking_the_file.close();
147
148 // Create the Reader object:
149 LHEF::Reader *inputReader = new LHEF::Reader(filename);
150
151 Long64_t allEntries = inputReader->getNumberOfEvents();
152 allEntries = (Nevt<0)?allEntries:min((int)allEntries,Nevt); // do not miss the "+1"
153
154 if(allEntries > 0) {
155 // Loop over all events
156 Long64_t entry = 0;
157 while(inputReader->readEvent()) {
158 treeWriter->Clear();
159 AnalyseEvent(inputReader, branchEvent, entry + 1);
160 AnalyseParticles(inputReader, branchParticle);
161 treeWriter->Fill();
162 ++entry;
163 if(allEntries<entry+1) break;
164 }
165 }
166 }
167 treeWriter->Write();
168
169
170 delete treeWriter;
171 //delete inputReader;
172}
173
Note: See TracBrowser for help on using the repository browser.