CMS 3D CMS Logo

ReclusterJets.cc
Go to the documentation of this file.
3 #include "fastjet/tools/Pruner.hh"
4 
5 using namespace std;
6 
7 //using namespace std;
8 using namespace fastjet;
9 
10 namespace heppy {
11 
12  ReclusterJets::ReclusterJets(const std::vector<LorentzVector> &objects, double ktpower, double rparam)
13  : ktpower_(ktpower), rparam_(rparam) {
14  // define jet inputs
15  fjInputs_.clear();
16  int index = 0;
17  for (const LorentzVector &o : objects) {
18  fastjet::PseudoJet j(o.Px(), o.Py(), o.Pz(), o.E());
19  j.set_user_index(index);
20  index++; // in case we want to know which piece ended where
21  fjInputs_.push_back(j);
22  }
23 
24  // choose a jet definition
25  fastjet::JetDefinition jet_def;
26 
27  // prepare jet def
28  if (ktpower_ == 1.0) {
29  jet_def = JetDefinition(kt_algorithm, rparam_);
30  } else if (ktpower_ == 0.0) {
31  jet_def = JetDefinition(cambridge_algorithm, rparam_);
32  } else if (ktpower_ == -1.0) {
33  jet_def = JetDefinition(antikt_algorithm, rparam_);
34  } else {
35  throw cms::Exception("InvalidArgument", "Unsupported ktpower value");
36  }
37 
38  // print out some infos
39  // cout << "Clustering with " << jet_def.description() << endl;
41  // define jet clustering sequence
42  fjClusterSeq_ = ClusterSequencePtr(new fastjet::ClusterSequence(fjInputs_, jet_def));
43  }
44 
45  std::vector<math::XYZTLorentzVector> ReclusterJets::makeP4s(const std::vector<fastjet::PseudoJet> &jets) {
46  std::vector<math::XYZTLorentzVector> JetObjectsAll;
47  for (const fastjet::PseudoJet &pj : jets) {
48  JetObjectsAll.push_back(LorentzVector(pj.px(), pj.py(), pj.pz(), pj.e()));
49  }
50  return JetObjectsAll;
51  }
52  std::vector<math::XYZTLorentzVector> ReclusterJets::getGrouping(double ptMin) {
53  // recluster jet
54  inclusiveJets_ = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(ptMin));
55  // return
56  return makeP4s(inclusiveJets_);
57  }
58 
59  std::vector<math::XYZTLorentzVector> ReclusterJets::getGroupingExclusive(double dcut) {
60  // recluster jet
61  exclusiveJets_ = fastjet::sorted_by_pt(fjClusterSeq_->exclusive_jets(dcut));
62  // return
63  return makeP4s(exclusiveJets_);
64  }
65 
66  std::vector<math::XYZTLorentzVector> ReclusterJets::getGroupingExclusive(int njets) {
67  // recluster jet
68  exclusiveJets_ = fastjet::sorted_by_pt(fjClusterSeq_->exclusive_jets(njets));
69  // return
70  return makeP4s(exclusiveJets_);
71  }
72 
74  // cluster everything first
75  exclusiveJets_ = fastjet::sorted_by_pt(fjClusterSeq_->exclusive_jets(1));
76  // get pruned exclusive
77  return getPrunedSubjetExclusive(0, zcut, rcutFactor);
78  }
79 
81  double zcut,
82  double rcutFactor) {
83  if (isubjet >= exclusiveJets_.size()) {
84  throw cms::Exception("InvalidArgument", "getPrunedSubjetExclusive called for non-existing exclusive subjet");
85  }
86  return getPruned(exclusiveJets_[isubjet], zcut, rcutFactor);
87  }
89  double zcut,
90  double rcutFactor) {
91  if (isubjet >= inclusiveJets_.size()) {
92  throw cms::Exception("InvalidArgument", "getPrunedSubjetInclusive called for non-existing inclusive subjet");
93  }
94  return getPruned(inclusiveJets_[isubjet], zcut, rcutFactor);
95  }
96 
97  math::XYZTLorentzVector ReclusterJets::getPruned(const fastjet::PseudoJet &jet, double zcut, double rcutFactor) {
98  // create pruner
99  fastjet::Pruner pruner(fastjet::cambridge_algorithm, zcut, rcutFactor);
100  // Prune the jet
101  fastjet::PseudoJet pruned_jet = pruner(jet);
102  return LorentzVector(pruned_jet.px(), pruned_jet.py(), pruned_jet.pz(), pruned_jet.e());
103  }
104 
105 } // namespace heppy
JetDefinition jet_def
math::XYZTLorentzVector LorentzVector
Definition: ReclusterJets.h:21
std::vector< fastjet::PseudoJet > fjInputs_
Definition: ReclusterJets.h:51
std::vector< LorentzVector > getGroupingExclusive(int njets)
get grouping (exclusive jets, until n are left)
LorentzVector getPrunedSubjetInclusive(unsigned int isubjet, double zcut, double rcutFactor)
get pruned 4-vector for a given subject (must be called after getGroupingInclusive) ...
std::vector< LorentzVector > getGrouping(double ptMin=0.0)
get grouping (inclusive jets)
std::vector< fastjet::PseudoJet > exclusiveJets_
Definition: ReclusterJets.h:60
ClusterSequencePtr fjClusterSeq_
Definition: ReclusterJets.h:58
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
TAKEN FROM http://cmssw.cvs.cern.ch/cgi-bin/cmssw.cgi/CMSSW/ElectroWeakAnalysis/Utilities/src/PdfWeig...
Definition: AlphaT.h:16
LorentzVector getPruned(double zcut, double rcutFactor)
get pruned 4-vector
LorentzVector getPrunedSubjetExclusive(unsigned int isubjet, double zcut, double rcutFactor)
get pruned 4-vector for a given subject (must be called after getGroupingExclusive) ...
std::shared_ptr< fastjet::ClusterSequence > ClusterSequencePtr
fastjet outputs
Definition: ReclusterJets.h:57
std::vector< fastjet::PseudoJet > inclusiveJets_
Definition: ReclusterJets.h:59
std::vector< LorentzVector > makeP4s(const std::vector< fastjet::PseudoJet > &jets)