CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Private Attributes
FastjetJetProducer Class Reference

#include <FastjetJetProducer.h>

Inheritance diagram for FastjetJetProducer:
VirtualJetProducer edm::stream::EDProducer<> cms::CATopJetProducer cms::HTTTopJetProducer

Public Types

typedef std::shared_ptr< DynamicRfiltDynamicRfiltPtr
 
typedef fastjet::Transformer transformer
 
typedef std::vector< transformer_ptrtransformer_coll
 
typedef std::unique_ptr< transformertransformer_ptr
 
- Public Types inherited from VirtualJetProducer
typedef std::shared_ptr< fastjet::GhostedAreaSpec > ActiveAreaSpecPtr
 
typedef std::shared_ptr< fastjet::AreaDefinition > AreaDefinitionPtr
 
typedef std::shared_ptr< fastjet::ClusterSequence > ClusterSequencePtr
 
typedef std::shared_ptr< fastjet::JetDefinition > JetDefPtr
 
typedef std::shared_ptr< fastjet::JetDefinition::Plugin > PluginPtr
 
typedef std::shared_ptr< fastjet::Selector > SelectorPtr
 
- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Public Member Functions

 FastjetJetProducer (const edm::ParameterSet &iConfig)
 
void produce (edm::Event &iEvent, const edm::EventSetup &iSetup) override
 
 ~FastjetJetProducer () override
 
- Public Member Functions inherited from VirtualJetProducer
std::string jetType () const
 
void produce (edm::Event &iEvent, const edm::EventSetup &iSetup) override
 
 VirtualJetProducer (const edm::ParameterSet &iConfig)
 
 ~VirtualJetProducer () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 
static void fillDescriptionsFromFastJetProducer (edm::ParameterSetDescription &desc)
 
- Static Public Member Functions inherited from VirtualJetProducer
static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 
static void fillDescriptionsFromVirtualJetProducer (edm::ParameterSetDescription &desc)
 

Protected Member Functions

virtual void produceTrackJets (edm::Event &iEvent, const edm::EventSetup &iSetup)
 
void runAlgorithm (edm::Event &iEvent, const edm::EventSetup &iSetup) override
 
- Protected Member Functions inherited from VirtualJetProducer
virtual void addHTTTopJetTagInfoCollection (edm::Event &iEvent, const edm::EventSetup &iSetup, edm::OrphanHandle< reco::BasicJetCollection > &oh)
 
virtual void copyConstituents (const std::vector< fastjet::PseudoJet > &fjConstituents, reco::Jet *jet)
 
virtual std::vector< reco::CandidatePtrgetConstituents (const std::vector< fastjet::PseudoJet > &fjConstituents)
 
CaloGeometry const & getGeometry (edm::EventSetup const &) const
 
HcalTopology const & getTopology (edm::EventSetup const &) const
 
virtual void inputTowers ()
 
virtual bool isAnomalousTower (reco::CandidatePtr input)
 
bool makeBasicJet (const JetType::Type &fTag)
 
bool makeCaloJet (const JetType::Type &fTag)
 
bool makeGenJet (const JetType::Type &fTag)
 
bool makePFClusterJet (const JetType::Type &fTag)
 
bool makePFJet (const JetType::Type &fTag)
 
virtual void makeProduces (std::string s, std::string tag="")
 
bool makeTrackJet (const JetType::Type &fTag)
 
void offsetCorrectJets (std::vector< fastjet::PseudoJet > &orphanInput)
 
virtual void output (edm::Event &iEvent, edm::EventSetup const &iSetup)
 
template<typename T >
void writeCompoundJets (edm::Event &iEvent, edm::EventSetup const &iSetup)
 function template to write out the outputs More...
 
template<typename T >
void writeJets (edm::Event &iEvent, edm::EventSetup const &iSetup)
 
template<typename T >
void writeJetsWithConstituents (edm::Event &iEvent, edm::EventSetup const &iSetup)
 function template to write out the outputs More...
 

Private Attributes

double beta_
 for constituent subtraction : R parameter for KT alg in jet median background estimator More...
 
bool correctShape_
 Soft drop. More...
 
double csRho_EtaMax_
 for pruning: constituent dR * pt/2m < rcut_factor More...
 
double csRParam_
 for constituent subtraction : maximum rapidity for ghosts More...
 
double dRMax_
 for CMSBoostedTauSeedingAlgorithm : min dR More...
 
double dRMin_
 for CMSBoostedTauSeedingAlgorithm : max asymmetry More...
 
float dxyTrVtxMax_
 
float dzTrVtxMax_
 
double gridMaxRapidity_
 for soft drop : R0 (angular distance normalization - should be set to jet radius in most cases) More...
 
double gridSpacing_
 for shape subtraction, get the fixed-grid rho More...
 
edm::EDGetTokenT< edm::View< reco::RecoChargedRefCandidate > > input_chrefcand_token_
 for CMSBoostedTauSeedingAlgorithm : max depth for descending into clustering sequence More...
 
int maxDepth_
 for CMSBoostedTauSeedingAlgorithm : max dR More...
 
float maxVtxZ_
 
int minVtxNdof_
 
double muCut_
 Correct the shape of the jets. More...
 
double muMax_
 for CMSBoostedTauSeedingAlgorithm : min mass-drop More...
 
double muMin_
 for CMSBoostedTauSeedingAlgorithm : subjet pt min More...
 
int nFilt_
 for dynamic filtering radius (as in arXiv:0802.2470) More...
 
double R0_
 for soft drop : beta (angular exponent) More...
 
double RcutFactor_
 for pruning OR soft drop: constituent minimum pt fraction of parent cluster More...
 
double rFilt_
 for mass-drop tagging, symmetry cut: min(pt1^2,pt2^2) * dR(1,2) / mjet > ycut More...
 
DynamicRfiltPtr rFiltDynamic_
 for filtering, trimming: dR scale of sub-clustering More...
 
