CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
JetClustering.cc
Go to the documentation of this file.
1 #include <vector>
2 #include <string>
3 
4 #include <fastjet/JetDefinition.hh>
5 #include <fastjet/PseudoJet.hh>
6 #include <fastjet/ClusterSequence.hh>
7 #include <fastjet/SISConePlugin.hh>
8 
11 
14 
15 namespace lhef {
16 
18  public:
21 
22  Algorithm(const edm::ParameterSet &params, double jetPtMin) :
23  jetPtMin(jetPtMin) {}
24  virtual ~Algorithm() {}
25 
26  virtual std::vector<Jet> operator () (
27  const ParticleVector &input) const = 0;
28 
29  double getJetPtMin() const { return jetPtMin; }
30 
31  private:
32  double jetPtMin;
33 };
34 
35 namespace {
36  class FastJetAlgorithmWrapper : public JetClustering::Algorithm {
37  public:
38  FastJetAlgorithmWrapper(const edm::ParameterSet &params,
39  double jetPtMin);
40  ~FastJetAlgorithmWrapper() {}
41 
42  protected:
43  std::vector<Jet> operator () (
44  const ParticleVector &input) const;
45 
46  std::auto_ptr<fastjet::JetDefinition::Plugin> plugin;
47  std::auto_ptr<fastjet::JetDefinition> jetDefinition;
48  };
49 
50  class KtAlgorithm : public FastJetAlgorithmWrapper {
51  public:
52  KtAlgorithm(const edm::ParameterSet &params,
53  double jetPtMin);
54  ~KtAlgorithm() {}
55  };
56 
57  class SISConeAlgorithm : public FastJetAlgorithmWrapper {
58  public:
59  SISConeAlgorithm(const edm::ParameterSet &params,
60  double jetPtMin);
61  ~SISConeAlgorithm() {}
62 
63  private:
64  static fastjet::SISConePlugin::SplitMergeScale
65  getScale(const std::string &name);
66  };
67 } // anonymous namespace
68 
69 FastJetAlgorithmWrapper::FastJetAlgorithmWrapper(
70  const edm::ParameterSet &params, double jetPtMin) :
71  JetClustering::Algorithm(params, jetPtMin)
72 {
73 }
74 
75 std::vector<JetClustering::Jet> FastJetAlgorithmWrapper::operator () (
76  const ParticleVector &input) const
77 {
78  if (input.empty())
79  return std::vector<JetClustering::Jet>();
80 
81  std::vector<fastjet::PseudoJet> jfInput;
82  jfInput.reserve(input.size());
83  for(ParticleVector::const_iterator iter = input.begin();
84  iter != input.end(); ++iter) {
85  jfInput.push_back(fastjet::PseudoJet(
86  (*iter)->momentum().px(), (*iter)->momentum().py(),
87  (*iter)->momentum().pz(), (*iter)->momentum().e()));
88  jfInput.back().set_user_index(iter - input.begin());
89  }
90 
91  fastjet::ClusterSequence sequence(jfInput, *jetDefinition);
92  std::vector<fastjet::PseudoJet> jets =
93  sequence.inclusive_jets(getJetPtMin());
94 
95  std::vector<Jet> result;
96  result.reserve(jets.size());
97  ParticleVector constituents;
98  for(std::vector<fastjet::PseudoJet>::const_iterator iter = jets.begin();
99  iter != jets.end(); ++iter) {
100  std::vector<fastjet::PseudoJet> fjConstituents =
101  sequence.constituents(*iter);
102  unsigned int size = fjConstituents.size();
103  constituents.resize(size);
104  for(unsigned int i = 0; i < size; i++)
105  constituents[i] =
106  input[fjConstituents[i].user_index()];
107 
108  result.push_back(
109  Jet(iter->px(), iter->py(), iter->pz(), iter->E(),
110  constituents));
111  }
112 
113  return result;
114 }
115 
116 KtAlgorithm::KtAlgorithm(const edm::ParameterSet &params, double jetPtMin) :
117  FastJetAlgorithmWrapper(params, jetPtMin)
118 {
119  jetDefinition.reset(new fastjet::JetDefinition(
120  fastjet::kt_algorithm,
121  params.getParameter<double>("ktRParam"),
122  fastjet::Best));
123 }
124 
125 SISConeAlgorithm::SISConeAlgorithm(
126  const edm::ParameterSet &params, double jetPtMin) :
127  FastJetAlgorithmWrapper(params, jetPtMin)
128 {
129  std::string splitMergeScale =
130  params.getParameter<std::string>("splitMergeScale");
131  fastjet::SISConePlugin::SplitMergeScale scale;
132 
133  if (splitMergeScale == "pt")
134  scale = fastjet::SISConePlugin::SM_pt;
135  else if (splitMergeScale == "Et")
136  scale = fastjet::SISConePlugin::SM_Et;
137  else if (splitMergeScale == "mt")
138  scale = fastjet::SISConePlugin::SM_mt;
139  else if (splitMergeScale == "pttilde")
140  scale = fastjet::SISConePlugin::SM_pttilde;
141  else
142  throw cms::Exception("Configuration")
143  << "JetClustering SISCone scale '" << splitMergeScale
144  << "' unknown." << std::endl;
145 
146  plugin.reset(new fastjet::SISConePlugin(
147  params.getParameter<double>("coneRadius"),
148  params.getParameter<double>("coneOverlapThreshold"),
149  params.getParameter<int>("maxPasses"),
150  params.getParameter<double>("protojetPtMin"),
151  params.getParameter<bool>("caching"),
152  scale));
153  jetDefinition.reset(new fastjet::JetDefinition(plugin.get()));
154 }
155 
157 {
158  double jetPtMin = params.getParameter<double>("jetPtMin");
159  init(params, jetPtMin);
160 }
161 
163  double jetPtMin)
164 {
165  init(params, jetPtMin);
166 }
167 
169 {
170 }
171 
172 void JetClustering::init(const edm::ParameterSet &params, double jetPtMin)
173 {
174  edm::ParameterSet algoParams =
175  params.getParameter<edm::ParameterSet>("algorithm");
176  std::string algoName = algoParams.getParameter<std::string>("name");
177 
178  if (algoName == "KT")
179  jetAlgo.reset(new KtAlgorithm(algoParams, jetPtMin));
180  else if (algoName == "SISCone")
181  jetAlgo.reset(new SISConeAlgorithm(algoParams, jetPtMin));
182  else
183  throw cms::Exception("Configuration")
184  << "JetClustering algorithm \"" << algoName
185  << "\" unknown." << std::endl;
186 }
187 
189 {
190  return jetAlgo->getJetPtMin();
191 }
192 
193 std::vector<JetClustering::Jet> JetClustering::operator () (
194  const ParticleVector &input) const
195 {
196  return (*jetAlgo)(input);
197 }
198 
199 } // namespace lhef
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
double getJetPtMin() const
std::auto_ptr< Algorithm > jetAlgo
Definition: JetClustering.h:70
JetClustering::ParticleVector ParticleVector
JetClustering(const edm::ParameterSet &params)
tuple plugin
Definition: pfTaus_cff.py:103
vector< PseudoJet > jets
tuple result
Definition: query.py:137
std::vector< Jet > operator()(const ParticleVector &input) const
JetInput::ParticleVector ParticleVector
Definition: JetClustering.h:17
virtual std::vector< Jet > operator()(const ParticleVector &input) const =0
void init(const edm::ParameterSet &params, double jetPtMin)
Algorithm(const edm::ParameterSet &params, double jetPtMin)
tuple size
Write out results.