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())
27  return _Rmax;
28 
29  std::vector<fastjet::PseudoJet> pieces = j.pieces();
30  if (pieces.size() != 2)
31  return _Rmax;
32 
33  double deltaR = pieces[0].delta_R(pieces[1]);
34  return std::min(_Rmax, _deltaR_factor * deltaR);
35  }
36 
37 private:
39 };
40 
42 public:
43  // typedefs
44  typedef fastjet::Transformer transformer;
45  typedef std::unique_ptr<transformer> transformer_ptr;
46  typedef std::vector<transformer_ptr> transformer_coll;
47  //
48  // construction/destruction
49  //
50  explicit FastjetJetProducer(const edm::ParameterSet& iConfig);
51  ~FastjetJetProducer() override;
52  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
53  static void fillDescriptionsFromFastJetProducer(edm::ParameterSetDescription& desc);
54 
55  void produce(edm::Event& iEvent, const edm::EventSetup& iSetup) override;
56 
57  // typedefs
58  typedef std::shared_ptr<DynamicRfilt> DynamicRfiltPtr;
59 
60 protected:
61  //
62  // member functions
63  //
64 
65  virtual void produceTrackJets(edm::Event& iEvent, const edm::EventSetup& iSetup);
66  void runAlgorithm(edm::Event& iEvent, const edm::EventSetup& iSetup) override;
67 
68 private:
69  // trackjet clustering parameters
72  float dzTrVtxMax_;
73  float dxyTrVtxMax_;
75  float maxVtxZ_;
76 
77  // jet trimming parameters
81  bool useTrimming_;
82  bool usePruning_;
86  bool useSoftDrop_;
88  double muCut_;
89  double yCut_;
90  double rFilt_;
91  DynamicRfiltPtr rFiltDynamic_;
92  double rFiltFactor_;
93  int nFilt_;
94  double trimPtFracMin_;
95  double zCut_;
96  double RcutFactor_;
97  double csRho_EtaMax_;
98  double csRParam_;
99  double beta_;
100  double R0_;
102  double gridSpacing_;
103 
104  double subjetPtMin_;
105  double muMin_;
106  double muMax_;
107  double yMin_;
108  double yMax_;
109  double dRMin_;
110  double dRMax_;
111  int maxDepth_;
112 
113  // tokens for the data access
115 };
116 
117 #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:224
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
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
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
std::shared_ptr< DynamicRfilt > DynamicRfiltPtr
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]