Fork me on GitHub

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

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

previously called JetUtils

File size: 4.7 KB
Line 
1/*
2 * ---- Delphes ----
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// \brief Trigger class, and some generic definitions
12
13#include "JetsUtil.h"
14
15using namespace std;
16
17JetsUtil::JetsUtil() : plugins(NULL) {
18 DET = new RESOLution();
19 init();
20}
21
22
23JetsUtil::JetsUtil(const string & DetDatacard) : plugins(NULL) {
24 DET = new RESOLution();
25 DET->ReadDataCard(DetDatacard);
26 init();
27}
28
29
30JetsUtil::JetsUtil(const RESOLution * DetDatacard) : plugins(NULL) {
31 DET = new RESOLution(*DetDatacard);
32 init();
33}
34
35
36JetsUtil::JetsUtil(const JetsUtil& jet){
37 DET = new RESOLution(*(jet.DET));
38}
39
40JetsUtil& JetsUtil::operator=(const JetsUtil& jet) {
41 if(this==&jet) return *this;
42 DET = new RESOLution(*(jet.DET));
43 return *this;
44}
45
46void JetsUtil::init() {
47 if(plugins) delete plugins;
48
49 switch(DET->JET_jetalgo) {
50 default:
51 case 1: {
52 // set up a CDF midpoint jet definition
53#ifdef ENABLE_PLUGIN_CDFCONES
54 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);
55 jet_def = fastjet::JetDefinition(plugins);
56#else
57 plugins = NULL;
58#endif
59 }
60 break;
61
62 case 2: {
63 // set up a CDF midpoint jet definition
64#ifdef ENABLE_PLUGIN_CDFCONES
65 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);
66 jet_def = fastjet::JetDefinition(plugins);
67#else
68 plugins = NULL;
69#endif
70 }
71 break;
72
73 case 3: {
74 // set up a siscone jet definition
75#ifdef ENABLE_PLUGIN_SISCONE
76 plugins = new fastjet::SISConePlugin (DET->JET_coneradius,DET->JET_overlap,DET->JET_S_npass, DET->JET_S_protojet_ptmin);
77 jet_def = fastjet::JetDefinition(plugins);
78#else
79 plugins = NULL;
80#endif
81 }
82 break;
83
84 case 4: {
85
86 jet_def = fastjet::JetDefinition(fastjet::kt_algorithm, DET->JET_coneradius);
87 }
88 break;
89
90 case 5: {
91 jet_def = fastjet::JetDefinition(fastjet::cambridge_algorithm,DET->JET_coneradius);
92 }
93 break;
94
95 case 6: {
96 jet_def = fastjet::JetDefinition(fastjet::antikt_algorithm,DET->JET_coneradius);
97 }
98 break;
99 }
100
101}
102
103vector<fastjet::PseudoJet> JetsUtil::RunJets(const vector<fastjet::PseudoJet>& input_particles)
104{
105 inclusive_jets.clear();
106 sorted_jets.clear();
107 // run the jet clustering with the above jet definition
108 if(input_particles.size()!=0)
109 {
110 fastjet::ClusterSequence clust_seq(input_particles, jet_def);
111 // extract the inclusive jets with pt > 5 GeV
112 double ptmin = 5.0;
113 inclusive_jets = clust_seq.inclusive_jets(ptmin);
114 // sort jets into increasing pt
115 sorted_jets = sorted_by_pt(inclusive_jets);
116 }
117
118 return sorted_jets;
119}
120
121void JetsUtil::RunJetBtagging(ExRootTreeWriter *treeWriter, ExRootTreeBranch *branchJet,const vector<fastjet::PseudoJet> & sorted_jets,const TSimpleArray<TRootGenParticle>& NFCentralQ)
122{
123 TRootJet *elementJet;
124 TLorentzVector JET;
125 for (unsigned int i = 0; i < sorted_jets.size(); i++) {
126 JET.SetPxPyPzE(sorted_jets[i].px(),sorted_jets[i].py(),sorted_jets[i].pz(),sorted_jets[i].E());
127 if(JET.Pt() > DET->PTCUT_jet)
128 {
129 elementJet = (TRootJet*) branchJet->NewEntry();
130 elementJet->Set(JET);
131 // b-jets
132 bool btag=false;
133 if((fabs(JET.Eta()) < DET->CEN_max_tracker && DET->Btaggedjet(JET, NFCentralQ)))btag=true;
134 elementJet->Btag = btag;
135 }
136 } // for itJet : loop on all jets
137
138}
139
140void JetsUtil::RunTauJets(ExRootTreeWriter *treeWriter, ExRootTreeBranch *branchTauJet,const vector<fastjet::PseudoJet> & sorted_jets,const vector<PhysicsTower> & towers, const vector<TLorentzVector> & TrackCentral)
141{
142 TRootTauJet *elementTauJet;
143 TLorentzVector JET;
144 for (unsigned int i = 0; i < sorted_jets.size(); i++) {
145 JET.SetPxPyPzE(sorted_jets[i].px(),sorted_jets[i].py(),sorted_jets[i].pz(),sorted_jets[i].E());
146 // Tau jet identification : 1! track and electromagnetic collimation
147 if(fabs(JET.Eta()) < (DET->CEN_max_tracker - DET->TAU_track_scone)) {
148 double Energie_tau_central = DET->EnergySmallCone(towers,JET.Eta(),JET.Phi());
149 if(
150 ( Energie_tau_central/JET.E() > DET->TAU_energy_frac ) &&
151 ( DET->NumTracks(TrackCentral,DET->TAU_track_pt,JET.Eta(),JET.Phi()) == 1 ) &&
152 ( JET.Pt() > DET->PTCUT_taujet)
153 ) {
154 elementTauJet = (TRootTauJet*) branchTauJet->NewEntry();
155 elementTauJet->Set(JET);
156 } // if tau jet
157 } // if JET.eta < tracker - tau_cone : Tau jet identification
158
159 } // for itJet : loop on all jets
160
161
162}
163
Note: See TracBrowser for help on using the repository browser.