double rFiltFactor_
 for dynamic filtering radius (as in arXiv:0802.2470) More...
 
double subjetPtMin_
 for shape subtraction, get the grid spacing More...
 
double trimPtFracMin_
 for filtering, pruning: number of subjets expected More...
 
bool useCMSBoostedTauSeedingAlgorithm_
 Jet pruning technique. More...
 
bool useConstituentSubtraction_
 Use Kt clustering algorithm for pruning (default is Cambridge/Aachen) More...
 
bool useDynamicFiltering_
 Jet filtering technique. More...
 
bool useFiltering_
 Mass-drop tagging for boosted Higgs. More...
 
bool useKtPruning_
 algorithm for seeding reconstruction of boosted Taus (similar to mass-drop tagging) More...
 
bool useMassDropTagger_
 
bool useOnlyOnePV_
 
bool useOnlyVertexTracks_
 
bool usePruning_
 Jet trimming technique. More...
 
bool useSoftDrop_
 constituent subtraction technique More...
 
bool useTrimming_
 Use dynamic filtering radius (as in arXiv:0802.2470) More...
 
double yCut_
 for mass-drop tagging, m0/mjet (m0 = mass of highest mass subjet) More...
 
double yMax_
 for CMSBoostedTauSeedingAlgorithm : min asymmetry More...
 
double yMin_
 for CMSBoostedTauSeedingAlgorithm : max mass-drop More...
 
double zCut_
 for trimming: constituent minimum pt fraction of full jet More...
 

Additional Inherited Members

- Protected Attributes inherited from VirtualJetProducer
int activeAreaRepeats_
 
bool applyWeight_
 
bool doAreaDiskApprox_
 
bool doAreaFastjet_
 
bool doFastJetNonUniform_
 
bool doPUOffsetCorr_
 
bool doPVCorrection_
 
bool doRhoFastjet_
 
ActiveAreaSpecPtr fjActiveArea_
 
AreaDefinitionPtr fjAreaDefinition_
 
ClusterSequencePtr fjClusterSeq_
 
std::vector< fastjet::PseudoJet > fjInputs_
 
JetDefPtr fjJetDefinition_
 
std::vector< fastjet::PseudoJet > fjJets_
 
PluginPtr fjPlugin_
 
SelectorPtr fjSelector_
 
bool fromHTTTopJetProducer_ = false
 
double ghostArea_
 
double ghostEtaMax_
 
edm::EDGetTokenT< reco::VertexCollectioninput_vertex_token_
 
edm::EDGetTokenT< edm::ValueMap< float > > input_weights_token_
 
double inputEMin_
 
double inputEtMin_
 
std::vector< edm::Ptr< reco::Candidate > > inputs_
 
std::string jetAlgorithm_
 
std::string jetCollInstanceName_
 
double jetPtMin_
 
std::string jetType_
 
JetType::Type jetTypeE
 
unsigned int maxInputs_
 
unsigned int minSeed_
 
std::string moduleLabel_
 
unsigned int nExclude_
 
std::vector< double > puCenters_
 
std::string puSubtractorName_
 
double puWidth_
 
bool restrictInputs_
 
double rhoEtaMax_
 
double rParam_
 
edm::InputTag src_
 
edm::InputTag srcPVs_
 
std::shared_ptr< PileUpSubtractorsubtractor_
 
bool useDeterministicSeed_
 
bool useExplicitGhosts_
 
int verbosity_
 
reco::Particle::Point vertex_
 
double voronoiRfact_
 
edm::ValueMap< float > weights_
 
bool writeCompound_
 
bool writeJetsWithConst_
 

Detailed Description

Definition at line 41 of file FastjetJetProducer.h.

Member Typedef Documentation

◆ DynamicRfiltPtr

Definition at line 58 of file FastjetJetProducer.h.

◆ transformer

typedef fastjet::Transformer FastjetJetProducer::transformer

Definition at line 44 of file FastjetJetProducer.h.

◆ transformer_coll

Definition at line 46 of file FastjetJetProducer.h.

◆ transformer_ptr

Definition at line 45 of file FastjetJetProducer.h.

Constructor & Destructor Documentation

◆ FastjetJetProducer()

FastjetJetProducer::FastjetJetProducer ( const edm::ParameterSet iConfig)
explicit

Definition at line 60 of file FastjetJetProducer.cc.

References beta_, correctShape_, csRho_EtaMax_, csRParam_, dRMax_, dRMin_, dxyTrVtxMax_, dzTrVtxMax_, Exception, VirtualJetProducer::fjAreaDefinition_, edm::ParameterSet::getParameter(), gridMaxRapidity_, gridSpacing_, input_chrefcand_token_, maxDepth_, maxVtxZ_, minVtxNdof_, muCut_, muMax_, muMin_, nFilt_, R0_, RcutFactor_, rFilt_, rFiltDynamic_, rFiltFactor_, subjetPtMin_, trimPtFracMin_, useCMSBoostedTauSeedingAlgorithm_, useConstituentSubtraction_, useDynamicFiltering_, VirtualJetProducer::useExplicitGhosts_, useFiltering_, useKtPruning_, useMassDropTagger_, useOnlyOnePV_, useOnlyVertexTracks_, usePruning_, useSoftDrop_, useTrimming_, yCut_, yMax_, yMin_, and zCut_.

