Fork me on GitHub

source: svn/trunk/src/STDHEPConverter.cc@ 390

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

new PDG table

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