Fork me on GitHub

source: svn/trunk/src/LHCOConverter.cc@ 313

Last change on this file since 313 was 313, checked in by severine ovyn, 16 years ago

OK etmis and plus joli

File size: 10.0 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// local includes
33#include "LHCOConverter.h"
34#include "ExRootTreeReader.h"
35#include "BlockClasses.h"
36
37// root includes
38#include "TFile.h"
39#include "TChain.h"
40#include "TClonesArray.h"
41
42// c++ includes
43#include <iostream>
44#include <string>
45#include <fstream>
46#include <sstream>
47#include <iomanip>
48using namespace std;
49
50
51//------------------------------------------------------------------------------
52extern const float UNDEFINED;
53
54LHCOConverter::LHCOConverter(const string& inputroot, const string& inputlog):
55 inputfilename_(inputroot), inputlogname_(inputlog), outputfilename_(inputroot), valid_(true) {
56 if (inputfilename_ == "") {valid_ = false;}
57 else {
58 TFile ftest(inputfilename_.c_str(),"READ");
59 if (!ftest.IsOpen()) { cout << "ERROR: file " << inputfilename_ << " could not be opened\n"; valid_=false;
60 }
61 else {
62 ftest.Close();
63 const unsigned int r_find = outputfilename_.rfind(".root");
64 if(r_find > outputfilename_.size()) {
65 cerr << "ERROR: the extension of the input file is not '.root'. Exiting...\n"; valid_=false;
66 }
67 else {
68 outputfilename_.replace(r_find,5,"_events.lhco");
69 ofstream outfile( outputfilename_.c_str());
70 outfile.close();
71 cout << left << setw(20) <<"** INFO: "<<""
72 << left << setw(29) << outputfilename_ <<""
73 << right << setw(20) <<"created **"<<endl;
74
75 }
76 }
77 }
78}
79
80
81
82void LHCOConverter::CopyRunLogFile() {
83 if(!valid_) {/*cout << "Nothing done\n";*/ return; }
84 if(inputlogname_ == "") {
85 // re-creates the logrun file name from the rootfile name
86 inputlogname_ = inputfilename_;
87 const unsigned int r_find = inputlogname_.rfind(".root");
88 inputlogname_.replace(r_find,5,"_run.log");
89 }
90 ifstream infile( inputlogname_.c_str() );
91 if (!infile.good()) { cout << "Warning: no input logfile found\n"; return; }
92
93 ofstream outfile( outputfilename_.c_str(),ios_base::app);
94
95 // else, if find is found
96 string linereader;
97 while ( getline(infile,linereader) ) {
98 outfile << " # " << linereader << endl;
99 }
100 outfile << " #" << endl;
101 infile.close();
102 outfile.close();
103}
104
105
106
107void LHCOConverter::ConvertExRootAnalysisToLHCO() {
108 if(!valid_) {/*cout << "Nothing done\n";*/ return; }
109
110 //output file
111 ofstream outfile( outputfilename_.c_str(),ios_base::app);
112 //# typ eta phi pt jmas ntrk btag had/em dum1 dum2
113 outfile << " ## More info on LHCO files: http://cp3wks05.fynu.ucl.ac.be/Manual/lhco.html\n #" << endl;
114 outfile << " # typ eta phi pt jmas ntrk btag had/em dum1 dum2" << endl;
115
116 // input files
117 TChain analysisChain("Analysis");
118 TChain triggerChain("Trigger");
119 analysisChain.Add(inputfilename_.c_str());
120 triggerChain.Add(inputfilename_.c_str());
121 ExRootTreeReader *analysisTree = new ExRootTreeReader(&analysisChain);
122 ExRootTreeReader *triggerTree = new ExRootTreeReader(&triggerChain);
123
124 TClonesArray *branchPhoton = analysisTree->UseBranch("Photon");
125 TClonesArray *branchElectron = analysisTree->UseBranch("Electron");
126 TClonesArray *branchMuon = analysisTree->UseBranch("Muon");
127 TClonesArray *branchTauJet = analysisTree->UseBranch("TauJet");
128 TClonesArray *branchJet = analysisTree->UseBranch("Jet");
129 TClonesArray *branchETmis = analysisTree->UseBranch("ETmis");
130
131 Long64_t Nevents = analysisTree->GetEntries();
132 Nevents = min(Nevents,triggerTree->GetEntries());
133 if (Nevents != analysisTree->GetEntries())
134 cout << "** WARNING: not the same number of entries in 'Analysis' and **\n** 'Trigger' trees\n";
135 for(Long64_t event = 0; event < Nevents; ++event) {
136 analysisTree->ReadEntry(event);
137 triggerTree->ReadEntry(event);
138 unsigned int triginfo = 0;
139 unsigned int line =0;
140
141 outfile << setw(3) << line++ << setw(4) << " " << setw(7) << event << setw(7) << triginfo << endl;
142
143 // 0 photon data
144 BranchReader(branchPhoton,line,lhcoPhotonID);
145
146 // 1 electron/positron data
147 BranchReader(branchElectron,line,lhcoElectronID);
148
149 // 2 muon data
150 BranchReader(branchMuon,line,lhcoMuonID);
151
152 // 3 tau-jets
153 BranchReader(branchTauJet,line,lhcoTauJetID);
154
155 // 4 jets
156 BranchReader(branchJet,line,lhcoJetID);
157
158 // 6 MET
159 BranchReaderETmis(branchETmis,line);
160
161 } // event loop
162
163 outfile.close();
164 delete triggerTree;
165 delete analysisTree;
166}
167
168
169//ostringstream LHCOConverter::BranchReader(const TClonesArray* branch, unsigned int& line, unsigned short int lhcoID) const {
170void LHCOConverter::BranchReader(const TClonesArray* branch, unsigned int& line, unsigned short int lhcoID) const {
171 //ostringstream outfile;
172 ofstream outfile( outputfilename_.c_str(),ios_base::app);
173 unsigned int branch_size = branch->GetEntries();
174 TRootParticle * particle = 0;
175 for (unsigned int i=0; i< branch_size; i++) {
176 particle = (TRootParticle*) branch->At(i);
177 TLorentzVector pmu;
178 pmu.SetPtEtaPhiE(particle->PT, particle->Eta, particle->Phi, particle->E);
179 float jmass = pmu.M();
180 float ntrk = 0.0;
181 float btag =0;
182 double ratioE = 0;
183
184 if(lhcoID == lhcoPhotonID){TRootPhoton gam(*((TRootPhoton*) branch->At(i))); jmass=0; ratioE = gam.EHoverEE;}
185 if(lhcoID == lhcoElectronID) { TRootElectron elec(*((TRootElectron*) branch->At(i))); ntrk = elec.Charge; ratioE = elec.EHoverEE;}
186 else if (lhcoID == lhcoMuonID) { TRootMuon muon(*((TRootMuon*) branch->At(i))); ntrk = muon.Charge; ratioE = muon.EHoverEE;}
187 else if (lhcoID == lhcoTauJetID) { TRootTauJet taujet(*((TRootTauJet*) branch->At(i))); ntrk = taujet.Charge; }
188 else if (lhcoID == lhcoJetID) { TRootJet jet(*((TRootJet*) branch->At(i))); ntrk = jet.NTracks; btag = jet.Btag; ratioE = jet.EHoverEE;}
189 if(lhcoID != lhcoETmisID) {
190 outfile << fixed << setprecision(3)
191 << setw(3) << line++ // line counter
192 << setw(4) << lhcoID << " " // object ID in lhco format
193 << setw(7) << particle->Eta <<" "
194 << setw(7) << particle->Phi <<" "
195 << setw(7) << particle->PT <<" "
196 << setw(7) << jmass <<" " // invariant mass
197 << setw(7) << ntrk <<" " // number of tracks
198 << setw(7) << btag <<" "
199 << setw(7) << ratioE <<" " // E_had/E_em
200 << setw(7) << 0.0 << setw(7) << 0.0 // dummy/dummy
201 << endl;
202 } else {
203 outfile << fixed << setprecision(3)
204 << setw(3) << line++ // line counter
205 << setw(4) << lhcoID << " " // object ID in lhco format
206 << setw(7) << 0.000 <<" "
207 << setw(7) << particle->Phi <<" "
208 << setw(7) << particle->PT <<" "
209 << setw(7) << 0. <<" " // invariant mass
210 << setw(7) << 0 <<" " // number of tracks
211 << setw(7) << 0. <<" "
212 << setw(7) << 0. <<" " // E_had/E_em
213 << setw(7) << 0.0 << setw(7) << 0.0 // dummy/dummy
214 << endl;
215 }
216 }
217 outfile.close();
218 // return outfile;
219}
220
221//ostringstream LHCOConverter::BranchReader(const TClonesArray* branch, unsigned int& line, unsigned short int lhcoID) const {
222void LHCOConverter::BranchReaderETmis(const TClonesArray* branch, unsigned int& line) const {
223 //ostringstream outfile;
224 ofstream outfile( outputfilename_.c_str(),ios_base::app);
225 TRootETmis * particle = 0;
226 particle = (TRootETmis*) branch->At(0);
227 outfile << fixed << setprecision(3)
228 << setw(3) << line++ // line counter
229 << setw(4) << 6 << " " // object ID in lhco format
230 << setw(7) << 0.000 <<" "
231 << setw(7) << particle->Phi <<" "
232 << setw(7) << particle->ET <<" "
233 << setw(7) << 0. <<" " // invariant mass
234 << setw(7) << 0 <<" " // number of tracks
235 << setw(7) << 0. <<" "
236 << setw(7) << 0. <<" " // E_had/E_em
237 << setw(7) << 0.0 << setw(7) << 0.0 // dummy/dummy
238 << endl;
239
240 outfile.close();
241 // return outfile;
242}
243
244
Note: See TracBrowser for help on using the repository browser.