CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TtJetPartonMatch.h
Go to the documentation of this file.
1 #ifndef TtJetPartonMatch_h
2 #define TtJetPartonMatch_h
3 
4 #include <memory>
5 #include <vector>
6 
11 
15 
16 // ----------------------------------------------------------------------
17 // template class for:
18 //
19 // * TtFullHadJetPartonMatch
20 // * TtSemiLepJetPartonMatch
21 // * TtFullLepJetPartonMatch
22 //
23 // the class provides plugins for jet-parton matching corresponding
24 // to the JetPartonMatching class; expected input is one of the
25 // classes in:
26 //
27 // AnalysisDataFormats/TopObjects/interface/TtFullHadEvtPartons.h
28 // AnalysisDataFormats/TopObjects/interface/TtSemiLepEvtPartons.h
29 // AnalysisDataFormats/TopObjects/interface/TtFullLepEvtPartons.h
30 //
31 // output is:
32 // * a vector of vectors containing the indices of the jets in the
33 // input collection matched to the partons in the order defined in
34 // the corresponding Tt*EvtPartons class
35 // * a set of vectors with quality parameters of the matching
36 //
37 // the matching can be performed on an arbitrary number of jet
38 // combinations; per default the combination which matches best
39 // according to the quality parameters will be stored; the length
40 // of the vectors will be 1 then
41 // ----------------------------------------------------------------------
42 
43 template <typename C>
45 
46  public:
47 
49  explicit TtJetPartonMatch(const edm::ParameterSet&);
53  virtual void produce(edm::Event&, const edm::EventSetup&);
54 
55  private:
56 
59 
68  int maxNJets_;
71  int maxNComb_;
75  bool useDeltaR_;
81  double maxDist_;
84 };
85 
86 template<typename C>
88  partons_ (cfg.getParameter<std::vector<std::string> >("partonsToIgnore")),
89  genEvt_ (consumes<TtGenEvent>(edm::InputTag("genEvt"))),
90  jets_ (consumes<edm::View<reco::Jet> >(cfg.getParameter<edm::InputTag>("jets"))),
91  maxNJets_ (cfg.getParameter<int> ("maxNJets" )),
92  maxNComb_ (cfg.getParameter<int> ("maxNComb" )),
93  algorithm_ (readAlgorithm(cfg.getParameter<std::string>("algorithm" ))),
94  useDeltaR_ (cfg.getParameter<bool> ("useDeltaR" )),
95  useMaxDist_(cfg.getParameter<bool> ("useMaxDist" )),
96  maxDist_ (cfg.getParameter<double> ("maxDist" )),
97  verbosity_ (cfg.getParameter<int> ("verbosity" ))
98 {
99  // produces a vector of jet/lepton indices in the order of
100  // * TtSemiLepEvtPartons
101  // * TtFullHadEvtPartons
102  // * TtFullLepEvtPartons
103  // and vectors of the corresponding quality parameters
104  produces<std::vector<std::vector<int> > >();
105  produces<std::vector<double> >("SumPt");
106  produces<std::vector<double> >("SumDR");
107  produces<int>("NumberOfConsideredJets");
108 }
109 
110 template<typename C>
112 {
113 }
114 
115 template<typename C>
116 void
118 {
119  // will write
120  // * parton match
121  // * sumPt
122  // * sumDR
123  // to the event
124  std::auto_ptr<std::vector<std::vector<int> > > match(new std::vector<std::vector<int> >);
125  std::auto_ptr<std::vector<double> > sumPt(new std::vector<double>);
126  std::auto_ptr<std::vector<double> > sumDR(new std::vector<double>);
127  std::auto_ptr<int> pJetsConsidered(new int);
128 
129  // get TtGenEvent and jet collection from the event
131  evt.getByToken(genEvt_, genEvt);
132 
134  evt.getByToken(jets_, topJets);
135 
136  // fill vector of partons in the order of
137  // * TtFullLepEvtPartons
138  // * TtSemiLepEvtPartons
139  // * TtFullHadEvtPartons
140  std::vector<const reco::Candidate*> partons = partons_.vec(*genEvt);
141 
142  // prepare vector of jets
143  std::vector<const reco::Candidate*> jets;
144  for(unsigned int ij=0; ij<topJets->size(); ++ij) {
145  // take all jets if maxNJets_ == -1; otherwise use
146  // maxNJets_ if maxNJets_ is big enough or use same
147  // number of jets as partons if maxNJets_ < number
148  // of partons
149  if(maxNJets_!=-1) {
150  if(maxNJets_>=(int)partons.size()) {
151  if((int)ij==maxNJets_) break;
152  }
153  else {
154  if(ij==partons.size()) break;
155  }
156  }
157  jets.push_back( (const reco::Candidate*) &(*topJets)[ij] );
158  }
159  *pJetsConsidered = jets.size();
160 
161  // do the matching with specified parameters
162  JetPartonMatching jetPartonMatch(partons, jets, algorithm_, useMaxDist_, useDeltaR_, maxDist_);
163 
164  // print some info for each event
165  // if corresponding verbosity level set
166  if(verbosity_>0)
167  jetPartonMatch.print();
168 
169  for(unsigned int ic=0; ic<jetPartonMatch.getNumberOfAvailableCombinations(); ++ic) {
170  if((int)ic>=maxNComb_ && maxNComb_>=0) break;
171  std::vector<int> matches = jetPartonMatch.getMatchesForPartons(ic);
172  partons_.expand(matches); // insert dummy indices for partons that were chosen to be ignored
173  match->push_back( matches );
174  sumPt->push_back( jetPartonMatch.getSumDeltaPt(ic) );
175  sumDR->push_back( jetPartonMatch.getSumDeltaR (ic) );
176  }
177  evt.put(match);
178  evt.put(sumPt, "SumPt");
179  evt.put(sumDR, "SumDR");
180  evt.put(pJetsConsidered, "NumberOfConsideredJets");
181 }
182 
183 template<typename C>
186 {
187  if (str == "totalMinDist" ) return JetPartonMatching::totalMinDist;
188  else if(str == "minSumDist" ) return JetPartonMatching::minSumDist;
189  else if(str == "ptOrderedMinDist") return JetPartonMatching::ptOrderedMinDist;
190  else if(str == "unambiguousOnly" ) return JetPartonMatching::unambiguousOnly;
191  else throw cms::Exception("Configuration")
192  << "Chosen algorithm is not supported: " << str << "\n";
193 }
194 
195 #endif
tuple cfg
Definition: looper.py:293
edm::EDGetTokenT< TtGenEvent > genEvt_
TtGenEvent collection input.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
unsigned int getNumberOfAvailableCombinations()
Class derived from the TopGenEvent for ttbar events.
Definition: TtGenEvent.h:18
JetPartonMatching::algorithms readAlgorithm(const std::string &str)
convert string for algorithm into corresponding enumerator type
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
vector< PseudoJet > jets
std::vector< int > getMatchesForPartons(const unsigned int comb=0)
bool useDeltaR_
switch to choose between deltaR/deltaTheta matching
JetPartonMatching::algorithms algorithm_
choice of algorithm
edm::EDGetTokenT< edm::View< reco::Jet > > jets_
jet collection input
int verbosity_
verbosity level
double getSumDeltaR(const unsigned int comb=0)
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
TtJetPartonMatch(const edm::ParameterSet &)
default conructor
genEvt_(cfg.getParameter< edm::InputTag >("genEvent"))
~TtJetPartonMatch()
default destructor
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
double getSumDeltaPt(const unsigned int comb=0)
virtual void produce(edm::Event &, const edm::EventSetup &)
write jet parton match objects into the event