CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Types | Public Member Functions | Private Types | Private Member Functions | Private Attributes
heppy::ReclusterJets Class Reference

#include <ReclusterJets.h>

Public Types

typedef math::XYZTLorentzVector LorentzVector
 

Public Member Functions

std::vector< LorentzVectorgetGrouping (double ptMin=0.0)
 get grouping (inclusive jets) More...
 
std::vector< LorentzVectorgetGroupingExclusive (int njets)
 get grouping (exclusive jets, until n are left) More...
 
std::vector< LorentzVectorgetGroupingExclusive (double dcut)
 get grouping (exclusive jets, up to cut dcut) More...
 
LorentzVector getPruned (double zcut, double rcutFactor)
 get pruned 4-vector More...
 
LorentzVector getPrunedSubjetExclusive (unsigned int isubjet, double zcut, double rcutFactor)
 get pruned 4-vector for a given subject (must be called after getGroupingExclusive) More...
 
LorentzVector getPrunedSubjetInclusive (unsigned int isubjet, double zcut, double rcutFactor)
 get pruned 4-vector for a given subject (must be called after getGroupingInclusive) More...
 
 ReclusterJets (const std::vector< LorentzVector > &objects, double ktpower, double rparam)
 

Private Types

typedef boost::shared_ptr
< fastjet::ClusterSequence > 
ClusterSequencePtr
 fastjet outputs More...
 

Private Member Functions

LorentzVector getPruned (const fastjet::PseudoJet &jet, double zcut, double rcutFactor)
 
std::vector< LorentzVectormakeP4s (const std::vector< fastjet::PseudoJet > &jets)
 

Private Attributes

std::vector< fastjet::PseudoJet > exclusiveJets_
 
ClusterSequencePtr fjClusterSeq_
 
std::vector< fastjet::PseudoJet > fjInputs_
 
std::vector< fastjet::PseudoJet > inclusiveJets_
 
double ktpower_
 
double rparam_
 

Detailed Description

Definition at line 20 of file ReclusterJets.h.

Member Typedef Documentation

typedef boost::shared_ptr<fastjet::ClusterSequence> heppy::ReclusterJets::ClusterSequencePtr
private

fastjet outputs

Definition at line 60 of file ReclusterJets.h.

Definition at line 23 of file ReclusterJets.h.

Constructor & Destructor Documentation

heppy::ReclusterJets::ReclusterJets ( const std::vector< LorentzVector > &  objects,
double  ktpower,
double  rparam 
)

Definition at line 12 of file ReclusterJets.cc.

References Exception, fjClusterSeq_, fjInputs_, cmsHarvester::index, j, fwrapper::jet_def, ktpower_, python.connectstrParser::o, and rparam_.

12  :
13  ktpower_(ktpower), rparam_(rparam)
14 {
15  // define jet inputs
16  fjInputs_.clear();
17  int index=0;
18  for (const LorentzVector &o : objects) {
19  fastjet::PseudoJet j(o.Px(),o.Py(),o.Pz(),o.E());
20  j.set_user_index(index); 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 }
JetDefinition jet_def
std::vector< fastjet::PseudoJet > fjInputs_
Definition: ReclusterJets.h:54
math::XYZTLorentzVector LorentzVector
ClusterSequencePtr fjClusterSeq_
Definition: ReclusterJets.h:61
boost::shared_ptr< fastjet::ClusterSequence > ClusterSequencePtr
fastjet outputs
Definition: ReclusterJets.h:60
int j
Definition: DBlmapReader.cc:9

Member Function Documentation

std::vector< math::XYZTLorentzVector > heppy::ReclusterJets::getGrouping ( double  ptMin = 0.0)

get grouping (inclusive jets)

Definition at line 52 of file ReclusterJets.cc.

References fjClusterSeq_, inclusiveJets_, and makeP4s().

52  {
53  // recluster jet
54  inclusiveJets_ = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(ptMin));
55  // return
56  return makeP4s(inclusiveJets_);
57 }
ClusterSequencePtr fjClusterSeq_
Definition: ReclusterJets.h:61
std::vector< fastjet::PseudoJet > inclusiveJets_
Definition: ReclusterJets.h:62
std::vector< LorentzVector > makeP4s(const std::vector< fastjet::PseudoJet > &jets)
std::vector< math::XYZTLorentzVector > heppy::ReclusterJets::getGroupingExclusive ( int  njets)

get grouping (exclusive jets, until n are left)

Definition at line 66 of file ReclusterJets.cc.

References exclusiveJets_, fjClusterSeq_, and makeP4s().

66  {
67  // recluster jet
68  exclusiveJets_ = fastjet::sorted_by_pt(fjClusterSeq_->exclusive_jets(njets));
69  // return
70  return makeP4s(exclusiveJets_);
71 }
std::vector< fastjet::PseudoJet > exclusiveJets_
Definition: ReclusterJets.h:63
ClusterSequencePtr fjClusterSeq_
Definition: ReclusterJets.h:61
std::vector< LorentzVector > makeP4s(const std::vector< fastjet::PseudoJet > &jets)
std::vector< math::XYZTLorentzVector > heppy::ReclusterJets::getGroupingExclusive ( double  dcut)

get grouping (exclusive jets, up to cut dcut)

Definition at line 59 of file ReclusterJets.cc.

References exclusiveJets_, fjClusterSeq_, and makeP4s().

