Fork me on GitHub

source: svn/trunk/Resolutions.cpp@ 22

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

plugins at the begining

File size: 12.1 KB
Line 
1/*
2 ---- FastSim ----
3 A Fast Simulator for general purpose LHC detector
4 S. Ovyn ~~~~ severine.ovyn@uclouvain.be
5
6 Center for Particle Physics and Phenomenology (CP3)
7 Universite Catholique de Louvain (UCL)
8 Louvain-la-Neuve, Belgium
9*/
10
11/// \file Smearing.cpp
12/// \brief executable for the FastSim
13
14#include "TChain.h"
15#include "TApplication.h"
16
17#include "Utilities/ExRootAnalysis/interface/ExRootTreeReader.h"
18#include "Utilities/ExRootAnalysis/interface/ExRootTreeWriter.h"
19#include "Utilities/ExRootAnalysis/interface/ExRootTreeBranch.h"
20
21#include "H_BeamParticle.h"
22#include "H_BeamLine.h"
23#include "H_RomanPot.h"
24
25#include "interface/DataConverter.h"
26#include "interface/HEPEVTConverter.h"
27#include "interface/LHEFConverter.h"
28#include "interface/STDHEPConverter.h"
29
30#include "interface/SmearUtil.h"
31#include "Utilities/Fastjet/include/fastjet/PseudoJet.hh"
32#include "Utilities/Fastjet/include/fastjet/ClusterSequence.hh"
33
34// get info on how fastjet was configured
35#include "Utilities/Fastjet/include/fastjet/config.h"
36
37// make sure we have what is needed
38#ifdef ENABLE_PLUGIN_SISCONE
39# include "Utilities/Fastjet/plugins/SISCone/SISConePlugin.hh"
40#endif
41#ifdef ENABLE_PLUGIN_CDFCONES
42# include "Utilities/Fastjet/plugins/CDFCones/CDFMidPointPlugin.hh"
43# include "Utilities/Fastjet/plugins/CDFCones/CDFJetCluPlugin.hh"
44#endif
45
46#include<vector>
47#include<iostream>
48
49#include "interface/TreeClasses.h"
50using namespace std;
51
52//------------------------------------------------------------------------------
53
54
55
56void PairingJet(TLorentzVector &JETSm, TLorentzVector JET, vector<fastjet::PseudoJet> sorted_jetsS)
57{
58 JETSm.SetPxPyPzE(0,0,0,0);
59 float deltaRtest=5000;
60 for (unsigned int i = 0; i < sorted_jetsS.size(); i++) {
61 TLorentzVector Att;
62 Att.SetPxPyPzE(sorted_jetsS[i].px(),sorted_jetsS[i].py(),sorted_jetsS[i].pz(),sorted_jetsS[i].E());
63 if(DeltaR(JET.Phi(),JET.Eta(),Att.Phi(),Att.Eta()) < deltaRtest)
64 {
65 deltaRtest = DeltaR(JET.Phi(),JET.Eta(),Att.Phi(),Att.Eta());
66 if(deltaRtest < 0.25)
67 {
68 JETSm = Att;
69 }
70 }
71 }
72}
73
74
75int main(int argc, char *argv[])
76{
77 int appargc = 2;
78 char *appName = "Smearing";
79 char *appargv[] = {appName, "-b"};
80 TApplication app(appName, &appargc, appargv);
81
82 if(argc != 4 && argc != 3) {
83 cout << " Usage: " << argv[0] << " input_file" << " output_file" << " data_card " << endl;
84 cout << " input_list - list of files in Ntpl, StdHep of LHEF format," << endl;
85 cout << " output_file - output file." << endl;
86 cout << " data_card - Datacard containing resolution variables for the detector simulation (optional) "<<endl;
87 exit(1);
88 }
89
90 srand (time (NULL)); /* Initialisation du générateur */
91
92 //read the input TROOT file
93 string inputFileList(argv[1]), outputfilename(argv[2]);
94 if(outputfilename.find(".root") > outputfilename.length() ) {
95 cout << "output_file should be a .root file!\n";
96 return -1;
97 }
98 TFile *outputFile = TFile::Open(outputfilename.c_str(), "RECREATE"); // Creates the file, but should be closed just after
99 outputFile->Close();
100
101 string line;
102 ifstream infile(inputFileList.c_str());
103 infile >> line; // the first line determines the type of input files
104
105 DataConverter *converter=0;
106
107 if(strstr(line.c_str(),".hep"))
108 {
109 cout<<"*************************************************************************"<<endl;
110 cout<<"************ StdHEP file format detected **************"<<endl;
111 cout<<"************ Starting convertion to TRoot format **************"<<endl;
112 cout<<"*************************************************************************"<<endl;
113 converter = new STDHEPConverter(inputFileList,outputfilename);//case ntpl file in input list
114 }
115 else if(strstr(line.c_str(),".lhe"))
116 {
117 cout<<"*************************************************************************"<<endl;
118 cout<<"************ LHEF file format detected **************"<<endl;
119 cout<<"************ Starting convertion to TRoot format **************"<<endl;
120 cout<<"*************************************************************************"<<endl;
121 converter = new LHEFConverter(inputFileList,outputfilename);//case ntpl file in input list
122 }
123 else if(strstr(line.c_str(),".root"))
124 {
125 cout<<"*************************************************************************"<<endl;
126 cout<<"************ h2root file format detected **************"<<endl;
127 cout<<"************ Starting convertion to TRoot format **************"<<endl;
128 cout<<"*************************************************************************"<<endl;
129 converter = new HEPEVTConverter(inputFileList,outputfilename);//case ntpl file in input list
130 }
131 else { cout << "*** " << line.c_str() << "\n*** file format not identified\n*** Exiting\n"; return -1;};
132
133 TChain chain("GEN");
134 chain.Add(outputfilename.c_str());
135 ExRootTreeReader *treeReader = new ExRootTreeReader(&chain);
136 const TClonesArray *branchGen = treeReader->UseBranch("Particle");
137 TIter itGen((TCollection*)branchGen);
138
139 //write the output root file
140 ExRootTreeWriter *treeWriter = new ExRootTreeWriter(outputfilename, "Analysis");
141 ExRootTreeBranch *branchjet = treeWriter->NewBranch("JetPTResol", RESOLJET::Class());
142 ExRootTreeBranch *branchelec = treeWriter->NewBranch("ElecEResol", RESOLELEC::Class());
143 ExRootTreeBranch *branchmuon = treeWriter->NewBranch("MuonPTResol", RESOLMUON::Class());
144 ExRootTreeBranch *branchtaujet = treeWriter->NewBranch("TauJetPTResol", TAUHAD::Class());
145 ExRootTreeBranch *branchetmis = treeWriter->NewBranch("ETmisResol",ETMIS::Class());
146
147 TRootGenParticle *particle;
148 TRootETmis *etmisc;
149
150 RESOLELEC *elementElec;
151 RESOLMUON *elementMuon;
152 RESOLJET *elementJet;
153 TAUHAD *elementTaujet;
154 ETMIS *elementEtmis;
155
156
157 //read the datacard input file
158 string DetDatacard("");
159 if(argc==4) DetDatacard =argv[3];
160 RESOLution *DET = new RESOLution();
161 DET->ReadDataCard(DetDatacard);
162
163 TLorentzVector genMomentum(0,0,0,0);
164 LorentzVector jetMomentum;
165 vector<TLorentzVector> TrackCentral;
166
167 vector<fastjet::PseudoJet> input_particles;//for FastJet algorithm
168 vector<fastjet::PseudoJet> inclusive_jets;
169 vector<fastjet::PseudoJet> sorted_jets;
170
171 vector<fastjet::PseudoJet> input_particlesS;//for FastJet algorithm
172 vector<fastjet::PseudoJet> inclusive_jetsS;
173 vector<fastjet::PseudoJet> sorted_jetsS;
174
175 fastjet::JetDefinition jet_def;
176 fastjet::JetDefinition jet_defS;
177 fastjet::JetDefinition::Plugin * plugins;
178 fastjet::JetDefinition::Plugin * pluginsS;
179
180 // set up a CDF midpoint jet definition
181 #ifdef ENABLE_PLUGIN_CDFCONES
182 plugins = new fastjet::CDFJetCluPlugin(DET->C_SEEDTHRESHOLD,DET->CONERADIUS,DET->C_ADJACENCYCUT,DET->C_MAXITERATIONS,DET->C_IRATCH,DET->C_OVERLAPTHRESHOLD);
183 jet_def = fastjet::JetDefinition(plugins);
184 #else
185 plugins = NULL;
186 #endif
187
188 // set up a CDF midpoint jet definition
189 #ifdef ENABLE_PLUGIN_CDFCONES
190 pluginsS = new fastjet::CDFJetCluPlugin(2,DET->CONERADIUS,DET->C_ADJACENCYCUT,DET->C_MAXITERATIONS,DET->C_IRATCH,DET->C_OVERLAPTHRESHOLD);
191 jet_defS = fastjet::JetDefinition(pluginsS);
192 #else
193 pluginsS = NULL;
194 #endif
195
196
197 // Loop over all events
198 Long64_t entry, allEntries = treeReader->GetEntries();
199 cout << "** Chain contains " << allEntries << " events" << endl;
200 for(entry = 0; entry < allEntries; ++entry)
201 {
202 TLorentzVector PTmis(0,0,0,0);
203 treeReader->ReadEntry(entry);
204 treeWriter->Clear();
205
206 if((entry % 100) == 0 && entry > 0 ) cout << "** Processing element # " << entry << endl;
207
208 TSimpleArray<TRootGenParticle> bGen;
209 itGen.Reset();
210 TrackCentral.clear();
211 TSimpleArray<TRootGenParticle> NFCentralQ;
212 input_particles.clear();
213 inclusive_jets.clear();
214 sorted_jets.clear();
215 input_particlesS.clear();
216 inclusive_jetsS.clear();
217 sorted_jetsS.clear();
218
219 // Loop over all particles in event
220 while( (particle = (TRootGenParticle*) itGen.Next()) )
221 {
222 genMomentum.SetPxPyPzE(particle->Px, particle->Py, particle->Pz, particle->E);
223
224 int pid = abs(particle->PID);
225 float eta = fabs(particle->Eta);
226
227 if(particle->Status == 1)
228 {
229 input_particles.push_back(fastjet::PseudoJet(genMomentum.Px(),genMomentum.Py(),genMomentum.Pz(), genMomentum.E()));
230 }
231
232 // keeps only final particles, visible by the central detector, including the fiducial volume
233 // the ordering of conditions have been optimised for speed : put first the STATUS condition
234 if( (particle->Status == 1) &&
235 (
236 (pid == pMU && eta < DET->MAX_MU) ||
237 (pid != pMU && (pid != pNU1) && (pid != pNU2) && (pid != pNU3) && eta < DET->MAX_CALO_FWD)
238 )
239 ) {
240 switch(pid) {
241
242 case pE: // all electrons with eta < DET->MAX_CALO_FWD
243 DET->SmearElectron(genMomentum);
244 break; // case pE
245
246 case pGAMMA: // all photons with eta < DET->MAX_CALO_FWD
247 DET->SmearElectron(genMomentum);
248 break; // case pGAMMA
249
250 case pMU: // all muons with eta < DET->MAX_MU
251 DET->SmearMu(genMomentum);
252 break; // case pMU
253
254 case pLAMBDA: // all lambdas with eta < DET->MAX_CALO_FWD
255 case pK0S: // all K0s with eta < DET->MAX_CALO_FWD
256 DET->SmearHadron(genMomentum, 0.7);
257 break; // case hadron
258
259 default: // all other final particles with eta < DET->MAX_CALO_FWD
260 DET->SmearHadron(genMomentum, 1.0);
261 break;
262 } // switch (pid)
263
264 // all final particles but muons and neutrinos
265 // for calorimetric towers and mission PT
266 if(genMomentum.E()!=0) PTmis = PTmis + genMomentum;//ptmis
267
268 if(pid != pMU && genMomentum.Pt() > DET->PT_TRACKS_MIN)
269 {
270 input_particlesS.push_back(fastjet::PseudoJet(genMomentum.Px(),genMomentum.Py(),genMomentum.Pz(), genMomentum.E()));
271 }
272
273 // all final charged particles
274 if(
275 ((rand()%100) < DET->TRACKING_EFF) &&
276 (genMomentum.E()!=0) &&
277 (fabs(particle->Eta) < DET->MAX_TRACKER) &&
278 (genMomentum.Pt() > DET->PT_TRACKS_MIN ) && // pt too small to be taken into account
279 (pid != pGAMMA) &&
280 (pid != pPI0) &&
281 (pid != pK0L) &&
282 (pid != pN) &&
283 (pid != pSIGMA0) &&
284 (pid != pDELTA0) &&
285 (pid != pK0S) // not charged particles : invisible by tracker
286 )
287 TrackCentral.push_back(genMomentum);
288 } // switch
289 } // while
290
291
292 //*****************************
293
294 double ptmin=1;
295 if(input_particles.size()!=0)
296 {
297 fastjet::ClusterSequence clust_seq(input_particles, jet_def);
298 inclusive_jets = clust_seq.inclusive_jets(ptmin);
299 sorted_jets = sorted_by_pt(inclusive_jets);
300 }
301
302 if(input_particlesS.size()!=0)
303 {
304 fastjet::ClusterSequence clust_seqS(input_particlesS, jet_defS);
305 inclusive_jetsS = clust_seqS.inclusive_jets(ptmin);
306 sorted_jetsS = sorted_by_pt(inclusive_jetsS);
307 }
308
309 TLorentzVector JETSm(0,0,0,0);
310 for (unsigned int i = 0; i < sorted_jets.size(); i++) {
311 TLorentzVector JET(0,0,0,0);
312 JET.SetPxPyPzE(sorted_jets[i].px(),sorted_jets[i].py(),sorted_jets[i].pz(),sorted_jets[i].E());
313 PairingJet(JETSm,JET,sorted_jetsS);
314 if(JETSm.Pt()>1)
315 {
316 elementJet= (RESOLJET*) branchjet->NewEntry();
317 elementJet->NonSmearePT = JET.Pt();
318 elementJet->SmearePT = (JETSm.Pt()/JET.Pt());
319 /*cout<<"valeur obtenue "<<JETSm.Pt()/JET.Pt()<<endl;
320 cout<<" pt non smeare "<<JET.Pt()<<" phi "<<JET.Phi()<<" eta "<<JET.Eta()<<endl;
321 cout<<"pt smeare "<<JETSm.Pt()<<" phi "<<JETSm.Phi()<<" eta "<<JETSm.Eta()<<endl;
322 cout<<"************"<<endl;
323*/
324 }
325
326 } // for itJet : loop on all jets
327
328 treeWriter->Fill();
329 } // Loop over all events
330 treeWriter->Write();
331
332 cout << "** Exiting..." << endl;
333
334 delete treeWriter;
335 delete treeReader;
336 delete DET;
337 if(converter) delete converter;
338}
339
Note: See TracBrowser for help on using the repository browser.