CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SubJetAlgorithm.cc
Go to the documentation of this file.
4 #include "fastjet/ClusterSequenceArea.hh"
5 
6 #include "boost/shared_ptr.hpp"
7 
8 using namespace std;
9 using namespace edm;
10 
12  zcut_ = z;
13 }
14 
16  rcut_factor_ = r;
17 }
18 
19 
20 // Run the algorithm
21 // ------------------
22 void SubJetAlgorithm::run( const vector<fastjet::PseudoJet> & cell_particles,
23  vector<CompoundPseudoJet> & hardjetsOutput ) {
24 
25  //for actual jet clustering, either the pruned or the original version is used.
26  //For the pruned version, a new jet definition using the PrunedRecombPlugin is required:
27  fastjet::FastPrunePlugin PRplugin(*fjJetDefinition_,
28  *fjJetDefinition_,
29  zcut_,
30  rcut_factor_);
31  fastjet::JetDefinition pjetdef(&PRplugin);
32 
33 
34 
35  // cluster the jets with the jet definition jetDef:
36  // run algorithm
37  boost::shared_ptr<fastjet::ClusterSequence> fjClusterSeq;
38  if ( !doAreaFastjet_ ) {
39  fjClusterSeq =
40  boost::shared_ptr<fastjet::ClusterSequence>( new fastjet::ClusterSequence( cell_particles,
41  pjetdef ) );
42  } else if (voronoiRfact_ <= 0) {
43  fjClusterSeq =
44  boost::shared_ptr<fastjet::ClusterSequence>( new fastjet::ClusterSequenceActiveArea( cell_particles,
45  pjetdef ,
46  *fjActiveArea_ ) );
47  } else {
48  fjClusterSeq =
49  boost::shared_ptr<fastjet::ClusterSequence>( new fastjet::ClusterSequenceVoronoiArea( cell_particles,
50  pjetdef ,
51  fastjet::VoronoiAreaSpec(voronoiRfact_) ) );
52  }
53 
54  vector<fastjet::PseudoJet> inclusiveJets = fjClusterSeq->inclusive_jets(ptMin_);
55 
56  // These will store the indices of each subjet that
57  // are present in each jet
58  vector<vector<int> > indices(inclusiveJets.size());
59  // Loop over inclusive jets, attempt to find substructure
60  vector<fastjet::PseudoJet>::iterator jetIt = inclusiveJets.begin();
61  for ( ; jetIt != inclusiveJets.end(); ++jetIt ) {
62  //decompose into requested number of subjets:
63  vector<fastjet::PseudoJet> subjets = fjClusterSeq->exclusive_subjets(*jetIt, nSubjets_);
64  //create the subjets objects to put into the "output" objects
65  vector<CompoundPseudoSubJet> subjetsOutput;
66  std::vector<fastjet::PseudoJet>::const_iterator itSubJetBegin = subjets.begin(),
67  itSubJet = itSubJetBegin, itSubJetEnd = subjets.end();
68  for (; itSubJet != itSubJetEnd; ++itSubJet ){
69  // Get the transient subjet constituents from fastjet
70  vector<fastjet::PseudoJet> subjetFastjetConstituents = fjClusterSeq->constituents( *itSubJet );
71  // Get the indices of the constituents:
72  vector<int> constituents;
73  vector<fastjet::PseudoJet>::const_iterator fastSubIt = subjetFastjetConstituents.begin(),
74  transConstEnd = subjetFastjetConstituents.end();
75  for ( ; fastSubIt != transConstEnd; ++fastSubIt ) {
76  if (fastSubIt->user_index() >= 0) {
77  constituents.push_back( fastSubIt->user_index() );
78  }
79  }
80 
81  double subJetArea = (doAreaFastjet_) ?
82  dynamic_cast<fastjet::ClusterSequenceActiveArea&>(*fjClusterSeq).area(*itSubJet) : 0.0;
83 
84 
85  // Make a CompoundPseudoSubJet object to hold this subjet and the indices of its constituents
86  subjetsOutput.push_back( CompoundPseudoSubJet( *itSubJet, subJetArea, constituents ) );
87 
88  }
89 
90 
91  double fatJetArea = (doAreaFastjet_) ?
92  dynamic_cast<fastjet::ClusterSequenceActiveArea&>(*fjClusterSeq).area(*jetIt) : 0.0;
93 
94  // Make a CompoundPseudoJet object to hold this hard jet, and the subjets that make it up
95  hardjetsOutput.push_back( CompoundPseudoJet( *jetIt,fatJetArea,subjetsOutput));
96  }
97 }
void set_rcut_factor(double r)
void set_zcut(double z)
CompoundPseudoJet holds an association of fastjet::PseudoJets that represent a &quot;hard&quot; top jet with su...
void run(const std::vector< fastjet::PseudoJet > &cell_particles, std::vector< CompoundPseudoJet > &hardjetsOutput)
Find the ProtoJets from the collection of input Candidates.