59  {
60  // recluster jet
61  exclusiveJets_ = fastjet::sorted_by_pt(fjClusterSeq_->exclusive_jets(dcut));
62  // return
63  return makeP4s(exclusiveJets_);
64 }
std::vector< fastjet::PseudoJet > exclusiveJets_
Definition: ReclusterJets.h:63
ClusterSequencePtr fjClusterSeq_
Definition: ReclusterJets.h:61
std::vector< LorentzVector > makeP4s(const std::vector< fastjet::PseudoJet > &jets)
math::XYZTLorentzVector heppy::ReclusterJets::getPruned ( double  zcut,
double  rcutFactor 
)

get pruned 4-vector

Definition at line 73 of file ReclusterJets.cc.

References exclusiveJets_, fjClusterSeq_, and getPrunedSubjetExclusive().

Referenced by getPrunedSubjetExclusive(), and getPrunedSubjetInclusive().

73  {
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 }
std::vector< fastjet::PseudoJet > exclusiveJets_
Definition: ReclusterJets.h:63
ClusterSequencePtr fjClusterSeq_
Definition: ReclusterJets.h:61
LorentzVector getPrunedSubjetExclusive(unsigned int isubjet, double zcut, double rcutFactor)
get pruned 4-vector for a given subject (must be called after getGroupingExclusive) ...
math::XYZTLorentzVector heppy::ReclusterJets::getPruned ( const fastjet::PseudoJet &  jet,
double  zcut,
double  rcutFactor 
)
private

Definition at line 93 of file ReclusterJets.cc.

93  {
94  // create pruner
95  fastjet::Pruner pruner(fastjet::cambridge_algorithm, zcut, rcutFactor);
96  // Prune the jet
97  fastjet::PseudoJet pruned_jet = pruner(jet);
98  return LorentzVector( pruned_jet.px(), pruned_jet.py(), pruned_jet.pz(), pruned_jet.e() );
99 }
math::XYZTLorentzVector LorentzVector
Definition: ReclusterJets.h:23
math::XYZTLorentzVector heppy::ReclusterJets::getPrunedSubjetExclusive ( unsigned int  isubjet,
double  zcut,
double  rcutFactor 
)

get pruned 4-vector for a given subject (must be called after getGroupingExclusive)

Definition at line 80 of file ReclusterJets.cc.

References Exception, exclusiveJets_, and getPruned().

Referenced by getPruned().

80  {
81  if (isubjet >= exclusiveJets_.size()) {
82  throw cms::Exception("InvalidArgument", "getPrunedSubjetExclusive called for non-existing exclusive subjet");
83  }
84  return getPruned(exclusiveJets_[isubjet], zcut, rcutFactor);
85 }
std::vector< fastjet::PseudoJet > exclusiveJets_
Definition: ReclusterJets.h:63
LorentzVector getPruned(double zcut, double rcutFactor)
get pruned 4-vector
math::XYZTLorentzVector heppy::ReclusterJets::getPrunedSubjetInclusive ( unsigned int  isubjet,
double  zcut,
double  rcutFactor 
)

get pruned 4-vector for a given subject (must be called after getGroupingInclusive)

Definition at line 86 of file ReclusterJets.cc.

References Exception, getPruned(), and inclusiveJets_.

86  {
87  if (isubjet >= inclusiveJets_.size()) {
88  throw cms::Exception("InvalidArgument", "getPrunedSubjetInclusive called for non-existing inclusive subjet");
89  }
90  return getPruned(inclusiveJets_[isubjet], zcut, rcutFactor);
91 }
LorentzVector getPruned(double zcut, double rcutFactor)
get pruned 4-vector
std::vector< fastjet::PseudoJet > inclusiveJets_
Definition: ReclusterJets.h:62
std::vector< math::XYZTLorentzVector > heppy::ReclusterJets::makeP4s ( const std::vector< fastjet::PseudoJet > &  jets)
private

Definition at line 45 of file ReclusterJets.cc.

Referenced by getGrouping(), and getGroupingExclusive().

45  {
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 }
math::XYZTLorentzVector LorentzVector
Definition: ReclusterJets.h:23
vector< PseudoJet > jets

Member Data Documentation

std::vector<fastjet::PseudoJet> heppy::ReclusterJets::exclusiveJets_
private

Definition at line 63 of file ReclusterJets.h.

Referenced by getGroupingExclusive(), getPruned(), and getPrunedSubjetExclusive().

ClusterSequencePtr heppy::ReclusterJets::fjClusterSeq_
private

Definition at line 61 of file ReclusterJets.h.

Referenced by getGrouping(), getGroupingExclusive(), getPruned(), and ReclusterJets().

std::vector<fastjet::PseudoJet> heppy::ReclusterJets::fjInputs_
private

Definition at line 54 of file ReclusterJets.h.

Referenced by ReclusterJets().

std::vector<fastjet::PseudoJet> heppy::ReclusterJets::inclusiveJets_
private

Definition at line 62 of file ReclusterJets.h.

Referenced by getGrouping(), and getPrunedSubjetInclusive().

double heppy::ReclusterJets::ktpower_
private

Definition at line 56 of file ReclusterJets.h.

Referenced by ReclusterJets().

double heppy::ReclusterJets::rparam_
private

Definition at line 57 of file ReclusterJets.h.

Referenced by ReclusterJets().