CMS 3D CMS Logo

CATopJetProducer.cc
Go to the documentation of this file.
3 
4 #include <memory>
5 
7 
8 using namespace edm;
9 using namespace cms;
10 using namespace reco;
11 using namespace std;
12 
13 CATopJetProducer::CATopJetProducer(edm::ParameterSet const& conf)
14  : FastjetJetProducer(conf),
15  tagAlgo_(conf.getParameter<int>("tagAlgo")),
16  ptMin_(conf.getParameter<double>("jetPtMin")),
17  centralEtaCut_(conf.getParameter<double>("centralEtaCut")),
18  verbose_(conf.getParameter<bool>("verbose")) {
19  if (tagAlgo_ == CA_TOPTAGGER) {
20  legacyCMSTopTagger_ = std::unique_ptr<CATopJetAlgorithm>(new CATopJetAlgorithm(
21  src_,
22  conf.getParameter<bool>("verbose"),
23  conf.getParameter<int>("algorithm"), // 0 = KT, 1 = CA, 2 = anti-KT
24  conf.getParameter<int>("useAdjacency"), // choose adjacency requirement:
25  // 0 = no adjacency
26  // 1 = deltar adjacency
27  // 2 = modified adjacency
28  // 3 = calotower neirest neigbor based adjacency (untested)
29  conf.getParameter<double>("centralEtaCut"), // eta for defining "central" jets
30  conf.getParameter<double>("jetPtMin"), // min jet pt
31  conf.getParameter<std::vector<double>>(
32  "sumEtBins"), // sumEt bins over which cuts may vary. vector={bin 0 lower bound, bin 1 lower bound, ...}
33  conf.getParameter<std::vector<double>>("rBins"), // Jet distance paramter R. R values depend on sumEt bins.
34  conf.getParameter<std::vector<double>>("ptFracBins"), // fraction of hard jet pt that subjet must have (deltap)
35  conf.getParameter<std::vector<double>>(
36  "deltarBins"), // Applicable only if useAdjacency=1. deltar adjacency values for each sumEtBin
37  conf.getParameter<std::vector<double>>(
38  "nCellBins"), // Applicable only if useAdjacency=3. number of cells to consider two subjets adjacent
39  conf.getParameter<double>("inputEtMin"), // seed threshold - NOT USED
40  conf.getParameter<bool>(
41  "useMaxTower"), // use max tower as adjacency criterion, otherwise use centroid - NOT USED
42  conf.getParameter<double>("sumEtEtaCut"), // eta for event SumEt - NOT USED
43  conf.getParameter<double>("etFrac") // fraction of event sumEt / 2 for a jet to be considered "hard" -NOT USED
44  ));
45  } else if (tagAlgo_ == FJ_CMS_TOPTAG) {
47  std::unique_ptr<fastjet::CMSTopTagger>(new fastjet::CMSTopTagger(conf.getParameter<double>("ptFrac"),
48  conf.getParameter<double>("rFrac"),
49  conf.getParameter<double>("adjacencyParam")));
50  } else if (tagAlgo_ == FJ_JHU_TOPTAG) {
52  std::unique_ptr<fastjet::JHTopTagger>(new fastjet::JHTopTagger(conf.getParameter<double>("ptFrac"),
53  conf.getParameter<double>("deltaRCut"),
54  conf.getParameter<double>("cosThetaWMax")));
55  } else if (tagAlgo_ == FJ_NSUB_TAG) {
56  fastjet::JetDefinition::Plugin* plugin = new fastjet::SISConePlugin(0.6, 0.75);
57  fastjet::JetDefinition NsubJetDef(plugin);
58  fjNSUBTagger_ = std::unique_ptr<fastjet::RestFrameNSubjettinessTagger>(
59  new fastjet::RestFrameNSubjettinessTagger(NsubJetDef,
60  conf.getParameter<double>("tau2Cut"),
61  conf.getParameter<double>("cosThetaSCut"),
62  conf.getParameter<bool>("useExclusive")));
63  }
64 }
65 
67 
69  if (!doAreaFastjet_ && !doRhoFastjet_) {
70  fjClusterSeq_ = std::make_shared<fastjet::ClusterSequence>(fjInputs_, *fjJetDefinition_);
71  } else if (voronoiRfact_ <= 0) {
73  ClusterSequencePtr(new fastjet::ClusterSequenceArea(fjInputs_, *fjJetDefinition_, *fjAreaDefinition_));
74  } else {
76  new fastjet::ClusterSequenceVoronoiArea(fjInputs_, *fjJetDefinition_, fastjet::VoronoiAreaSpec(voronoiRfact_)));
77  }
78 
79  if (tagAlgo_ == CA_TOPTAGGER) {
80  (*legacyCMSTopTagger_).run(fjInputs_, fjJets_, fjClusterSeq_);
81 
82  } else {
83  //Run the jet clustering
84  vector<fastjet::PseudoJet> inclusiveJets = fjClusterSeq_->inclusive_jets(ptMin_);
85 
86  if (verbose_)
87  cout << "Getting central jets" << endl;
88  // Find the transient central jets
89  vector<fastjet::PseudoJet> centralJets;
90  for (unsigned int i = 0; i < inclusiveJets.size(); i++) {
91  if (inclusiveJets[i].perp() > ptMin_ && fabs(inclusiveJets[i].rapidity()) < centralEtaCut_) {
92  centralJets.push_back(inclusiveJets[i]);
93  }
94  }
95 
96  fastjet::CMSTopTagger& CMSTagger = *fjCMSTopTagger_;
97  fastjet::JHTopTagger& JHUTagger = *fjJHUTopTagger_;
98  fastjet::RestFrameNSubjettinessTagger& NSUBTagger = *fjNSUBTagger_;
99 
100  vector<fastjet::PseudoJet>::iterator jetIt = centralJets.begin(), centralJetsEnd = centralJets.end();
101  if (verbose_)
102  cout << "Loop over jets" << endl;
103  for (; jetIt != centralJetsEnd; ++jetIt) {
104  if (verbose_)
105  cout << "CMS FJ jet pt: " << (*jetIt).perp() << endl;
106 
107  fastjet::PseudoJet taggedJet;
108  if (tagAlgo_ == FJ_CMS_TOPTAG)
109  taggedJet = CMSTagger.result(*jetIt);
110  else if (tagAlgo_ == FJ_JHU_TOPTAG)
111  taggedJet = JHUTagger.result(*jetIt);
112  else if (tagAlgo_ == FJ_NSUB_TAG)
113  taggedJet = NSUBTagger.result(*jetIt);
114  else
115  cout << "NOT A VALID TAGGING ALGORITHM CHOICE!" << endl;
116 
117  if (taggedJet != 0)
118  fjJets_.push_back(taggedJet);
119  }
120  }
121 }
122 
126  desc.add<int>("tagAlgo", 0); // choice of top tagging algorithm
127  desc.add<double>("centralEtaCut", 2.5); // eta for defining "central" jets
128  desc.add<bool>("verbose", false);
129  desc.add<string>("jetCollInstanceName", "caTopSubJets"); // subjet collection
130  desc.add<int>("algorithm", 1); // 0 = KT, 1 = CA, 2 = anti-KT
131  desc.add<int>("useAdjacency", 2); // veto adjacent subjets:
132  // 0, no adjacency
133  // 1, deltar adjacency
134  // 2, modified adjacency
135  // 3, calotower neirest neigbor based adjacency (untested)
136  vector<double> sumEtBinsDefault = {0., 1600., 2600.};
137  desc.add<vector<double>>(
138  "sumEtBins",
139  sumEtBinsDefault); // sumEt bins over which cuts vary. vector={bin 0 lower bound, bin 1 lower bound, ...}
140  vector<double> rBinsDefault(3, 0.8);
141  desc.add<vector<double>>("rBins", rBinsDefault); // Jet distance paramter R. R values depend on sumEt bins.
142  vector<double> ptFracBinsDefault(3, 0.05);
143  desc.add<vector<double>>("ptFracBins", ptFracBinsDefault); // minimum fraction of central jet pt for subjets (deltap)
144  vector<double> deltarBinsDefault(3, 0.019);
145  desc.add<vector<double>>(
146  "deltarBins", deltarBinsDefault); // Applicable only if useAdjacency=1. deltar adjacency values for each sumEtBin
147  vector<double> nCellBinsDefault(3, 1.9);
148  desc.add<vector<double>>(
149  "nCellBins",
150  nCellBinsDefault); // Applicable only if useAdjacency=3. number of cells apart for two subjets to be considered "independent"
151  desc.add<bool>("useMaxTower", false); // use max tower in adjacency criterion, otherwise use centroid - NOT USED
152  desc.add<double>("sumEtEtaCut", 3.0); // eta for event SumEt - NOT USED
153  desc.add<double>("etFrac", 0.7); // fraction of event sumEt / 2 for a jet to be considered "hard" - NOT USED
154  desc.add<double>("ptFrac", 0.05);
155  desc.add<double>("rFrac", 0.);
156  desc.add<double>("adjacencyParam", 0.);
157  desc.add<double>("deltaRCut", 0.19);
158  desc.add<double>("cosThetaWMax", 0.7);
159  desc.add<double>("tau2Cut", 0.);
160  desc.add<double>("cosThetaSCut", 0.);
161  desc.add<bool>("useExclusive", false);
167  descriptions.add("CATopJetProducer", desc);
168 }
169 //define this as a plug-in
cms::CATopJetProducer::tagAlgo_
int tagAlgo_
Definition: CATopJetProducer.h:88
cms::CATopJetProducer::FJ_CMS_TOPTAG
Definition: CATopJetProducer.h:92
cms::CATopJetProducer::produce
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
Definition: CATopJetProducer.cc:66
electrons_cff.bool
bool
Definition: electrons_cff.py:372
mps_fire.i
i
Definition: mps_fire.py:355
edm::ParameterSetDescription::add
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:95
cms::CATopJetProducer::ptMin_
double ptMin_
Definition: CATopJetProducer.h:89
VirtualJetProducer::ClusterSequencePtr
std::shared_ptr< fastjet::ClusterSequence > ClusterSequencePtr
Definition: VirtualJetProducer.h:73
edm
HLT enums.
Definition: AlignableModifier.h:19
gather_cfg.cout
cout
Definition: gather_cfg.py:144
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
cms::CATopJetProducer::fjNSUBTagger_
std::unique_ptr< fastjet::RestFrameNSubjettinessTagger > fjNSUBTagger_
Definition: CATopJetProducer.h:86
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
perp
T perp() const
Magnitude of transverse component.
Definition: Basic3DVectorLD.h:133
VirtualJetProducer::voronoiRfact_
double voronoiRfact_
Definition: VirtualJetProducer.h:170
VirtualJetProducer::fjClusterSeq_
ClusterSequencePtr fjClusterSeq_
Definition: VirtualJetProducer.h:185
MakerMacros.h
cms::CATopJetProducer::legacyCMSTopTagger_
std::unique_ptr< CATopJetAlgorithm > legacyCMSTopTagger_
Definition: CATopJetProducer.h:83
cms::CATopJetProducer::FJ_JHU_TOPTAG
Definition: CATopJetProducer.h:92
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
cms::CATopJetProducer::fjJHUTopTagger_
std::unique_ptr< fastjet::JHTopTagger > fjJHUTopTagger_
Definition: CATopJetProducer.h:85
cms::CATopJetProducer::runAlgorithm
void runAlgorithm(edm::Event &iEvent, const edm::EventSetup &iSetup) override
Definition: CATopJetProducer.cc:68
cms::CATopJetProducer::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: CATopJetProducer.cc:123
FastjetJetProducer::fillDescriptionsFromFastJetProducer
static void fillDescriptionsFromFastJetProducer(edm::ParameterSetDescription &desc)
Definition: FastjetJetProducer.cc:481
VirtualJetProducer::doRhoFastjet_
bool doRhoFastjet_
Definition: VirtualJetProducer.h:167
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
cms::CATopJetProducer::verbose_
bool verbose_
Definition: CATopJetProducer.h:91
VirtualJetProducer::src_
edm::InputTag src_
Definition: VirtualJetProducer.h:148
edm::ParameterSet
Definition: ParameterSet.h:36
createfilelist.int
int
Definition: createfilelist.py:10
FastjetJetProducer::produce
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
Definition: FastjetJetProducer.cc:176
iEvent
int iEvent
Definition: GenABIO.cc:224
fwrapper::plugin
unique_ptr< JetDefinition::Plugin > plugin
Definition: fastjetfortran_madfks.cc:43
edm::EventSetup
Definition: EventSetup.h:57
VirtualJetProducer::fjAreaDefinition_
AreaDefinitionPtr fjAreaDefinition_
Definition: VirtualJetProducer.h:189
HltBtagPostValidation_cff.c
c
Definition: HltBtagPostValidation_cff.py:31
cms::CATopJetProducer::FJ_NSUB_TAG
Definition: CATopJetProducer.h:92
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
FastjetJetProducer.h
std
Definition: JetResolutionObject.h:76
cms::CATopJetProducer::centralEtaCut_
double centralEtaCut_
Definition: CATopJetProducer.h:90
cms::CATopJetProducer::CA_TOPTAGGER
Definition: CATopJetProducer.h:92
VirtualJetProducer::doAreaFastjet_
bool doAreaFastjet_
Definition: VirtualJetProducer.h:163
CATopJetAlgorithm
Definition: CATopJetAlgorithm.h:41
VirtualJetProducer::fjJets_
std::vector< fastjet::PseudoJet > fjJets_
Definition: VirtualJetProducer.h:192
VirtualJetProducer::fjJetDefinition_
JetDefPtr fjJetDefinition_
Definition: VirtualJetProducer.h:186
CATopJetProducer.h
VirtualJetProducer::fjInputs_
std::vector< fastjet::PseudoJet > fjInputs_
Definition: VirtualJetProducer.h:191
edm::Event
Definition: Event.h:73
FastjetJetProducer
Definition: FastjetJetProducer.h:41
cms::CATopJetProducer::fjCMSTopTagger_
std::unique_ptr< fastjet::CMSTopTagger > fjCMSTopTagger_
The algorithm to do the work.
Definition: CATopJetProducer.h:84
VirtualJetProducer::fillDescriptionsFromVirtualJetProducer
static void fillDescriptionsFromVirtualJetProducer(edm::ParameterSetDescription &desc)
Definition: VirtualJetProducer.cc:1016
cms
Namespace of DDCMS conversion namespace.
Definition: ProducerAnalyzer.cc:21
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37
cms::CATopJetProducer
Definition: CATopJetProducer.h:71