CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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_HEP_TOPTAG ) {
50  fjHEPTopTagger_ = std::auto_ptr<fastjet::HEPTopTagger>(
51  new fastjet::HEPTopTagger(conf.getParameter<double>("muCut"),
52  conf.getParameter<double>("maxSubjetMass"),
53  conf.getParameter<bool>("useSubjetMass")
54  )
55  );
56  }
57  else if (tagAlgo_ == FJ_JHU_TOPTAG ) {
58  fjJHUTopTagger_ = std::auto_ptr<fastjet::JHTopTagger>(
59  new fastjet::JHTopTagger(conf.getParameter<double>("ptFrac"),
60  conf.getParameter<double>("deltaRCut"),
61  conf.getParameter<double>("cosThetaWMax")
62  )
63  );
64  }
65  else if (tagAlgo_ == FJ_NSUB_TAG ) {
66 
67  fastjet::JetDefinition::Plugin *plugin = new fastjet::SISConePlugin(0.6, 0.75);
68  fastjet::JetDefinition NsubJetDef(plugin);
69  fjNSUBTagger_ = std::auto_ptr<fastjet::RestFrameNSubjettinessTagger>(
70  new fastjet::RestFrameNSubjettinessTagger(NsubJetDef,
71  conf.getParameter<double>("tau2Cut"),
72  conf.getParameter<double>("cosThetaSCut"),
73  conf.getParameter<bool>("useExclusive")
74  )
75  );
76  }
77 
78 
79 
80 
81 
82 }
83 
85 {
87 }
88 
90 {
91  if ( !doAreaFastjet_ && !doRhoFastjet_) {
92  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequence( fjInputs_, *fjJetDefinition_ ) );
93  } else if (voronoiRfact_ <= 0) {
94  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequenceArea( fjInputs_, *fjJetDefinition_ , *fjAreaDefinition_ ) );
95  } else {
96  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequenceVoronoiArea( fjInputs_, *fjJetDefinition_ , fastjet::VoronoiAreaSpec(voronoiRfact_) ) );
97  }
98 
99  if (tagAlgo_ == CA_TOPTAGGER){
100  (*legacyCMSTopTagger_).run( fjInputs_, fjJets_, fjClusterSeq_ );
101 
102  }
103  else {
104 
105  //Run the jet clustering
106  vector<fastjet::PseudoJet> inclusiveJets = fjClusterSeq_->inclusive_jets(ptMin_);
107 
108  if ( verbose_ ) cout << "Getting central jets" << endl;
109  // Find the transient central jets
110  vector<fastjet::PseudoJet> centralJets;
111  for (unsigned int i = 0; i < inclusiveJets.size(); i++) {
112 
113  if (inclusiveJets[i].perp() > ptMin_ && fabs(inclusiveJets[i].rapidity()) < centralEtaCut_) {
114  centralJets.push_back(inclusiveJets[i]);
115  }
116  }
117 
118  fastjet::CMSTopTagger & CMSTagger = *fjCMSTopTagger_;
119  fastjet::HEPTopTagger & HEPTagger = *fjHEPTopTagger_;
120  fastjet::JHTopTagger & JHUTagger = *fjJHUTopTagger_;
121  fastjet::RestFrameNSubjettinessTagger & NSUBTagger = *fjNSUBTagger_;
122 
123 
124  vector<fastjet::PseudoJet>::iterator jetIt = centralJets.begin(), centralJetsEnd = centralJets.end();
125  if ( verbose_ )cout<<"Loop over jets"<<endl;
126  for ( ; jetIt != centralJetsEnd; ++jetIt ) {
127 
128  if (verbose_) cout << "CMS FJ jet pt: " << (*jetIt).perp() << endl;
129 
130  fastjet::PseudoJet taggedJet;
131  if (tagAlgo_ == FJ_CMS_TOPTAG) taggedJet = CMSTagger.result(*jetIt);
132  else if (tagAlgo_ == FJ_HEP_TOPTAG) taggedJet = HEPTagger.result(*jetIt);
133  else if (tagAlgo_ == FJ_JHU_TOPTAG) taggedJet = JHUTagger.result(*jetIt);
134  else if (tagAlgo_ == FJ_NSUB_TAG) taggedJet = NSUBTagger.result(*jetIt);
135  else cout << "NOT A VALID TAGGING ALGORITHM CHOICE!" << endl;
136 
137  if (taggedJet != 0) fjJets_.push_back(taggedJet);
138  }
139  }
140 }
141 //define this as a plug-in
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
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::HEPTopTagger > fjHEPTopTagger_
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)
tuple conf
Definition: dbtoconf.py:185
virtual void runAlgorithm(edm::Event &iEvent, const edm::EventSetup &iSetup)
T perp() const
Magnitude of transverse component.
tuple cout
Definition: gather_cfg.py:121
std::auto_ptr< fastjet::RestFrameNSubjettinessTagger > fjNSUBTagger_
AreaDefinitionPtr fjAreaDefinition_
boost::shared_ptr< fastjet::ClusterSequence > ClusterSequencePtr