CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SubJetAlgorithm.h
Go to the documentation of this file.
1 #ifndef RecoJets_JetAlgorithms_CATopJetAlgorithm2_h
2 #define RecoJets_JetAlgorithms_CATopJetAlgorithm2_h
3 
4 #include <vector>
5 
6 #include "boost/shared_ptr.hpp"
7 
10 #include "RecoJets/JetAlgorithms/interface/FastPrunePlugin.hh"
11 #include <fastjet/JetDefinition.hh>
12 #include <fastjet/PseudoJet.hh>
13 #include <fastjet/ClusterSequence.hh>
14 #include <fastjet/GhostedAreaSpec.hh>
15 
16 
18  public:
20  unsigned int subjets,
21  double zcut,
22  double rcut_factor,
23  boost::shared_ptr<fastjet::JetDefinition> fjJetDefinition,
24  bool doAreaFastjet,
25  boost::shared_ptr<fastjet::GhostedAreaSpec> fjActiveArea,
26  double voronoiRfact
27  ) :
28  ptMin_ (ptMin ),
29  nSubjets_ (subjets ),
30  zcut_ (zcut ),
31  rcut_factor_ (rcut_factor ),
32  fjJetDefinition_(fjJetDefinition),
33  doAreaFastjet_ (doAreaFastjet),
34  fjActiveArea_ (fjActiveArea),
35  voronoiRfact_ (voronoiRfact)
36  {
37 
38  }
39 
40  void set_zcut(double z);
41  void set_rcut_factor(double r);
42  double zcut() const { return zcut_;}
43  double rcut_factor() const { return rcut_factor_; }
44 
46  void run( const std::vector<fastjet::PseudoJet> & cell_particles,
47  std::vector<CompoundPseudoJet> & hardjetsOutput);
48 
49 
50  private:
51 
52  double ptMin_; //<! lower pt cut on which jets to reco
53  int nSubjets_; //<! number of subjets to produce.
54  double zcut_; //<! zcut parameter (see arXiv:0903.5081). Only relevant if pruning is enabled.
55  double rcut_factor_; //<! r-cut factor (see arXiv:0903.5081).
56  boost::shared_ptr<fastjet::JetDefinition> fjJetDefinition_; //<! jet definition to use
57  bool doAreaFastjet_; //<! whether or not to use the fastjet area
58  boost::shared_ptr<fastjet::GhostedAreaSpec> fjActiveArea_; //<! fastjet area spec
59  double voronoiRfact_; //<! fastjet voronoi area R factor
60 };
61 
62 #endif
tuple voronoiRfact
double zcut() const
void set_rcut_factor(double r)
void set_zcut(double z)
tuple doAreaFastjet
SubJetAlgorithm(double ptMin, unsigned int subjets, double zcut, double rcut_factor, boost::shared_ptr< fastjet::JetDefinition > fjJetDefinition, bool doAreaFastjet, boost::shared_ptr< fastjet::GhostedAreaSpec > fjActiveArea, double voronoiRfact)
void run(const std::vector< fastjet::PseudoJet > &cell_particles, std::vector< CompoundPseudoJet > &hardjetsOutput)
Find the ProtoJets from the collection of input Candidates.
boost::shared_ptr< fastjet::GhostedAreaSpec > fjActiveArea_
double rcut_factor() const
boost::shared_ptr< fastjet::JetDefinition > fjJetDefinition_