Fork me on GitHub

source: svn/trunk/src/JetsUtil.cc@ 307

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

LHCO advanced

File size: 7.8 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#include "JetsUtil.h"
33
34using namespace std;
35
36JetsUtil::JetsUtil() : plugins(NULL) {
37 DET = new RESOLution();
38 init();
39}
40
41
42JetsUtil::JetsUtil(const string & DetDatacard) : plugins(NULL) {
43 DET = new RESOLution();
44 DET->ReadDataCard(DetDatacard);
45 init();
46}
47
48
49JetsUtil::JetsUtil(const RESOLution * DetDatacard) : plugins(NULL) {
50 DET = new RESOLution(*DetDatacard);
51 init();
52}
53
54
55JetsUtil::JetsUtil(const JetsUtil& jet){
56 DET = new RESOLution(*(jet.DET));
57}
58
59JetsUtil& JetsUtil::operator=(const JetsUtil& jet) {
60 if(this==&jet) return *this;
61 DET = new RESOLution(*(jet.DET));
62 return *this;
63}
64
65void JetsUtil::init() {
66 if(plugins) delete plugins;
67
68 switch(DET->JET_jetalgo) {
69 default:
70 case 1: {
71 // set up a CDF midpoint jet definition
72#ifdef ENABLE_PLUGIN_CDFCONES
73 plugins = new fastjet::CDFJetCluPlugin(DET->JET_seed,DET->JET_coneradius,DET->JET_C_adjacencycut,DET->JET_C_maxiterations,DET->JET_C_iratch,DET->JET_overlap);
74 jet_def = fastjet::JetDefinition(plugins);
75#else
76 plugins = NULL;
77#endif
78 }
79 break;
80
81 case 2: {
82 // set up a CDF midpoint jet definition
83#ifdef ENABLE_PLUGIN_CDFCONES
84 plugins = new fastjet::CDFMidPointPlugin (DET->JET_seed,DET->JET_coneradius,DET->JET_M_coneareafraction,DET->JET_M_maxpairsize,DET->JET_M_maxiterations,DET->JET_overlap);
85 jet_def = fastjet::JetDefinition(plugins);
86#else
87 plugins = NULL;
88#endif
89 }
90 break;
91
92 case 3: {
93 // set up a siscone jet definition
94#ifdef ENABLE_PLUGIN_SISCONE
95 plugins = new fastjet::SISConePlugin (DET->JET_coneradius,DET->JET_overlap,DET->JET_S_npass, DET->JET_S_protojet_ptmin);
96 jet_def = fastjet::JetDefinition(plugins);
97#else
98 plugins = NULL;
99#endif
100 }
101 break;
102
103 case 4: {
104
105 jet_def = fastjet::JetDefinition(fastjet::kt_algorithm, DET->JET_coneradius);
106 }
107 break;
108
109 case 5: {
110 jet_def = fastjet::JetDefinition(fastjet::cambridge_algorithm,DET->JET_coneradius);
111 }
112 break;
113
114 case 6: {
115 jet_def = fastjet::JetDefinition(fastjet::antikt_algorithm,DET->JET_coneradius);
116 }
117 break;
118 }
119
120}
121
122vector<fastjet::PseudoJet> JetsUtil::RunJets(const vector<fastjet::PseudoJet>& input_particles, const vector<TRootTracks> & TrackCentral, vector<int> &NTrackJet)
123{
124 inclusive_jets.clear();
125 sorted_jets.clear();
126 // run the jet clustering with the above jet definition
127 if(input_particles.size()!=0)
128 {
129 fastjet::ClusterSequence clust_seq(input_particles, jet_def);
130 // extract the inclusive jets with pt > 5 GeV
131 double ptmin = 5.0;
132 inclusive_jets = clust_seq.inclusive_jets(ptmin);
133
134 // sort jets into increasing pt
135 sorted_jets = sorted_by_pt(inclusive_jets);
136
137 //check number of tracks in jets
138 float iEtaTrack[TrackCentral.size()];
139 float iPhiTrack[TrackCentral.size()];
140 for(unsigned int t = 0; t < TrackCentral.size(); t++)
141 {
142 DET->BinEtaPhi(TrackCentral[t].PhiOuter,TrackCentral[t].EtaOuter,iPhiTrack[t],iEtaTrack[t]);
143 }
144 int numTrackJet;
145 for (unsigned int i = 0; i < sorted_jets.size(); i++)
146 {
147 numTrackJet=0;
148 vector<fastjet::PseudoJet> constituents = clust_seq.constituents(sorted_jets[i]);
149 for(unsigned int it = 0; it < TrackCentral.size(); it++)
150 {
151 for (unsigned int j = 0; j < constituents.size(); j++)
152 {
153 if(DeltaR(iPhiTrack[it], iEtaTrack[it], constituents[j].phi(), constituents[j].eta())<0.001)numTrackJet++;
154 }
155 }
156 NTrackJet.push_back(numTrackJet);
157 }
158 }
159
160 return sorted_jets;
161}
162
163
164
165void JetsUtil::RunJetBtagging(ExRootTreeWriter *treeWriter, ExRootTreeBranch *branchJet,const vector<fastjet::PseudoJet> & sorted_jets,const TSimpleArray<TRootGenParticle>& NFCentralQ, const vector<int> &NTrackJet)
166{
167 TRootJet *elementJet;
168 TLorentzVector JET;
169 for (unsigned int i = 0; i < sorted_jets.size(); i++) {
170
171 JET.SetPxPyPzE(sorted_jets[i].px(),sorted_jets[i].py(),sorted_jets[i].pz(),sorted_jets[i].E());
172 if(JET.Pt() > DET->PTCUT_jet)
173 {
174 elementJet = (TRootJet*) branchJet->NewEntry();
175 elementJet->Set(JET);
176 elementJet->NTracks = NTrackJet[i];
177
178 // b-jets
179 bool btag=false;
180 if((fabs(JET.Eta()) < DET->CEN_max_tracker && DET->Btaggedjet(JET, NFCentralQ)))btag=true;
181 elementJet->Btag = btag;
182 }
183 } // for itJet : loop on all jets
184
185}
186
187void JetsUtil::RunTauJets(ExRootTreeWriter *treeWriter, ExRootTreeBranch *branchTauJet,const vector<fastjet::PseudoJet> & sorted_jets,const vector<PhysicsTower> & towers, const vector<TRootTracks> & TrackCentral, const vector<int> &NTrackJet)
188{
189 TRootTauJet *elementTauJet;
190 TLorentzVector JET;
191 float charge=0;
192 for (unsigned int i = 0; i < sorted_jets.size(); i++) {
193 JET.SetPxPyPzE(sorted_jets[i].px(),sorted_jets[i].py(),sorted_jets[i].pz(),sorted_jets[i].E());
194 // Tau jet identification : 1! track and electromagnetic collimation
195 if(fabs(JET.Eta()) < (DET->CEN_max_tracker - DET->TAU_track_scone)) {
196 double Energie_tau_central = DET->EnergySmallCone(towers,JET.Eta(),JET.Phi());
197 if(
198 ( Energie_tau_central/JET.E() > DET->TAU_energy_frac ) &&
199 //( DET->NumTracks(charge,TrackCentral,DET->TAU_track_pt,JET.Eta(),JET.Phi(),DET->TAU_track_scone) == 1 ) &&
200 ( DET->NumTracks(charge,TrackCentral,DET->TAU_track_pt,JET.Eta(),JET.Phi()) == 1 ) &&
201 ( JET.Pt() > DET->PTCUT_taujet)
202 ) {
203 elementTauJet = (TRootTauJet*) branchTauJet->NewEntry();
204 elementTauJet->Set(JET);
205 elementTauJet->NTracks=NTrackJet[i];
206 elementTauJet->Charge = charge;
207 } // if tau jet
208 } // if JET.eta < tracker - tau_cone : Tau jet identification
209
210 } // for itJet : loop on all jets
211
212
213}
214
Note: See TracBrowser for help on using the repository browser.