60  : VirtualJetProducer(iConfig) {
61  useOnlyVertexTracks_ = iConfig.getParameter<bool>("UseOnlyVertexTracks");
62  useOnlyOnePV_ = iConfig.getParameter<bool>("UseOnlyOnePV");
63  dzTrVtxMax_ = iConfig.getParameter<double>("DzTrVtxMax");
64  dxyTrVtxMax_ = iConfig.getParameter<double>("DxyTrVtxMax");
65  minVtxNdof_ = iConfig.getParameter<int>("MinVtxNdof");
66  maxVtxZ_ = iConfig.getParameter<double>("MaxVtxZ");
67 
68  useMassDropTagger_ = iConfig.getParameter<bool>("useMassDropTagger");
69  muCut_ = iConfig.getParameter<double>("muCut");
70  yCut_ = iConfig.getParameter<double>("yCut");
71 
72  useFiltering_ = iConfig.getParameter<bool>("useFiltering");
73  rFilt_ = iConfig.getParameter<double>("rFilt");
74  nFilt_ = iConfig.getParameter<int>("nFilt");
75  useDynamicFiltering_ = iConfig.getParameter<bool>("useDynamicFiltering");
77  rFiltDynamic_ = std::make_shared<DynamicRfilt>(rFilt_, rFiltFactor_);
78  rFiltFactor_ = iConfig.getParameter<double>("rFiltFactor");
79 
80  useTrimming_ = iConfig.getParameter<bool>("useTrimming");
81  trimPtFracMin_ = iConfig.getParameter<double>("trimPtFracMin");
82 
83  usePruning_ = iConfig.getParameter<bool>("usePruning");
84  zCut_ = iConfig.getParameter<double>("zcut");
85  RcutFactor_ = iConfig.getParameter<double>("rcut_factor");
86  useKtPruning_ = iConfig.getParameter<bool>("useKtPruning");
87 
88  useCMSBoostedTauSeedingAlgorithm_ = iConfig.getParameter<bool>("useCMSBoostedTauSeedingAlgorithm");
89  subjetPtMin_ = iConfig.getParameter<double>("subjetPtMin");
90  muMin_ = iConfig.getParameter<double>("muMin");
91  muMax_ = iConfig.getParameter<double>("muMax");
92  yMin_ = iConfig.getParameter<double>("yMin");
93  yMax_ = iConfig.getParameter<double>("yMax");
94  dRMin_ = iConfig.getParameter<double>("dRMin");
95  dRMax_ = iConfig.getParameter<double>("dRMax");
96  maxDepth_ = iConfig.getParameter<int>("maxDepth");
97 
98  useConstituentSubtraction_ = iConfig.getParameter<bool>("useConstituentSubtraction");
99  csRho_EtaMax_ = iConfig.getParameter<double>("csRho_EtaMax");
100  csRParam_ = iConfig.getParameter<double>("csRParam");
101 
102  useSoftDrop_ = iConfig.getParameter<bool>("useSoftDrop");
103  zCut_ = iConfig.getParameter<double>("zcut");
104  beta_ = iConfig.getParameter<double>("beta");
105  R0_ = iConfig.getParameter<double>("R0");
106 
107  correctShape_ = iConfig.getParameter<bool>("correctShape");
108  gridMaxRapidity_ = iConfig.getParameter<double>("gridMaxRapidity");
109  gridSpacing_ = iConfig.getParameter<double>("gridSpacing");
110 
112  consumes<edm::View<reco::RecoChargedRefCandidate>>(iConfig.getParameter<edm::InputTag>("src"));
113 
116  useExplicitGhosts_ = true;
117 
119 
120  if ((useMassDropTagger_) && ((muCut_ == -1) || (yCut_ == -1)))
121  throw cms::Exception("useMassDropTagger")
122  << "Parameters muCut and/or yCut for Mass Drop are not defined." << std::endl;
123 
124  if ((useFiltering_) && ((rFilt_ == -1) || (nFilt_ == -1))) {
125  throw cms::Exception("useFiltering") << "Parameters rFilt and/or nFilt for Filtering are not defined." << std::endl;
126  if ((useDynamicFiltering_) && (rFiltFactor_ == -1))
127  throw cms::Exception("useDynamicFiltering")
128  << "Parameters rFiltFactor for DynamicFiltering is not defined." << std::endl;
129  }
130 
131  if ((useTrimming_) && ((rFilt_ == -1) || (trimPtFracMin_ == -1)))
132  throw cms::Exception("useTrimming") << "Parameters rFilt and/or trimPtFracMin for Trimming are not defined."
133  << std::endl;
134 
135  if ((usePruning_) && ((zCut_ == -1) || (RcutFactor_ == -1) || (nFilt_ == -1)))
136  throw cms::Exception("usePruning") << "Parameters zCut and/or RcutFactor and/or nFilt for Pruning are not defined."
137  << std::endl;
138 
140  ((subjetPtMin_ == -1) || (maxDepth_ == -1) || (muMin_ == -1) || (muMax_ == -1) || (yMin_ == -1) ||
141  (yMax_ == -1) || (dRMin_ == -1) || (dRMax_ == -1)))
142  throw cms::Exception("useCMSBoostedTauSeedingAlgorithm")
143  << "Parameters subjetPtMin, muMin, muMax, yMin, yMax, dRmin, dRmax, maxDepth for CMSBoostedTauSeedingAlgorithm "
144  "are not defined."
145  << std::endl;
146 
147  if (useConstituentSubtraction_ && (fjAreaDefinition_.get() == nullptr))
148  throw cms::Exception("AreaMustBeSet")
149  << "Logic error. The area definition must be set if you use constituent subtraction." << std::endl;
150 
151  if ((useConstituentSubtraction_) && ((csRho_EtaMax_ == -1) || (csRParam_ == -1)))
152  throw cms::Exception("useConstituentSubtraction")
153  << "Parameters csRho_EtaMax and/or csRParam for ConstituentSubtraction are not defined." << std::endl;
154 
155  if (useSoftDrop_ && usePruning_)
156  throw cms::Exception("PruningAndSoftDrop")
157  << "Logic error. Soft drop is a generalized pruning, do not run them together." << std::endl;
158 
159  if ((useSoftDrop_) && ((zCut_ == -1) || (beta_ == -1) || (R0_ == -1)))
160  throw cms::Exception("useSoftDrop") << "Parameters zCut and/or beta and/or R0 for SoftDrop are not defined."
161  << std::endl;
162 
163  if ((correctShape_) && ((gridMaxRapidity_ == -1) || (gridSpacing_ == -1)))
164  throw cms::Exception("correctShape")
165  << "Parameters gridMaxRapidity and/or gridSpacing for SoftDrop are not defined." << std::endl;
166 }
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
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
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
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) ...
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
double beta_
for constituent subtraction : R parameter for KT alg in jet median background estimator ...
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
double R0_
for soft drop : beta (angular exponent)
double yCut_
for mass-drop tagging, m0/mjet (m0 = mass of highest mass subjet)
VirtualJetProducer(const edm::ParameterSet &iConfig)
bool useConstituentSubtraction_
Use Kt clustering algorithm for pruning (default is Cambridge/Aachen)
AreaDefinitionPtr fjAreaDefinition_

