CMS 3D CMS Logo

CATopJetProducer.cc
Go to the documentation of this file.
2 #include "CATopJetProducer.h"
3 
4 
5 using namespace edm;
6 using namespace cms;
7 using namespace reco;
8 using namespace std;
9 
10 CATopJetProducer::CATopJetProducer(edm::ParameterSet const& conf):
11  FastjetJetProducer( conf ),
12  tagAlgo_(conf.getParameter<int>("tagAlgo")),
13  ptMin_(conf.getParameter<double>("jetPtMin")),
14  centralEtaCut_(conf.getParameter<double>("centralEtaCut")),
15  verbose_(conf.getParameter<bool>("verbose"))
16 {
17 
18  if (tagAlgo_ == CA_TOPTAGGER ) {
19 
20  legacyCMSTopTagger_ = std::auto_ptr<CATopJetAlgorithm>(
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> >("sumEtBins"), // sumEt bins over which cuts may vary. vector={bin 0 lower bound, bin 1 lower bound, ...}
32  conf.getParameter<std::vector<double> >("rBins"), // Jet distance paramter R. R values depend on sumEt bins.
33  conf.getParameter<std::vector<double> >("ptFracBins"), // fraction of hard jet pt that subjet must have (deltap)
34  conf.getParameter<std::vector<double> >("deltarBins"), // Applicable only if useAdjacency=1. deltar adjacency values for each sumEtBin
35  conf.getParameter<std::vector<double> >("nCellBins"), // Applicable only if useAdjacency=3. number of cells to consider two subjets adjacent
36  conf.getParameter<double>("inputEtMin"), // seed threshold - NOT USED
37  conf.getParameter<bool> ("useMaxTower"), // use max tower as adjacency criterion, otherwise use centroid - NOT USED
38  conf.getParameter<double>("sumEtEtaCut"), // eta for event SumEt - NOT USED
39  conf.getParameter<double>("etFrac") // fraction of event sumEt / 2 for a jet to be considered "hard" -NOT USED
40  ));
41  }
42  else if (tagAlgo_ == FJ_CMS_TOPTAG ) {
43  fjCMSTopTagger_ = std::auto_ptr<fastjet::CMSTopTagger>(
44  new fastjet::CMSTopTagger(conf.getParameter<double> ("ptFrac"),
45  conf.getParameter<double> ("rFrac"),
46  conf.getParameter<double> ("adjacencyParam"))
47  );
48  }
49  else if (tagAlgo_ == FJ_JHU_TOPTAG ) {
50  fjJHUTopTagger_ = std::auto_ptr<fastjet::JHTopTagger>(
51  new fastjet::JHTopTagger(conf.getParameter<double>("ptFrac"),
52  conf.getParameter<double>("deltaRCut"),
53  conf.getParameter<double>("cosThetaWMax")
54  )
55  );
56  }
57  else if (tagAlgo_ == FJ_NSUB_TAG ) {
58 
59  fastjet::JetDefinition::Plugin *plugin = new fastjet::SISConePlugin(0.6, 0.75);
60  fastjet::JetDefinition NsubJetDef(plugin);
61  fjNSUBTagger_ = std::auto_ptr<fastjet::RestFrameNSubjettinessTagger>(
62  new fastjet::RestFrameNSubjettinessTagger(NsubJetDef,
63  conf.getParameter<double>("tau2Cut"),
64  conf.getParameter<double>("cosThetaSCut"),
65  conf.getParameter<bool>("useExclusive")
66  )
67  );
68  }
69 
70 
71 
72 
73 
74 }
75 
77 {
79 }
80 
82 {
83  if ( !doAreaFastjet_ && !doRhoFastjet_) {
84  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequence( fjInputs_, *fjJetDefinition_ ) );
85  } else if (voronoiRfact_ <= 0) {
86  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequenceArea( fjInputs_, *fjJetDefinition_ , *fjAreaDefinition_ ) );
87  } else {
88  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequenceVoronoiArea( fjInputs_, *fjJetDefinition_ , fastjet::VoronoiAreaSpec(voronoiRfact_) ) );
89  }
90 
91  if (tagAlgo_ == CA_TOPTAGGER){
92  (*legacyCMSTopTagger_).run( fjInputs_, fjJets_, fjClusterSeq_ );
93 
94  }
95  else {
96 
97  //Run the jet clustering
98  vector<fastjet::PseudoJet> inclusiveJets = fjClusterSeq_->inclusive_jets(ptMin_);
99 
100  if ( verbose_ ) cout << "Getting central jets" << endl;
101  // Find the transient central jets
102  vector<fastjet::PseudoJet> centralJets;
103  for (unsigned int i = 0; i < inclusiveJets.size(); i++) {
104 
105  if (inclusiveJets[i].perp() > ptMin_ && fabs(inclusiveJets[i].rapidity()) < centralEtaCut_) {
106  centralJets.push_back(inclusiveJets[i]);
107  }
108  }
109 
110  fastjet::CMSTopTagger & CMSTagger = *fjCMSTopTagger_;
111  fastjet::JHTopTagger & JHUTagger = *fjJHUTopTagger_;
112  fastjet::RestFrameNSubjettinessTagger & NSUBTagger = *fjNSUBTagger_;
113 
114 
115  vector<fastjet::PseudoJet>::iterator jetIt = centralJets.begin(), centralJetsEnd = centralJets.end();
116  if ( verbose_ )cout<<"Loop over jets"<<endl;
117  for ( ; jetIt != centralJetsEnd; ++jetIt ) {
118 
119  if (verbose_) cout << "CMS FJ jet pt: " << (*jetIt).perp() << endl;
120 
121  fastjet::PseudoJet taggedJet;
122  if (tagAlgo_ == FJ_CMS_TOPTAG) taggedJet = CMSTagger.result(*jetIt);
123  else if (tagAlgo_ == FJ_JHU_TOPTAG) taggedJet = JHUTagger.result(*jetIt);
124  else if (tagAlgo_ == FJ_NSUB_TAG) taggedJet = NSUBTagger.result(*jetIt);
125  else cout << "NOT A VALID TAGGING ALGORITHM CHOICE!" << endl;
126 
127  if (taggedJet != 0) fjJets_.push_back(taggedJet);
128  }
129  }
130 }
131 //define this as a plug-in
T getParameter(std::string const &) const
std::auto_ptr< CATopJetAlgorithm > legacyCMSTopTagger_
auto_ptr< JetDefinition::Plugin > plugin
std::vector< fastjet::PseudoJet > fjJets_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::auto_ptr< fastjet::CMSTopTagger > fjCMSTopTagger_
The algorithm to do the work.
std::auto_ptr< fastjet::JHTopTagger > fjJHUTopTagger_
std::vector< fastjet::PseudoJet > fjInputs_
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup)
int iEvent
Definition: GenABIO.cc:230
ClusterSequencePtr fjClusterSeq_
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup)
virtual void runAlgorithm(edm::Event &iEvent, const edm::EventSetup &iSetup)
T perp() const
Magnitude of transverse component.
fixed size matrix
HLT enums.
std::auto_ptr< fastjet::RestFrameNSubjettinessTagger > fjNSUBTagger_
AreaDefinitionPtr fjAreaDefinition_
boost::shared_ptr< fastjet::ClusterSequence > ClusterSequencePtr