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 
66  int maxNJets_;
69  int maxNComb_;
73  bool useDeltaR_;
79  double maxDist_;
82 };
83 
84 template<typename C>
86  partons_ (cfg.getParameter<std::vector<std::string> >("partonsToIgnore")),
87  jets_ (cfg.getParameter<edm::InputTag> ("jets" )),
88  maxNJets_ (cfg.getParameter<int> ("maxNJets" )),
89  maxNComb_ (cfg.getParameter<int> ("maxNComb" )),
90  algorithm_ (readAlgorithm(cfg.getParameter<std::string>("algorithm" ))),
91  useDeltaR_ (cfg.getParameter<bool> ("useDeltaR" )),
92  useMaxDist_(cfg.getParameter<bool> ("useMaxDist" )),
93  maxDist_ (cfg.getParameter<double> ("maxDist" )),
94  verbosity_ (cfg.getParameter<int> ("verbosity" ))
95 {
96  // produces a vector of jet/lepton indices in the order of
97  // * TtSemiLepEvtPartons
98  // * TtFullHadEvtPartons
99  // * TtFullLepEvtPartons
100  // and vectors of the corresponding quality parameters
101  produces<std::vector<std::vector<int> > >();
102  produces<std::vector<double> >("SumPt");
103  produces<std::vector<double> >("SumDR");
104  produces<int>("NumberOfConsideredJets");
105 }
106 
107 template<typename C>
109 {
110 }
111 
112 template<typename C>
113 void
115 {
116  // will write
117  // * parton match
118  // * sumPt
119  // * sumDR
120  // to the event
121  std::auto_ptr<std::vector<std::vector<int> > > match(new std::vector<std::vector<int> >);
122  std::auto_ptr<std::vector<double> > sumPt(new std::vector<double>);
123  std::auto_ptr<std::vector<double> > sumDR(new std::vector<double>);
124  std::auto_ptr<int> pJetsConsidered(new int);
125 
126  // get TtGenEvent and jet collection from the event
128  evt.getByLabel("genEvt", genEvt);
129 
131  evt.getByLabel(jets_, topJets);
132 
133  // fill vector of partons in the order of
134  // * TtFullLepEvtPartons
135  // * TtSemiLepEvtPartons
136  // * TtFullHadEvtPartons
137  std::vector<const reco::Candidate*> partons = partons_.vec(*genEvt);
138 
139  // prepare vector of jets
140  std::vector<const reco::Candidate*> jets;
141  for(unsigned int ij=0; ij<topJets->size(); ++ij) {
142  // take all jets if maxNJets_ == -1; otherwise use
143  // maxNJets_ if maxNJets_ is big enough or use same
144  // number of jets as partons if maxNJets_ < number
145  // of partons
146  if(maxNJets_!=-1) {
147  if(maxNJets_>=(int)partons.size()) {
148  if((int)ij==maxNJets_) break;
149  }
150  else {
151  if(ij==partons.size()) break;
152  }
153  }
154  jets.push_back( (const reco::Candidate*) &(*topJets)[ij] );
155  }
156  *pJetsConsidered = jets.size();
157 
158  // do the matching with specified parameters
159  JetPartonMatching jetPartonMatch(partons, jets, algorithm_, useMaxDist_, useDeltaR_, maxDist_);
160 
161  // print some info for each event
162  // if corresponding verbosity level set
163  if(verbosity_>0)
164  jetPartonMatch.print();
165 
166  for(unsigned int ic=0; ic<jetPartonMatch.getNumberOfAvailableCombinations(); ++ic) {
167  if((int)ic>=maxNComb_ && maxNComb_>=0) break;
168  std::vector<int> matches = jetPartonMatch.getMatchesForPartons(ic);
169  partons_.expand(matches); // insert dummy indices for partons that were chosen to be ignored
170  match->push_back( matches );
171  sumPt->push_back( jetPartonMatch.getSumDeltaPt(ic) );
172  sumDR->push_back( jetPartonMatch.getSumDeltaR (ic) );
173  }
174  evt.put(match);
175  evt.put(sumPt, "SumPt");
176  evt.put(sumDR, "SumDR");
177  evt.put(pJetsConsidered, "NumberOfConsideredJets");
178 }
179 
180 template<typename C>
183 {
184  if (str == "totalMinDist" ) return JetPartonMatching::totalMinDist;
185  else if(str == "minSumDist" ) return JetPartonMatching::minSumDist;
186  else if(str == "ptOrderedMinDist") return JetPartonMatching::ptOrderedMinDist;
187  else if(str == "unambiguousOnly" ) return JetPartonMatching::unambiguousOnly;
188  else throw cms::Exception("Configuration")
189  << "Chosen algorithm is not supported: " << str << "\n";
190 }
191 
192 #endif
edm::InputTag jets_
jet collection input
unsigned int getNumberOfAvailableCombinations()
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:94
vector< PseudoJet > jets
std::vector< int > getMatchesForPartons(const unsigned int comb=0)
bool useDeltaR_
switch to choose between deltaR/deltaTheta matching
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
JetPartonMatching::algorithms algorithm_
choice of algorithm
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:6
TtJetPartonMatch(const edm::ParameterSet &)
default conructor
~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