◆ ~FastjetJetProducer()

FastjetJetProducer::~FastjetJetProducer ( )
override

Definition at line 169 of file FastjetJetProducer.cc.

169 {}

Member Function Documentation

◆ fillDescriptions()

void FastjetJetProducer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 464 of file FastjetJetProducer.cc.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), fillDescriptionsFromFastJetProducer(), and VirtualJetProducer::fillDescriptionsFromVirtualJetProducer().

464  {
465  edm::ParameterSetDescription descFastjetJetProducer;
467  fillDescriptionsFromFastJetProducer(descFastjetJetProducer);
471  descFastjetJetProducer.add<string>("jetCollInstanceName", "");
472  descFastjetJetProducer.add<bool>("sumRecHits", false);
473 
475  descriptions.add("FastjetJetProducer", descFastjetJetProducer);
476 }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillDescriptionsFromFastJetProducer(edm::ParameterSetDescription &desc)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static void fillDescriptionsFromVirtualJetProducer(edm::ParameterSetDescription &desc)

◆ fillDescriptionsFromFastJetProducer()

void FastjetJetProducer::fillDescriptionsFromFastJetProducer ( edm::ParameterSetDescription desc)
static

Definition at line 478 of file FastjetJetProducer.cc.

References submitPVResolutionJobs::desc.

Referenced by fillDescriptions(), cms::CATopJetProducer::fillDescriptions(), and cms::HTTTopJetProducer::fillDescriptions().

478  {
479  desc.add<bool>("useMassDropTagger", false);
480  desc.add<bool>("useFiltering", false);
481  desc.add<bool>("useDynamicFiltering", false);
482  desc.add<bool>("useTrimming", false);
483  desc.add<bool>("usePruning", false);
484  desc.add<bool>("useCMSBoostedTauSeedingAlgorithm", false);
485  desc.add<bool>("useKtPruning", false);
486  desc.add<bool>("useConstituentSubtraction", false);
487  desc.add<bool>("useSoftDrop", false);
488  desc.add<bool>("correctShape", false);
489  desc.add<bool>("UseOnlyVertexTracks", false);
490  desc.add<bool>("UseOnlyOnePV", false);
491  desc.add<double>("muCut", -1.0);
492  desc.add<double>("yCut", -1.0);
493  desc.add<double>("rFilt", -1.0);
494  desc.add<double>("rFiltFactor", -1.0);
495  desc.add<double>("trimPtFracMin", -1.0);
496  desc.add<double>("zcut", -1.0);
497  desc.add<double>("rcut_factor", -1.0);
498  desc.add<double>("csRho_EtaMax", -1.0);
499  desc.add<double>("csRParam", -1.0);
500  desc.add<double>("beta", -1.0);
501  desc.add<double>("R0", -1.0);
502  desc.add<double>("gridMaxRapidity", -1.0); // For fixed-grid rho
503  desc.add<double>("gridSpacing", -1.0); // For fixed-grid rho
504  desc.add<double>("DzTrVtxMax", 999999.);
505  desc.add<double>("DxyTrVtxMax", 999999.);
506  desc.add<double>("MaxVtxZ", 15.0);
507  desc.add<double>("subjetPtMin", -1.0);
508  desc.add<double>("muMin", -1.0);
509  desc.add<double>("muMax", -1.0);
510  desc.add<double>("yMin", -1.0);
511  desc.add<double>("yMax", -1.0);
512  desc.add<double>("dRMin", -1.0);
513  desc.add<double>("dRMax", -1.0);
514  desc.add<int>("maxDepth", -1);
515  desc.add<int>("nFilt", -1);
516  desc.add<int>("MinVtxNdof", 5);
517 }

◆ produce()

void FastjetJetProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 175 of file FastjetJetProducer.cc.

References VirtualJetProducer::fjClusterSeq_, iEvent, VirtualJetProducer::jetTypeE, VirtualJetProducer::makeTrackJet(), VirtualJetProducer::produce(), and produceTrackJets().

Referenced by cms::CATopJetProducer::produce(), and cms::HTTTopJetProducer::produce().

175  {
176  // for everything but track jets
177  if (!makeTrackJet(jetTypeE)) {
178  // use the default production from one collection
180 
181  } else { // produce trackjets from tracks grouped per primary vertex
182 
183  produceTrackJets(iEvent, iSetup);
184  }
185 
186  // fjClusterSeq_ retains quite a lot of memory - about 1 to 7Mb at 200 pileup
187  // depending on the exact configuration; and there are 24 FastjetJetProducers in the
188  // sequence so this adds up to about 60 Mb. It's allocated every time runAlgorithm
189  // is called, so safe to delete here.
190  fjClusterSeq_.reset();
191 }
JetType::Type jetTypeE
bool makeTrackJet(const JetType::Type &fTag)
int iEvent
Definition: GenABIO.cc:224
ClusterSequencePtr fjClusterSeq_
virtual void produceTrackJets(edm::Event &iEvent, const edm::EventSetup &iSetup)
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override

◆ produceTrackJets()

void FastjetJetProducer::produceTrackJets ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
protectedvirtual

Definition at line 193 of file FastjetJetProducer.cc.

