Fork me on GitHub

source: svn/trunk/src/STDHEPConverter.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: 7.4 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#include "SmearUtil.h"
38#include "ExRootTreeWriter.h"
39#include "ExRootTreeBranch.h"
40#include "stdhep_mcfio.h"
41#include "stdhep_declarations.h"
42#include "STDHEPConverter.h"
43
44using namespace std;
45
46//---------------------------------------------------------------------------
47
48void STDHEPConverter::AnalyseEvent(ExRootTreeBranch *branch, const Long64_t eventNumber)
49{
50 TRootGenEvent *element;
51
52 element = static_cast<TRootGenEvent*>(branch->NewEntry());
53
54 element->Number = eventNumber;
55}
56
57//---------------------------------------------------------------------------
58
59void STDHEPConverter::AnalyseParticles(ExRootTreeBranch *branch)
60{
61 int Nmax = -1; // useless for the moment -- should be changed in order
62 if(Nmax>0) Nmax = min(Nmax,myhepevt.nhep); // to use Nevt in DataConverter
63 else Nmax = myhepevt.nhep;
64
65 TRootC::GenParticle *element;
66
67 Double_t signPz;
68 TLorentzVector momentum;
69 Int_t number;
70
71 for(number = 0; number < Nmax; ++number)
72 {
73 element = static_cast<TRootC::GenParticle*>(branch->NewEntry());
74
75 element->PID = myhepevt.idhep[number];
76 element->Status = myhepevt.isthep[number];
77 element->M1 = myhepevt.jmohep[number][0] - 1;
78 element->M2 = myhepevt.jmohep[number][1] - 1;
79 element->D1 = myhepevt.jdahep[number][0] - 1;
80 element->D2 = myhepevt.jdahep[number][1] - 1;
81// element->Charge = myhepevt.hepchg(element->PID);
82 element->Charge = ChargeVal(element->PID);
83
84 element->E = myhepevt.phep[number][3];
85 element->Px = myhepevt.phep[number][0];
86 element->Py = myhepevt.phep[number][1];
87 element->Pz = myhepevt.phep[number][2];
88 element->M = myhepevt.phep[number][4];
89
90 momentum.SetPxPyPzE(element->Px, element->Py, element->Pz, element->E);
91 element->PT = momentum.Perp();
92 signPz = (element->Pz >= 0.0) ? 1.0 : -1.0;
93 //element->Eta = element->PT == 0.0 ? signPz*999.9 : momentum.Eta(); to avoid a warning from ROOT, replace the "==0" by "< 1e-6"
94 element->Eta = element->PT <1e-6 ? signPz*999.9 : momentum.Eta();
95 element->Phi = momentum.Phi();
96
97 element->T = myhepevt.vhep[number][3];
98 element->X = myhepevt.vhep[number][0];
99 element->Y = myhepevt.vhep[number][1];
100 element->Z = myhepevt.vhep[number][2];
101 }
102}
103
104
105//------------------------------------------------------------------------------
106
107STDHEPConverter::~STDHEPConverter()
108{
109}
110
111//------------------------------------------------------------------------------
112
113STDHEPConverter::STDHEPConverter(const string& inputFileList, const string& outputFileName, const int& Nevents) : DataConverter(Nevents)
114{
115 int ierr, entryType;
116 int istr = 0;
117 int nevt = 0;
118
119 ExRootTreeWriter *treeWriter = new ExRootTreeWriter(outputFileName, "GEN");
120
121 // information about generated event
122 ExRootTreeBranch *branchGenEvent = treeWriter->NewBranch("Event", TRootGenEvent::Class());
123 // generated particles from HEPEVT
124 ExRootTreeBranch *branchGenParticle = treeWriter->NewBranch("Particle", TRootC::GenParticle::Class());
125
126 // Open a stream connected to an event file:
127 ifstream infile(inputFileList.c_str());
128 string filename;
129 if(!infile.is_open())
130 {
131 cerr << left << setw(30) <<"** ERROR: Can't open "<<""
132 << left << setw(20) << inputFileList <<""
133 << right << setw(19) <<"for input **"<<""<<endl;
134
135 exit(1);
136 }
137 while(1)
138 {
139 infile >> filename;
140 if(!infile.good()) break;
141 ifstream checking_the_file(filename.c_str());
142 if(!checking_the_file.good())
143 {
144 cerr << left << setw(30) <<"** ERROR: Can't find file "<<""
145 << left << setw(20) << filename <<""
146 << right << setw(19) <<"for input **"<<""<<endl;
147 continue;
148 }
149 else checking_the_file.close();
150
151 ierr = StdHepXdrReadInit(const_cast<char*>(filename.c_str()), &nevt, istr);
152
153 if(ierr != 0)
154 {
155 cerr <<"** ERROR: Can't open file for input **"<<endl;
156 break;
157 }
158
159 Long64_t allEntries = nevt;
160
161 if(allEntries > 0)
162 {
163 // Loop over all objects
164 Long64_t entry = 0;
165 Long64_t recordNumber = 1;
166 for(entry = 0; entry < allEntries; ++entry)
167 {
168 ierr = StdHepXdrRead(&entryType, istr);
169
170 if(ierr != 0)
171 {
172 cerr << left << setw(49) <<"** ERROR: Unexpected end of file after"<<""
173 << left << setw(10) << entry <<""
174 << right << setw(10) <<"entries **"<<""<<endl;
175 break;
176 }
177
178 // analyse only entries with standard HEPEVT common block
179 if(entryType == 1)
180 {
181 // add empty events for missing event numbers
182 while(recordNumber < myhepevt.nevhep)
183 {
184 treeWriter->Clear();
185 AnalyseEvent(branchGenEvent, recordNumber);
186 treeWriter->Fill();
187 ++recordNumber;
188 }
189
190 treeWriter->Clear();
191
192 AnalyseEvent(branchGenEvent, myhepevt.nevhep);
193 AnalyseParticles(branchGenParticle);
194
195 treeWriter->Fill();
196
197 ++recordNumber;
198
199 }
200 }
201 }
202 StdHepXdrEnd(istr);
203 }
204 treeWriter->Write();
205
206 delete treeWriter;
207
208}
209
210
211
212
Note: See TracBrowser for help on using the repository browser.