Fork me on GitHub

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

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

checks that the input file is present

File size: 4.8 KB
RevLine 
[2]1#include <iostream>
2#include <utility>
3#include <deque>
4#include <sstream>
5#include <fstream>
6
7#include "TROOT.h"
8#include "TApplication.h"
9#include "TFile.h"
10#include "TTree.h"
11#include "TBranch.h"
12#include "TLeaf.h"
13#include "TLorentzVector.h"
14#include "TChain.h"
15
16#include "Utilities/ExRootAnalysis/interface/BlockClasses.h"
17#include "Utilities/ExRootAnalysis/interface/ExRootTreeReader.h"
18#include "Utilities/ExRootAnalysis/interface/ExRootTreeWriter.h"
19#include "Utilities/ExRootAnalysis/interface/ExRootTreeBranch.h"
20#include "Utilities/ExRootAnalysis/interface/LHEF.h"
21#include "Utilities/stdhep/stdhep_mcfio.h"
22#include "Utilities/stdhep/stdhep_declarations.h"
23
24#include "interface/STDHEPConverter.h"
25
26
27
28using namespace std;
29
30//---------------------------------------------------------------------------
31
32void STDHEPConverter::AnalyseEvent(ExRootTreeBranch *branch, const Long64_t eventNumber)
33{
34 TRootGenEvent *element;
35
36 element = static_cast<TRootGenEvent*>(branch->NewEntry());
37
38 element->Number = eventNumber;
39}
40
41//---------------------------------------------------------------------------
42
43void STDHEPConverter::AnalyseParticles(ExRootTreeBranch *branch)
44{
45 TRootGenParticle *element;
46
47 Double_t signPz;
48 TLorentzVector momentum;
49 Int_t number;
50
51 for(number = 0; number < myhepevt.nhep; ++number)
52 {
53
54 element = static_cast<TRootGenParticle*>(branch->NewEntry());
55
56 element->PID = myhepevt.idhep[number];
57 element->Status = myhepevt.isthep[number];
58 element->M1 = myhepevt.jmohep[number][0] - 1;
59 element->M2 = myhepevt.jmohep[number][1] - 1;
60 element->D1 = myhepevt.jdahep[number][0] - 1;
61 element->D2 = myhepevt.jdahep[number][1] - 1;
62// element->Charge = myhepevt.hepchg(element->PID);
63
64 element->E = myhepevt.phep[number][3];
65 element->Px = myhepevt.phep[number][0];
66 element->Py = myhepevt.phep[number][1];
67 element->Pz = myhepevt.phep[number][2];
68
69 momentum.SetPxPyPzE(element->Px, element->Py, element->Pz, element->E);
70 element->PT = momentum.Perp();
71 signPz = (element->Pz >= 0.0) ? 1.0 : -1.0;
72 element->Eta = element->PT == 0.0 ? signPz*999.9 : momentum.Eta();
73 element->Phi = momentum.Phi();
74
75 element->T = myhepevt.vhep[number][3];
76 element->X = myhepevt.vhep[number][0];
77 element->Y = myhepevt.vhep[number][1];
78 element->Z = myhepevt.vhep[number][2];
79 }
80}
81
82
83//------------------------------------------------------------------------------
84
85STDHEPConverter::~STDHEPConverter()
86{
87}
88
89//------------------------------------------------------------------------------
90
91STDHEPConverter::STDHEPConverter(const string& inputFileList, const string& outputFileName)
92{
93
94 int ierr, entryType;
95 int istr = 0;
96 int nevt = 0;
97
98 ExRootTreeWriter *treeWriter = new ExRootTreeWriter(outputFileName, "GEN");
99
100 // information about generated event
101 ExRootTreeBranch *branchGenEvent = treeWriter->NewBranch("Event", TRootGenEvent::Class());
102 // generated particles from HEPEVT
103 ExRootTreeBranch *branchGenParticle = treeWriter->NewBranch("Particle", TRootGenParticle::Class());
104
105 // Open a stream connected to an event file:
106 ifstream infile(inputFileList.c_str());
107 string filename;
108 if(!infile.is_open())
109 {
110 cerr << "** ERROR: Can't open '" << inputFileList << "' for input" << endl;
111 exit(1);
112 }
113 while(1)
114 {
115 infile >> filename;
116 if(!infile.good()) break;
[18]117 ifstream checking_the_file(filename.c_str());
118 if(!checking_the_file.good()) { cout << filename << ": file not found\n"; continue;}
119 else checking_the_file.close();
[2]120
121 ierr = StdHepXdrReadInit(const_cast<char*>(filename.c_str()), &nevt, istr);
122
123 if(ierr != 0)
124 {
125 cerr << "** ERROR: Can't open file for input" << endl;
126 break;
127 }
128
129 Long64_t allEntries = nevt;
130 cout << "** Input file contains " << allEntries << " entries" << endl;
131
132 if(allEntries > 0)
133 {
134 // Loop over all objects
135 Long64_t entry = 0;
136 Long64_t recordNumber = 1;
137 for(entry = 0; entry < allEntries; ++entry)
138 {
139 ierr = StdHepXdrRead(&entryType, istr);
140
141 if(ierr != 0)
142 {
143 cerr << "** ERROR: Unexpected end of file after " << entry << " entries" << endl;
144 break;
145 }
146
147 // analyse only entries with standard HEPEVT common block
148 if(entryType == 1)
149 {
150 // add empty events for missing event numbers
151 while(recordNumber < myhepevt.nevhep)
152 {
153 treeWriter->Clear();
154 AnalyseEvent(branchGenEvent, recordNumber);
155 treeWriter->Fill();
156 ++recordNumber;
157 }
158
159 treeWriter->Clear();
160
161 AnalyseEvent(branchGenEvent, myhepevt.nevhep);
162 AnalyseParticles(branchGenParticle);
163
164 treeWriter->Fill();
165
166 ++recordNumber;
167
168 }
169 }
170 }
171 }
172 treeWriter->Write();
173 cout << "** Exiting..." << endl;
174
175 delete treeWriter;
176 StdHepXdrEnd(istr);
177
178}
179
180
181
182
Note: See TracBrowser for help on using the repository browser.