References PVValHelper::dxy, dxyTrVtxMax_, PVValHelper::dz, dzTrVtxMax_, VirtualJetProducer::fjClusterSeq_, VirtualJetProducer::fjInputs_, VirtualJetProducer::fjJets_, newFWLiteAna::found, VirtualJetProducer::getConstituents(), mps_fire::i, iEvent, input_chrefcand_token_, VirtualJetProducer::input_vertex_token_, VirtualJetProducer::inputs_, VirtualJetProducer::inputTowers(), metsig::jet, PDWG_EXODelayedJetMET_cff::jets, LogDebug, maxVtxZ_, minVtxNdof_, eostools::move(), heavyFlavorDQMFirstStep_cff::pvCollection, multPhiCorr_741_25nsDY_cfi::px, multPhiCorr_741_25nsDY_cfi::py, runAlgorithm(), edm::swap(), useOnlyOnePV_, useOnlyVertexTracks_, trackerHitRTTI::vector, VirtualJetProducer::vertex_, and reco::writeSpecific().

Referenced by produce().

193  {
194  // read in the track candidates
196  iEvent.getByToken(input_chrefcand_token_, inputsHandle);
197 
198  // make collection with pointers so we can play around with it
199  std::vector<edm::Ptr<reco::RecoChargedRefCandidate>> allInputs;
200  std::vector<edm::Ptr<reco::Candidate>> origInputs;
201  for (size_t i = 0; i < inputsHandle->size(); ++i) {
202  allInputs.push_back(inputsHandle->ptrAt(i));
203  origInputs.push_back(inputsHandle->ptrAt(i));
204  }
205 
206  // read in the PV collection
209  // define the overall output jet container
210  auto jets = std::make_unique<std::vector<reco::TrackJet>>();
211 
212  // loop over the good vertices, clustering for each vertex separately
213  for (reco::VertexCollection::const_iterator itVtx = pvCollection->begin(); itVtx != pvCollection->end(); ++itVtx) {
214  if (itVtx->isFake() || itVtx->ndof() < minVtxNdof_ || fabs(itVtx->z()) > maxVtxZ_)
215  continue;
216 
217  // clear the intermediate containers
218  inputs_.clear();
219  fjInputs_.clear();
220  fjJets_.clear();
221 
222  // if only vertex-associated tracks should be used
223  if (useOnlyVertexTracks_) {
224  // loop over the tracks associated to the vertex
225  for (reco::Vertex::trackRef_iterator itTr = itVtx->tracks_begin(); itTr != itVtx->tracks_end(); ++itTr) {
226  // whether a match was found in the track candidate input
227  bool found = false;
228  // loop over input track candidates
229  for (std::vector<edm::Ptr<reco::RecoChargedRefCandidate>>::iterator itIn = allInputs.begin();
230  itIn != allInputs.end();
231  ++itIn) {
232  // match the input track candidate to the track from the vertex
233  reco::TrackRef trref(itTr->castTo<reco::TrackRef>());
234  // check if the tracks match
235  if ((*itIn)->track() == trref) {
236  found = true;
237  // add this track candidate to the input for clustering
238  inputs_.push_back(*itIn);
239  // erase the track candidate from the total list of input, so we don't reuse it later
240  allInputs.erase(itIn);
241  // found the candidate track corresponding to the vertex track, so stop the loop
242  break;
243  } // end if match found
244  } // end loop over input tracks
245  // give an info message in case no match is found (can happen if candidates are subset of tracks used for clustering)
246  if (!found)
247  edm::LogInfo("FastjetTrackJetProducer")
248  << "Ignoring a track at vertex which is not in input track collection!";
249  } // end loop over tracks associated to vertex
250  // if all inpt track candidates should be used
251  } else {
252  // loop over input track candidates
253  for (std::vector<edm::Ptr<reco::RecoChargedRefCandidate>>::iterator itIn = allInputs.begin();
254  itIn != allInputs.end();
255  ++itIn) {
256  // check if the track is close enough to the vertex
257  float dz = (*itIn)->track()->dz(itVtx->position());
258  float dxy = (*itIn)->track()->dxy(itVtx->position());
259  if (fabs(dz) > dzTrVtxMax_)
260  continue;
261  if (fabs(dxy) > dxyTrVtxMax_)
262  continue;
263  bool closervtx = false;
264  // now loop over the good vertices a second time
265  for (reco::VertexCollection::const_iterator itVtx2 = pvCollection->begin(); itVtx2 != pvCollection->end();
266  ++itVtx2) {
267  if (itVtx->isFake() || itVtx->ndof() < minVtxNdof_ || fabs(itVtx->z()) > maxVtxZ_)
268  continue;
269  // and check this track is closer to any other vertex (if more than 1 vertex considered)
270  if (!useOnlyOnePV_ && itVtx != itVtx2 && fabs((*itIn)->track()->dz(itVtx2->position())) < fabs(dz)) {
271  closervtx = true;
272  break; // 1 closer vertex makes the track already not matched, so break
273  }
274  }
275  // don't add this track if another vertex is found closer
276  if (closervtx)
277  continue;
278  // add this track candidate to the input for clustering
279  inputs_.push_back(*itIn);
280  // erase the track candidate from the total list of input, so we don't reuse it later
281  allInputs.erase(itIn);
282  // take a step back in the loop since we just erased
283  --itIn;
284  }
285  }
286 
287  // convert candidates in inputs_ to fastjet::PseudoJets in fjInputs_
288  fjInputs_.reserve(inputs_.size());
289  inputTowers();
290  LogDebug("FastjetTrackJetProducer") << "Inputted towers\n";
291 
292  // run algorithm, using fjInputs_, modifying fjJets_ and allocating fjClusterSeq_
293  runAlgorithm(iEvent, iSetup);
294  LogDebug("FastjetTrackJetProducer") << "Ran algorithm\n";
295 
296  // convert our jets and add to the overall jet vector
297  for (unsigned int ijet = 0; ijet < fjJets_.size(); ++ijet) {
298  // get the constituents from fastjet
299  std::vector<fastjet::PseudoJet> fjConstituents = sorted_by_pt(fjClusterSeq_->constituents(fjJets_[ijet]));
300  // convert them to CandidatePtr vector
301  std::vector<reco::CandidatePtr> constituents = getConstituents(fjConstituents);
302  // fill the trackjet
304  // write the specifics to the jet (simultaneously sets 4-vector, vertex).
306  jet,
307  reco::Particle::LorentzVector(fjJets_[ijet].px(), fjJets_[ijet].py(), fjJets_[ijet].pz(), fjJets_[ijet].E()),
308  vertex_,
309  constituents);
310  jet.setJetArea(0);
311  jet.setPileup(0);
312  jet.setPrimaryVertex(edm::Ref<reco::VertexCollection>(pvCollection, (int)(itVtx - pvCollection->begin())));
313  jet.setVertex(itVtx->position());
314  jets->push_back(jet);
315  }
316 
317  if (useOnlyOnePV_)
318  break; // stop vertex loop if only one vertex asked for
319  } // end loop over vertices
320 
321  // put the jets in the collection
322  LogDebug("FastjetTrackJetProducer") << "Put " << jets->size() << " jets in the event.\n";
323  iEvent.put(std::move(jets));
324 
325  // Clear the work vectors so that memory is free for other modules.
326  // Use the trick of swapping with an empty vector so that the memory
327  // is actually given back rather than silently kept.
328  decltype(fjInputs_)().swap(fjInputs_);
329  decltype(fjJets_)().swap(fjJets_);
330  decltype(inputs_)().swap(inputs_);
331 }
reco::Particle::Point vertex_
virtual std::vector< reco::CandidatePtr > getConstituents(const std::vector< fastjet::PseudoJet > &fjConstituents)
std::vector< fastjet::PseudoJet > fjJets_
virtual void inputTowers()
void writeSpecific(reco::CaloJet &jet, reco::Particle::LorentzVector const &p4, reco::Particle::Point const &point, std::vector< reco::CandidatePtr > const &constituents, CaloGeometry const &geometry, HcalTopology const &topology)
Definition: JetSpecific.cc:32
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:112
std::vector< fastjet::PseudoJet > fjInputs_
int iEvent
Definition: GenABIO.cc:224
std::vector< edm::Ptr< reco::Candidate > > inputs_
ClusterSequencePtr fjClusterSeq_
edm::EDGetTokenT< edm::View< reco::RecoChargedRefCandidate > > input_chrefcand_token_
for CMSBoostedTauSeedingAlgorithm : max depth for descending into clustering sequence ...
Jets made out of tracks.
Definition: TrackJet.h:24
Log< level::Info, false > LogInfo
edm::EDGetTokenT< reco::VertexCollection > input_vertex_token_
void runAlgorithm(edm::Event &iEvent, const edm::EventSetup &iSetup) override
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
def move(src, dest)
Definition: eostools.py:511
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:38
#define LogDebug(id)

