CMS 3D CMS Logo

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