CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ReclusterJets.cc
Go to the documentation of this file.
3 
4 using namespace std;
5 
6 //using namespace std;
7 using namespace fastjet;
8 
9 namespace heppy{
10 
11 ReclusterJets::ReclusterJets(const std::vector<LorentzVector> & objects, double ktpower, double rparam) :
12  ktpower_(ktpower), rparam_(rparam)
13 {
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); index++; // in case we want to know which piece ended where
20  fjInputs_.push_back(j);
21  }
22 
23  // choose a jet definition
24  fastjet::JetDefinition jet_def;
25 
26  // prepare jet def
27  if (ktpower_ == 1.0) {
28  jet_def = JetDefinition(kt_algorithm, rparam_);
29  } else if (ktpower_ == 0.0) {
30  jet_def = JetDefinition(cambridge_algorithm, rparam_);
31  } else if (ktpower_ == -1.0) {
32  jet_def = JetDefinition(antikt_algorithm, rparam_);
33  } else {
34  throw cms::Exception("InvalidArgument", "Unsupported ktpower value");
35  }
36 
37  // print out some infos
38  // cout << "Clustering with " << jet_def.description() << endl;
40  // define jet clustering sequence
41  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequence( fjInputs_, jet_def));
42 }
43 
44 std::vector<math::XYZTLorentzVector> ReclusterJets::makeP4s(const std::vector<fastjet::PseudoJet> &jets) {
45  std::vector<math::XYZTLorentzVector> JetObjectsAll;
46  for (const fastjet::PseudoJet & pj : jets) {
47  JetObjectsAll.push_back( LorentzVector( pj.px(), pj.py(), pj.pz(), pj.e() ) );
48  }
49  return JetObjectsAll;
50 }
51 std::vector<math::XYZTLorentzVector> ReclusterJets::getGrouping(double ptMin) {
52  // recluster jet
53  inclusiveJets_ = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(ptMin));
54  // return
55  return makeP4s(inclusiveJets_);
56 }
57 
58 std::vector<math::XYZTLorentzVector> ReclusterJets::getGroupingExclusive(double dcut) {
59  // recluster jet
60  exclusiveJets_ = fastjet::sorted_by_pt(fjClusterSeq_->exclusive_jets(dcut));
61  // return
62  return makeP4s(exclusiveJets_);
63 }
64 
65 std::vector<math::XYZTLorentzVector> ReclusterJets::getGroupingExclusive(int njets) {
66  // recluster jet
67  exclusiveJets_ = fastjet::sorted_by_pt(fjClusterSeq_->exclusive_jets(njets));
68  // return
69  return makeP4s(exclusiveJets_);
70 }
71 }
JetDefinition jet_def
math::XYZTLorentzVector LorentzVector
Definition: ReclusterJets.h:23
std::vector< fastjet::PseudoJet > fjInputs_
Definition: ReclusterJets.h:40
std::vector< LorentzVector > getGroupingExclusive(int njets)
get grouping (exclusive jets, until n are left)
std::vector< LorentzVector > getGrouping(double ptMin=0.0)
get grouping (inclusive jets)
std::vector< fastjet::PseudoJet > exclusiveJets_
Definition: ReclusterJets.h:49
ClusterSequencePtr fjClusterSeq_
Definition: ReclusterJets.h:47
boost::shared_ptr< fastjet::ClusterSequence > ClusterSequencePtr
fastjet outputs
Definition: ReclusterJets.h:46
vector< PseudoJet > jets
int j
Definition: DBlmapReader.cc:9
std::vector< fastjet::PseudoJet > inclusiveJets_
Definition: ReclusterJets.h:48
std::vector< LorentzVector > makeP4s(const std::vector< fastjet::PseudoJet > &jets)