◆ runAlgorithm()

void FastjetJetProducer::runAlgorithm ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprotectedvirtual

Implements VirtualJetProducer.

Reimplemented in cms::HTTTopJetProducer.

Definition at line 334 of file FastjetJetProducer.cc.

References beta_, correctShape_, csRho_EtaMax_, csRParam_, VirtualJetProducer::doAreaFastjet_, VirtualJetProducer::doRhoFastjet_, dRMax_, dRMin_, ALCARECOTkAlBeamHalo_cff::filter, VirtualJetProducer::fjAreaDefinition_, VirtualJetProducer::fjClusterSeq_, VirtualJetProducer::fjInputs_, VirtualJetProducer::fjJetDefinition_, VirtualJetProducer::fjJets_, gridMaxRapidity_, gridSpacing_, VirtualJetProducer::jetPtMin_, maxDepth_, muCut_, muMax_, muMin_, nFilt_, TriggerAnalyzer::passed, R0_, RcutFactor_, rFilt_, rFiltDynamic_, ALCARECOPromptCalibProdSiPixelAli0T_cff::Selector, subjetPtMin_, trimPtFracMin_, useCMSBoostedTauSeedingAlgorithm_, useConstituentSubtraction_, useDynamicFiltering_, useFiltering_, useKtPruning_, useMassDropTagger_, usePruning_, useSoftDrop_, useTrimming_, VirtualJetProducer::verbosity_, VirtualJetProducer::voronoiRfact_, yCut_, yMax_, yMin_, and zCut_.

Referenced by produceTrackJets().

