test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CATopJetAlgorithm.h
Go to the documentation of this file.
1 #ifndef RecoJets_JetAlgorithms_CATopJetAlgorithm_h
2 #define RecoJets_JetAlgorithms_CATopJetAlgorithm_h
3 
4 
5 
6 /* *********************************************************
7  * \class CATopJetAlgorithm
8  * Jet producer to produce top jets using the C-A algorithm to break
9  * jets into subjets as described here:
10  * "Top-tagging: A Method for Identifying Boosted Hadronic Tops"
11  * David E. Kaplan, Keith Rehermann, Matthew D. Schwartz, Brock Tweedie
12  * arXiv:0806.0848v1 [hep-ph]
13  *
14  ************************************************************/
15 
16 #include <vector>
17 #include <list>
18 #include <functional>
19 #include <TMath.h>
20 #include <iostream>
21 
22 #include "boost/shared_ptr.hpp"
23 
36 
38 
39 
40 #include <fastjet/JetDefinition.hh>
41 #include <fastjet/PseudoJet.hh>
42 #include <fastjet/ClusterSequence.hh>
43 #include <fastjet/GhostedAreaSpec.hh>
44 #include <fastjet/ClusterSequenceArea.hh>
45 
47  public:
51  bool verbose,
52  int algorithm,
53  int useAdjacency,
54  double centralEtaCut,
55  double ptMin,
56  const std::vector<double> & sumEtBins,
57  const std::vector<double> & rBins,
58  const std::vector<double> & ptFracBins,
59  const std::vector<double> & deltarBins,
60  const std::vector<double> & nCellBins,
61  double seedThreshold,
62  bool useMaxTower,
63  double sumEtEtaCut,
64  double etFrac) :
65  mSrc_ (mSrc ),
66  verbose_ (verbose ),
67  algorithm_ (algorithm ),
68  useAdjacency_ (useAdjacency ),
69  centralEtaCut_ (centralEtaCut ),
70  ptMin_ (ptMin ),
71  sumEtBins_ (sumEtBins ),
72  rBins_ (rBins ),
73  ptFracBins_ (ptFracBins ),
74  deltarBins_ (deltarBins ),
75  nCellBins_ (nCellBins ),
76  seedThreshold_ (seedThreshold ),
77  useMaxTower_ (useMaxTower ),
78  sumEtEtaCut_ (sumEtEtaCut ),
79  etFrac_ (etFrac )
80 
81  { }
82 
84  void run( const std::vector<fastjet::PseudoJet> & cell_particles,
85  std::vector<fastjet::PseudoJet> & hardjetsOutput ,
86  boost::shared_ptr<fastjet::ClusterSequence> & fjClusterSeq
87  );
88 
89  private:
90 
91  edm::InputTag mSrc_; //<! calo tower input source
92  bool verbose_; //<!
93  int algorithm_; //<! 0 = KT, 1 = CA, 2 = anti-KT
94  int useAdjacency_; //<! choose adjacency requirement:
95  //<! 0 = no adjacency
96  //<! 1 = deltar adjacency
97  //<! 2 = modified adjacency
98  //<! 3 = calotower neirest neigbor based adjacency (untested)
99  double centralEtaCut_; //<! eta for defining "central" jets
100  double ptMin_; //<! lower pt cut on which jets to reco
101  std::vector<double> sumEtBins_; //<! sumEt bins over which cuts vary. vector={bin 0 lower bound, bin 1 lower bound, ...}
102  std::vector<double> rBins_; //<! Jet distance paramter R. R values depend on sumEt bins.
103  std::vector<double> ptFracBins_; //<! deltap = fraction of full jet pt for a subjet to be consider "hard".
104  std::vector<double> deltarBins_; //<! Applicable only if useAdjacency=1. deltar adjacency values for each sumEtBin
105  std::vector<double> nCellBins_; //<! Applicable only if useAdjacency=3. number of cells apart for two subjets to be considered "independent"
106  // NOT USED:
107  double seedThreshold_; //<! calo tower seed threshold - NOT USED
108  bool useMaxTower_; //<! use max tower for jet adjacency criterion, false is to use the centroid - NOT USED
109  double sumEtEtaCut_; //<! eta for event SumEt - NOT USED
110  double etFrac_; //<! fraction of event sumEt / 2 for a jet to be considered "hard" - NOT USED
111  std::string jetType_; //<! CaloJets or GenJets - NOT USED
112 
113 
114  // Decide if the two jets are in adjacent cells
115  bool adjacentCells(const fastjet::PseudoJet & jet1, const fastjet::PseudoJet & jet2,
116  const std::vector<fastjet::PseudoJet> & cell_particles,
117  const fastjet::ClusterSequence & theClusterSequence,
118  double nCellMin ) const;
119 
120 
121  // Attempt to break up one "hard" jet into two "soft" jets
122 
123  bool decomposeJet(const fastjet::PseudoJet & theJet,
124  const fastjet::ClusterSequence & theClusterSequence,
125  const std::vector<fastjet::PseudoJet> & cell_particles,
126  double ptHard, double nCellMin, double deltarcut,
127  fastjet::PseudoJet & ja, fastjet::PseudoJet & jb,
128  std::vector<fastjet::PseudoJet> & leftovers) const;
129 
130 };
131 
132 
133 
134 #endif
std::vector< double > sumEtBins_
std::vector< double > rBins_
edm::InputTag mSrc_
bool adjacentCells(const fastjet::PseudoJet &jet1, const fastjet::PseudoJet &jet2, const std::vector< fastjet::PseudoJet > &cell_particles, const fastjet::ClusterSequence &theClusterSequence, double nCellMin) const
bool decomposeJet(const fastjet::PseudoJet &theJet, const fastjet::ClusterSequence &theClusterSequence, const std::vector< fastjet::PseudoJet > &cell_particles, double ptHard, double nCellMin, double deltarcut, fastjet::PseudoJet &ja, fastjet::PseudoJet &jb, std::vector< fastjet::PseudoJet > &leftovers) const
CATopJetAlgorithm(const edm::InputTag &mSrc, bool verbose, int algorithm, int useAdjacency, double centralEtaCut, double ptMin, const std::vector< double > &sumEtBins, const std::vector< double > &rBins, const std::vector< double > &ptFracBins, const std::vector< double > &deltarBins, const std::vector< double > &nCellBins, double seedThreshold, bool useMaxTower, double sumEtEtaCut, double etFrac)
std::vector< double > ptFracBins_
std::vector< double > deltarBins_
std::vector< double > nCellBins_
void run(const std::vector< fastjet::PseudoJet > &cell_particles, std::vector< fastjet::PseudoJet > &hardjetsOutput, boost::shared_ptr< fastjet::ClusterSequence > &fjClusterSeq)
Find the ProtoJets from the collection of input Candidates.