Fork me on GitHub

source: svn/trunk/src/JetUtils.cc@ 73

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

remove fill bug

File size: 4.0 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 "interface/JetUtils.h"
14
15using namespace std;
16
17JetsUtil::JetsUtil() {
18
19 switch(JETALGO) {
20 default:
21 case 1: {
22
23 // set up a CDF midpoint jet definition
24#ifdef ENABLE_PLUGIN_CDFCONES
25 plugins = new fastjet::CDFJetCluPlugin(SEEDTHRESHOLD,CONERADIUS,C_ADJACENCYCUT,C_MAXITERATIONS,C_IRATCH,OVERLAPTHRESHOLD);
26 jet_def = fastjet::JetDefinition(plugins);
27#else
28 plugins = NULL;
29#endif
30 }
31 break;
32
33 case 2: {
34 // set up a CDF midpoint jet definition
35#ifdef ENABLE_PLUGIN_CDFCONES
36 plugins = new fastjet::CDFMidPointPlugin (SEEDTHRESHOLD,CONERADIUS,M_CONEAREAFRACTION,M_MAXPAIRSIZE,M_MAXITERATIONS,OVERLAPTHRESHOLD);
37 jet_def = fastjet::JetDefinition(plugins);
38#else
39 plugins = NULL;
40#endif
41 }
42 break;
43
44 case 3: {
45 // set up a siscone jet definition
46#ifdef ENABLE_PLUGIN_SISCONE
47 plugins = new fastjet::SISConePlugin (CONERADIUS,OVERLAPTHRESHOLD,NPASS, PROTOJET_PTMIN);
48 jet_def = fastjet::JetDefinition(plugins);
49#else
50 plugins = NULL;
51#endif
52 }
53 break;
54
55 case 4: {
56 jet_def = fastjet::JetDefinition(fastjet::kt_algorithm, CONERADIUS);
57 }
58 break;
59
60 case 5: {
61 jet_def = fastjet::JetDefinition(fastjet::cambridge_algorithm,CONERADIUS);
62 }
63 break;
64
65 case 6: {
66 jet_def = fastjet::JetDefinition(fastjet::antikt_algorithm,CONERADIUS);
67 }
68 break;
69 }
70
71}
72
73vector<fastjet::PseudoJet> JetsUtil::RunJets(const vector<fastjet::PseudoJet> & input_particles)
74{
75 inclusive_jets.clear();
76 sorted_jets.clear();
77
78 // run the jet clustering with the above jet definition
79 if(input_particles.size()!=0)
80 {
81 fastjet::ClusterSequence clust_seq(input_particles, jet_def);
82 // extract the inclusive jets with pt > 5 GeV
83 double ptmin = 5.0;
84 inclusive_jets = clust_seq.inclusive_jets(ptmin);
85 // sort jets into increasing pt
86 sorted_jets = sorted_by_pt(inclusive_jets);
87 }
88
89 return sorted_jets;
90}
91
92void JetsUtil::RunJetBtagging(ExRootTreeWriter *treeWriter, ExRootTreeBranch *branchJet,const vector<fastjet::PseudoJet> & sorted_jets,const TSimpleArray<TRootGenParticle>& NFCentralQ)
93{
94 TRootJet *elementJet;
95 TLorentzVector JET;
96 for (unsigned int i = 0; i < sorted_jets.size(); i++) {
97 JET.SetPxPyPzE(sorted_jets[i].px(),sorted_jets[i].py(),sorted_jets[i].pz(),sorted_jets[i].E());
98 if(JET.Pt() > JET_pt)
99 {
100 elementJet = (TRootJet*) branchJet->NewEntry();
101 elementJet->Set(JET);
102 // b-jets
103 bool btag=false;
104 if((fabs(JET.Eta()) < MAX_TRACKER && Btaggedjet(JET, NFCentralQ)))btag=true;
105 elementJet->Btag = btag;
106 }
107 } // for itJet : loop on all jets
108
109}
110
111void JetsUtil::RunTauJets(ExRootTreeWriter *treeWriter, ExRootTreeBranch *branchTauJet,const vector<fastjet::PseudoJet> & sorted_jets,const vector<PhysicsTower> & towers, const vector<TLorentzVector> & TrackCentral)
112{
113 TRootTauJet *elementTauJet;
114 TLorentzVector JET;
115 for (unsigned int i = 0; i < sorted_jets.size(); i++) {
116 JET.SetPxPyPzE(sorted_jets[i].px(),sorted_jets[i].py(),sorted_jets[i].pz(),sorted_jets[i].E());
117 // Tau jet identification : 1! track and electromagnetic collimation
118 if(fabs(JET.Eta()) < (MAX_TRACKER - TAU_CONE_TRACKS)) {
119 double Energie_tau_central = EnergySmallCone(towers,JET.Eta(),JET.Phi());
120 if(
121 ( Energie_tau_central/JET.E() > TAU_EM_COLLIMATION ) &&
122 ( NumTracks(TrackCentral,PT_TRACK_TAU,JET.Eta(),JET.Phi()) == 1 ) &&
123 ( JET.Pt() > TAUJET_pt)
124 ) {
125 elementTauJet = (TRootTauJet*) branchTauJet->NewEntry();
126 elementTauJet->Set(JET);
127 } // if tau jet
128 } // if JET.eta < tracker - tau_cone : Tau jet identification
129
130 } // for itJet : loop on all jets
131
132
133}
134
Note: See TracBrowser for help on using the repository browser.