334  {
335  // run algorithm
336  /*
337  fjInputs_.clear();
338  double px, py , pz, E;
339  string line;
340  std::ifstream fin("dump3.txt");
341  while (getline(fin, line)){
342  if (line == "#END") break;
343  if (line.substr(0,1) == "#") {continue;}
344  istringstream istr(line);
345  istr >> px >> py >> pz >> E;
346  // create a fastjet::PseudoJet with these components and put it onto
347  // back of the input_particles vector
348  fastjet::PseudoJet j(px,py,pz,E);
349  //if ( fabs(j.rap()) < inputEtaMax )
350  fjInputs_.push_back(fastjet::PseudoJet(px,py,pz,E));
351  }
352  fin.close();
353  */
354 
355  if (!doAreaFastjet_ && !doRhoFastjet_) {
356  fjClusterSeq_ = std::make_shared<fastjet::ClusterSequence>(fjInputs_, *fjJetDefinition_);
357  } else if (voronoiRfact_ <= 0) {
358  fjClusterSeq_ =
359  ClusterSequencePtr(new fastjet::ClusterSequenceArea(fjInputs_, *fjJetDefinition_, *fjAreaDefinition_));
360  } else {
362  new fastjet::ClusterSequenceVoronoiArea(fjInputs_, *fjJetDefinition_, fastjet::VoronoiAreaSpec(voronoiRfact_)));
363  }
364 
367  fjJets_ = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
368  } else {
369  fjJets_.clear();
370 
371  transformer_coll transformers;
372 
373  std::vector<fastjet::PseudoJet> tempJets = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
374 
375  unique_ptr<fastjet::JetMedianBackgroundEstimator> bge_rho;
377  fastjet::Selector rho_range = fastjet::SelectorAbsRapMax(csRho_EtaMax_);
378  bge_rho = std::make_unique<fastjet::JetMedianBackgroundEstimator>(
379  rho_range, fastjet::JetDefinition(fastjet::kt_algorithm, csRParam_), *fjAreaDefinition_);
380  bge_rho->set_particles(fjInputs_);
381  fastjet::contrib::ConstituentSubtractor* constituentSubtractor =
382  new fastjet::contrib::ConstituentSubtractor(bge_rho.get());
383 
384  transformers.push_back(transformer_ptr(constituentSubtractor));
385  };
386  if (useMassDropTagger_) {
387  fastjet::MassDropTagger* md_tagger = new fastjet::MassDropTagger(muCut_, yCut_);
388  transformers.push_back(transformer_ptr(md_tagger));
389  }
391  fastjet::contrib::CMSBoostedTauSeedingAlgorithm* tau_tagger = new fastjet::contrib::CMSBoostedTauSeedingAlgorithm(
393  transformers.push_back(transformer_ptr(tau_tagger));
394  }
395  if (useTrimming_) {
396  fastjet::Filter* trimmer = new fastjet::Filter(fastjet::JetDefinition(fastjet::kt_algorithm, rFilt_),
397  fastjet::SelectorPtFractionMin(trimPtFracMin_));
398  transformers.push_back(transformer_ptr(trimmer));
399  }
400  if ((useFiltering_) && (!useDynamicFiltering_)) {
401  fastjet::Filter* filter = new fastjet::Filter(fastjet::JetDefinition(fastjet::cambridge_algorithm, rFilt_),
402  fastjet::SelectorNHardest(nFilt_));
403  transformers.push_back(transformer_ptr(filter));
404  }
405 
406  if ((usePruning_) && (!useKtPruning_)) {
407  fastjet::Pruner* pruner = new fastjet::Pruner(fastjet::cambridge_algorithm, zCut_, RcutFactor_);
408  transformers.push_back(transformer_ptr(pruner));
409  }
410 
411  if (useDynamicFiltering_) {
413  new fastjet::Filter(fastjet::Filter(&*rFiltDynamic_, fastjet::SelectorNHardest(nFilt_)));
414  transformers.push_back(transformer_ptr(filter));
415  }
416 
417  if (useKtPruning_) {
418  fastjet::Pruner* pruner = new fastjet::Pruner(fastjet::kt_algorithm, zCut_, RcutFactor_);
419  transformers.push_back(transformer_ptr(pruner));
420  }
421 
422  if (useSoftDrop_) {
423  fastjet::contrib::SoftDrop* sd = new fastjet::contrib::SoftDrop(beta_, zCut_, R0_);
424  transformers.push_back(transformer_ptr(sd));
425  }
426 
427  unique_ptr<fastjet::Subtractor> subtractor;
428  unique_ptr<fastjet::GridMedianBackgroundEstimator> bge_rho_grid;
429  if (correctShape_) {
430  bge_rho_grid = std::make_unique<fastjet::GridMedianBackgroundEstimator>(gridMaxRapidity_, gridSpacing_);
431  bge_rho_grid->set_particles(fjInputs_);
432  subtractor = std::make_unique<fastjet::Subtractor>(bge_rho_grid.get());
433  subtractor->set_use_rho_m();
434  //subtractor->use_common_bge_for_rho_and_rhom(true);
435  }
436 
437  for (std::vector<fastjet::PseudoJet>::const_iterator ijet = tempJets.begin(), ijetEnd = tempJets.end();
438  ijet != ijetEnd;
439  ++ijet) {
440  fastjet::PseudoJet transformedJet = *ijet;
441  bool passed = true;
442  for (transformer_coll::const_iterator itransf = transformers.begin(), itransfEnd = transformers.end();
443  itransf != itransfEnd;
444  ++itransf) {
445  if (transformedJet != 0) {
446  transformedJet = (**itransf)(transformedJet);
447  } else {
448  passed = false;
449  }
450  }
451 
452  if (correctShape_) {
453  transformedJet = (*subtractor)(transformedJet);
454  }
455 
456  if (passed) {
457  fjJets_.push_back(transformedJet);
458  }
459  }
460  }
461 }
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
std::vector< fastjet::PseudoJet > fjJets_
bool useDynamicFiltering_
Jet filtering technique.
std::shared_ptr< fastjet::ClusterSequence > ClusterSequencePtr
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
std::vector< fastjet::PseudoJet > fjInputs_
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) ...
bool usePruning_
Jet trimming technique.
bool useCMSBoostedTauSeedingAlgorithm_
Jet pruning technique.
ClusterSequencePtr fjClusterSeq_
dd4hep::Filter Filter
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.
DynamicRfiltPtr rFiltDynamic_
for filtering, trimming: dR scale of sub-clustering
double dRMin_
for CMSBoostedTauSeedingAlgorithm : max asymmetry
double beta_
for constituent subtraction : R parameter for KT alg in jet median background estimator ...
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
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)
AreaDefinitionPtr fjAreaDefinition_

Member Data Documentation

◆ beta_

double FastjetJetProducer::beta_
private

for constituent subtraction : R parameter for KT alg in jet median background estimator

Definition at line 99 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ correctShape_

bool FastjetJetProducer::correctShape_
private

Soft drop.

Definition at line 87 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ csRho_EtaMax_

double FastjetJetProducer::csRho_EtaMax_
private

