CMS 3D CMS Logo

FastjetJetProducer.h
Go to the documentation of this file.
1 #ifndef RecoJets_JetProducers_plugins_FastjetJetProducer_h
2 #define RecoJets_JetProducers_plugins_FastjetJetProducer_h
3 
7 
10 
12 
13 #include <fastjet/tools/Transformer.hh>
14 
15 // a function returning
16 // min(Rmax, deltaR_factor * deltaR(j1,j2))
17 // where j1 and j2 are the 2 subjets of j
18 // if the jet does not have exactly 2 pieces, Rmax is used.
19 class DynamicRfilt : public fastjet::FunctionOfPseudoJet<double>{
20 public:
21  // default ctor
22  DynamicRfilt(double Rmax, double deltaR_factor) : _Rmax(Rmax), _deltaR_factor(deltaR_factor){}
23 
24  // action of the function
25  double result(const fastjet::PseudoJet &j) const override{
26  if (! j.has_pieces()) return _Rmax;
27 
28  std::vector<fastjet::PseudoJet> pieces = j.pieces();
29  if (pieces.size() != 2) return _Rmax;
30 
31  double deltaR = pieces[0].delta_R(pieces[1]);
32  return std::min(_Rmax, _deltaR_factor * deltaR);
33  }
34 
35 private:
37 };
38 
39 
40 
42 {
43 
44 public:
45  // typedefs
46  typedef fastjet::Transformer transformer;
47  typedef std::unique_ptr<transformer> transformer_ptr;
48  typedef std::vector<transformer_ptr> transformer_coll;
49  //
50  // construction/destruction
51  //
52  explicit FastjetJetProducer(const edm::ParameterSet& iConfig);
53  ~FastjetJetProducer() override;
54  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
55  static void fillDescriptionsFromFastJetProducer(edm::ParameterSetDescription& desc);
56 
57  void produce( edm::Event & iEvent, const edm::EventSetup & iSetup ) override;
58 
59  // typedefs
60  typedef boost::shared_ptr<DynamicRfilt> DynamicRfiltPtr;
61 
62 protected:
63 
64  //
65  // member functions
66  //
67 
68  virtual void produceTrackJets( edm::Event & iEvent, const edm::EventSetup & iSetup );
69  void runAlgorithm( edm::Event& iEvent, const edm::EventSetup& iSetup ) override;
70 
71  private:
72 
73  // trackjet clustering parameters
76  float dzTrVtxMax_;
77  float dxyTrVtxMax_;
79  float maxVtxZ_;
80 
81  // jet trimming parameters
85  bool useTrimming_;
86  bool usePruning_;
90  bool useSoftDrop_;
92  double muCut_;
93  double yCut_;
94  double rFilt_;
95  DynamicRfiltPtr rFiltDynamic_;
96  double rFiltFactor_;
97  int nFilt_;
98  double trimPtFracMin_;
99  double zCut_;
100  double RcutFactor_;
101  double csRho_EtaMax_;
102  double csRParam_;
103  double beta_;
104  double R0_;
106  double gridSpacing_;
107 
108 
109  double subjetPtMin_;
110  double muMin_;
111  double muMax_;
112  double yMin_;
113  double yMax_;
114  double dRMin_;
115  double dRMax_;
116  int maxDepth_;
117 
118 
119  // tokens for the data access
121 
122 };
123 
124 
125 #endif
double muMin_
for CMSBoostedTauSeedingAlgorithm : subjet pt min
double yMin_
for CMSBoostedTauSeedingAlgorithm : max mass-drop
bool useFiltering_
Mass-drop tagging for boosted Higgs.
bool correctShape_
Soft drop.
double rFilt_
for mass-drop tagging, symmetry cut: min(pt1^2,pt2^2) * dR(1,2) / mjet > ycut
double dRMax_
for CMSBoostedTauSeedingAlgorithm : min dR
double subjetPtMin_
for shape subtraction, get the grid spacing
bool useDynamicFiltering_
Jet filtering technique.
double RcutFactor_
for pruning OR soft drop: constituent minimum pt fraction of parent cluster
bool useSoftDrop_
constituent subtraction technique
double yMax_
for CMSBoostedTauSeedingAlgorithm : min asymmetry
double gridSpacing_
for shape subtraction, get the fixed-grid rho
bool useKtPruning_
algorithm for seeding reconstruction of boosted Taus (similar to mass-drop tagging) ...
bool useTrimming_
Use dynamic filtering radius (as in arXiv:0802.2470)
double csRParam_
for constituent subtraction : maximum rapidity for ghosts
DynamicRfilt(double Rmax, double deltaR_factor)
std::vector< transformer_ptr > transformer_coll
double trimPtFracMin_
for filtering, pruning: number of subjets expected
double gridMaxRapidity_
for soft drop : R0 (angular distance normalization - should be set to jet radius in most cases) ...
int iEvent
Definition: GenABIO.cc:230
bool usePruning_
Jet trimming technique.
bool useCMSBoostedTauSeedingAlgorithm_
Jet pruning technique.
edm::EDGetTokenT< edm::View< reco::RecoChargedRefCandidate > > input_chrefcand_token_
for CMSBoostedTauSeedingAlgorithm : max depth for descending into clustering sequence ...
double csRho_EtaMax_
for pruning: constituent dR * pt/2m < rcut_factor
double muMax_
for CMSBoostedTauSeedingAlgorithm : min mass-drop
double muCut_
Correct the shape of the jets.
double rFiltFactor_
for dynamic filtering radius (as in arXiv:0802.2470)
DynamicRfiltPtr rFiltDynamic_
for filtering, trimming: dR scale of sub-clustering
double dRMin_
for CMSBoostedTauSeedingAlgorithm : max asymmetry
boost::shared_ptr< DynamicRfilt > DynamicRfiltPtr
T min(T a, T b)
Definition: MathUtil.h:58
double beta_
for constituent subtraction : R parameter for KT alg in jet median background estimator ...
double result(const fastjet::PseudoJet &j) const override
double zCut_
for trimming: constituent minimum pt fraction of full jet
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
int nFilt_
for dynamic filtering radius (as in arXiv:0802.2470)
int maxDepth_
for CMSBoostedTauSeedingAlgorithm : max dR
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
fastjet::Transformer transformer
double R0_
for soft drop : beta (angular exponent)
std::unique_ptr< transformer > transformer_ptr
double yCut_
for mass-drop tagging, m0/mjet (m0 = mass of highest mass subjet)
bool useConstituentSubtraction_
Use Kt clustering algorithm for pruning (default is Cambridge/Aachen)
const double Rmax[kNumberCalorimeter]