CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups 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 public:
47  explicit TtJetPartonMatch(const edm::ParameterSet&);
49  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
50 
51 private:
54 
63  int maxNJets_;
66  int maxNComb_;
70  bool useDeltaR_;
76  double maxDist_;
79 };
80 
81 template <typename C>
83  : partons_(cfg.getParameter<std::vector<std::string>>("partonsToIgnore")),
84  genEvt_(consumes<TtGenEvent>(edm::InputTag("genEvt"))),
85  jets_(consumes<edm::View<reco::Jet>>(cfg.getParameter<edm::InputTag>("jets"))),
86  maxNJets_(cfg.getParameter<int>("maxNJets")),
87  maxNComb_(cfg.getParameter<int>("maxNComb")),
88  algorithm_(readAlgorithm(cfg.getParameter<std::string>("algorithm"))),
89  useDeltaR_(cfg.getParameter<bool>("useDeltaR")),
90  useMaxDist_(cfg.getParameter<bool>("useMaxDist")),
91  maxDist_(cfg.getParameter<double>("maxDist")),
92  verbosity_(cfg.getParameter<int>("verbosity")) {
93  // produces a vector of jet/lepton indices in the order of
94  // * TtSemiLepEvtPartons
95  // * TtFullHadEvtPartons
96  // * TtFullLepEvtPartons
97  // and vectors of the corresponding quality parameters
98  produces<std::vector<std::vector<int>>>();
99  produces<std::vector<double>>("SumPt");
100  produces<std::vector<double>>("SumDR");
101  produces<int>("NumberOfConsideredJets");
102 }
103 
104 template <typename C>
106  // will write
107  // * parton match
108  // * sumPt
109  // * sumDR
110  // to the event
111  auto match = std::make_unique<std::vector<std::vector<int>>>();
112  auto sumPt = std::make_unique<std::vector<double>>();
113  auto sumDR = std::make_unique<std::vector<double>>();
114  std::unique_ptr<int> pJetsConsidered(new int);
115 
116  // get TtGenEvent and jet collection from the event
117  const TtGenEvent& genEvt = evt.get(genEvt_);
118 
119  const edm::View<reco::Jet>& topJets = evt.get(jets_);
120 
121  // fill vector of partons in the order of
122  // * TtFullLepEvtPartons
123  // * TtSemiLepEvtPartons
124  // * TtFullHadEvtPartons
125  std::vector<const reco::Candidate*> partons = partons_.vec(genEvt);
126 
127  // prepare vector of jets
128  std::vector<const reco::Candidate*> jets;
129  for (unsigned int ij = 0; ij < topJets.size(); ++ij) {
130  // take all jets if maxNJets_ == -1; otherwise use
131  // maxNJets_ if maxNJets_ is big enough or use same
132  // number of jets as partons if maxNJets_ < number
133  // of partons
134  if (maxNJets_ != -1) {
135  if (maxNJets_ >= (int)partons.size()) {
136  if ((int)ij == maxNJets_)
137  break;
138  } else {
139  if (ij == partons.size())
140  break;
141  }
142  }
143  jets.push_back(&topJets[ij]);
144  }
145  *pJetsConsidered = jets.size();
146 
147  // do the matching with specified parameters
148  JetPartonMatching jetPartonMatch(partons, jets, algorithm_, useMaxDist_, useDeltaR_, maxDist_);
149 
150  // print some info for each event
151  // if corresponding verbosity level set
152  if (verbosity_ > 0)
153  jetPartonMatch.print();
154 
155  for (unsigned int ic = 0; ic < jetPartonMatch.getNumberOfAvailableCombinations(); ++ic) {
156  if ((int)ic >= maxNComb_ && maxNComb_ >= 0)
157  break;
158  std::vector<int> matches = jetPartonMatch.getMatchesForPartons(ic);
159  partons_.expand(matches); // insert dummy indices for partons that were chosen to be ignored
160  match->push_back(matches);
161  sumPt->push_back(jetPartonMatch.getSumDeltaPt(ic));
162  sumDR->push_back(jetPartonMatch.getSumDeltaR(ic));
163  }
164  evt.put(std::move(match));
165  evt.put(std::move(sumPt), "SumPt");
166  evt.put(std::move(sumDR), "SumDR");
167  evt.put(std::move(pJetsConsidered), "NumberOfConsideredJets");
168 }
169 
170 template <typename C>
172  if (str == "totalMinDist")
174  else if (str == "minSumDist")
176  else if (str == "ptOrderedMinDist")
178  else if (str == "unambiguousOnly")
180  else
181  throw cms::Exception("Configuration") << "Chosen algorithm is not supported: " << str << "\n";
182 }
183 
184 #endif
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
tuple cfg
Definition: looper.py:296
edm::EDGetTokenT< TtGenEvent > genEvt_
TtGenEvent collection input.
size_type size() const
unsigned int getNumberOfAvailableCombinations()
Class derived from the TopGenEvent for ttbar events.
Definition: TtGenEvent.h:18
vector< PseudoJet > jets
std::vector< int > getMatchesForPartons(const unsigned int comb=0)
def move
Definition: eostools.py:511
JetPartonMatching::algorithms readAlgorithm(const std::string &str) const
convert string for algorithm into corresponding enumerator type
bool useDeltaR_
switch to choose between deltaR/deltaTheta matching
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:346
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
write jet parton match objects into the event
JetPartonMatching::algorithms algorithm_
choice of algorithm
constexpr char Jet[]
Definition: modules.cc:9
genEvt_(cfg.getParameter< edm::InputTag >("genEvent"))
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
#define str(s)
double getSumDeltaPt(const unsigned int comb=0)