for pruning: constituent dR * pt/2m < rcut_factor

Definition at line 97 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ csRParam_

double FastjetJetProducer::csRParam_
private

for constituent subtraction : maximum rapidity for ghosts

Definition at line 98 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ dRMax_

double FastjetJetProducer::dRMax_
private

for CMSBoostedTauSeedingAlgorithm : min dR

Definition at line 110 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ dRMin_

double FastjetJetProducer::dRMin_
private

for CMSBoostedTauSeedingAlgorithm : max asymmetry

Definition at line 109 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ dxyTrVtxMax_

float FastjetJetProducer::dxyTrVtxMax_
private

Definition at line 73 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and produceTrackJets().

◆ dzTrVtxMax_

float FastjetJetProducer::dzTrVtxMax_
private

Definition at line 72 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and produceTrackJets().

◆ gridMaxRapidity_

double FastjetJetProducer::gridMaxRapidity_
private

for soft drop : R0 (angular distance normalization - should be set to jet radius in most cases)

Definition at line 101 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ gridSpacing_

double FastjetJetProducer::gridSpacing_
private

for shape subtraction, get the fixed-grid rho

Definition at line 102 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ input_chrefcand_token_

edm::EDGetTokenT<edm::View<reco::RecoChargedRefCandidate> > FastjetJetProducer::input_chrefcand_token_
private

for CMSBoostedTauSeedingAlgorithm : max depth for descending into clustering sequence

Definition at line 114 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and produceTrackJets().

◆ maxDepth_

int FastjetJetProducer::maxDepth_
private

for CMSBoostedTauSeedingAlgorithm : max dR

Definition at line 111 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ maxVtxZ_

float FastjetJetProducer::maxVtxZ_
private

Definition at line 75 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and produceTrackJets().

◆ minVtxNdof_

int FastjetJetProducer::minVtxNdof_
private

Definition at line 74 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and produceTrackJets().

◆ muCut_

double FastjetJetProducer::muCut_
private

Correct the shape of the jets.

Definition at line 88 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ muMax_

double FastjetJetProducer::muMax_
private

for CMSBoostedTauSeedingAlgorithm : min mass-drop

Definition at line 106 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ muMin_

double FastjetJetProducer::muMin_
private

for CMSBoostedTauSeedingAlgorithm : subjet pt min

Definition at line 105 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ nFilt_

int FastjetJetProducer::nFilt_
private

for dynamic filtering radius (as in arXiv:0802.2470)

Definition at line 93 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ R0_

double FastjetJetProducer::R0_
private

for soft drop : beta (angular exponent)

Definition at line 100 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ RcutFactor_

double FastjetJetProducer::RcutFactor_
private

for pruning OR soft drop: constituent minimum pt fraction of parent cluster

Definition at line 96 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ rFilt_

double FastjetJetProducer::rFilt_
private

for mass-drop tagging, symmetry cut: min(pt1^2,pt2^2) * dR(1,2) / mjet > ycut

Definition at line 90 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ rFiltDynamic_

DynamicRfiltPtr FastjetJetProducer::rFiltDynamic_
private

for filtering, trimming: dR scale of sub-clustering

Definition at line 91 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ rFiltFactor_

double FastjetJetProducer::rFiltFactor_
private

for dynamic filtering radius (as in arXiv:0802.2470)

Definition at line 92 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer().

◆ subjetPtMin_

double FastjetJetProducer::subjetPtMin_
private

for shape subtraction, get the grid spacing

Definition at line 104 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ trimPtFracMin_

double FastjetJetProducer::trimPtFracMin_
private

for filtering, pruning: number of subjets expected

Definition at line 94 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ useCMSBoostedTauSeedingAlgorithm_

bool FastjetJetProducer::useCMSBoostedTauSeedingAlgorithm_
private

Jet pruning technique.

Definition at line 83 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ useConstituentSubtraction_

bool FastjetJetProducer::useConstituentSubtraction_
private

Use Kt clustering algorithm for pruning (default is Cambridge/Aachen)

Definition at line 85 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ useDynamicFiltering_

bool FastjetJetProducer::useDynamicFiltering_
private

Jet filtering technique.

Definition at line 80 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ useFiltering_

bool FastjetJetProducer::useFiltering_
private

Mass-drop tagging for boosted Higgs.

Definition at line 79 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ useKtPruning_

bool FastjetJetProducer::useKtPruning_
private

algorithm for seeding reconstruction of boosted Taus (similar to mass-drop tagging)

Definition at line 84 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ useMassDropTagger_

bool FastjetJetProducer::useMassDropTagger_
private

Definition at line 78 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ useOnlyOnePV_

bool FastjetJetProducer::useOnlyOnePV_
private

Definition at line 71 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and produceTrackJets().

◆ useOnlyVertexTracks_

bool FastjetJetProducer::useOnlyVertexTracks_
private

Definition at line 70 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and produceTrackJets().

◆ usePruning_

bool FastjetJetProducer::usePruning_
private

Jet trimming technique.

Definition at line 82 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ useSoftDrop_

bool FastjetJetProducer::useSoftDrop_
private

constituent subtraction technique

Definition at line 86 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ useTrimming_

bool FastjetJetProducer::useTrimming_
private

Use dynamic filtering radius (as in arXiv:0802.2470)

Definition at line 81 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ yCut_

double FastjetJetProducer::yCut_
private

for mass-drop tagging, m0/mjet (m0 = mass of highest mass subjet)

Definition at line 89 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ yMax_

double FastjetJetProducer::yMax_
private

for CMSBoostedTauSeedingAlgorithm : min asymmetry

Definition at line 108 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ yMin_

double FastjetJetProducer::yMin_
private

for CMSBoostedTauSeedingAlgorithm : max mass-drop

Definition at line 107 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().

◆ zCut_

double FastjetJetProducer::zCut_
private

for trimming: constituent minimum pt fraction of full jet

Definition at line 95 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and runAlgorithm().