CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TtEvtBuilder.cc
Go to the documentation of this file.
1 #include <vector>
2 
13 
30 template <typename C>
32 public:
34  explicit TtEvtBuilder(const edm::ParameterSet&);
35 
36 private:
39  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
41  void fillSpecific(C&, const edm::Event&) const;
42 
43 private:
47  std::vector<edm::EDGetTokenT<int> > hypKeyTokens_;
48  std::vector<edm::EDGetTokenT<std::vector<TtEvent::HypoCombPair> > > hypTokens_;
49  std::vector<edm::EDGetTokenT<int> > hypNeutrTokens_;
50  std::vector<edm::EDGetTokenT<int> > hypJetTokens_;
51  typedef std::vector<edm::EDGetTokenT<int> >::const_iterator EventHypoIntToken;
52  typedef std::vector<edm::EDGetTokenT<std::vector<TtEvent::HypoCombPair> > >::const_iterator EventHypoToken;
60 
88 
90 };
91 
92 template <typename C>
94  : verbosity_(cfg.getParameter<int>("verbosity")),
95  hypKeyTokens_(edm::vector_transform(
96  cfg.getParameter<std::vector<edm::InputTag> >("hypotheses"),
97  [this](edm::InputTag const& tag) { return consumes<int>(edm::InputTag(tag.label(), "Key")); })),
99  cfg.getParameter<std::vector<edm::InputTag> >("hypotheses"),
100  [this](edm::InputTag const& tag) { return consumes<std::vector<TtEvent::HypoCombPair> >(tag); })),
101  hypNeutrTokens_(edm::vector_transform(cfg.getParameter<std::vector<edm::InputTag> >("hypotheses"),
102  [this](edm::InputTag const& tag) {
103  return consumes<int>(
104  edm::InputTag(tag.label(), "NumberOfRealNeutrinoSolutions"));
105  })),
106  hypJetTokens_(edm::vector_transform(cfg.getParameter<std::vector<edm::InputTag> >("hypotheses"),
107  [this](edm::InputTag const& tag) {
108  return consumes<int>(edm::InputTag(tag.label(), "NumberOfConsideredJets"));
109  })),
110  genEvt_(cfg.getParameter<edm::InputTag>("genEvent")),
111  genEvtToken_(mayConsume<TtGenEvent>(genEvt_)),
112  decayChnTop1_(cfg.getParameter<int>("decayChannel1")),
113  decayChnTop2_(cfg.getParameter<int>("decayChannel2")) {
114  // parameter subsets for kKinFit
115  if (cfg.exists("kinFit")) {
116  kinFit_ = cfg.getParameter<edm::ParameterSet>("kinFit");
117  fitChi2Token_ = mayConsume<std::vector<double> >(kinFit_.getParameter<edm::InputTag>("chi2"));
118  fitProbToken_ = mayConsume<std::vector<double> >(kinFit_.getParameter<edm::InputTag>("prob"));
119  }
120  // parameter subsets for kHitFit
121  if (cfg.exists("hitFit")) {
122  hitFit_ = cfg.getParameter<edm::ParameterSet>("hitFit");
123  hitFitChi2Token_ = mayConsume<std::vector<double> >(hitFit_.getParameter<edm::InputTag>("chi2"));
124  hitFitProbToken_ = mayConsume<std::vector<double> >(hitFit_.getParameter<edm::InputTag>("prob"));
125  hitFitMTToken_ = mayConsume<std::vector<double> >(hitFit_.getParameter<edm::InputTag>("mt"));
126  hitFitSigMTToken_ = mayConsume<std::vector<double> >(hitFit_.getParameter<edm::InputTag>("sigmt"));
127  }
128  // parameter subsets for kKinSolution
129  if (cfg.exists("kinSolution")) {
130  kinSolution_ = cfg.getParameter<edm::ParameterSet>("kinSolution");
131  solWeightToken_ = mayConsume<std::vector<double> >(kinSolution_.getParameter<edm::InputTag>("solWeight"));
132  wrongChargeToken_ = mayConsume<bool>(kinSolution_.getParameter<edm::InputTag>("wrongCharge"));
133  }
134  // parameter subsets for kGenMatch
135  if (cfg.exists("genMatch")) {
136  genMatch_ = cfg.getParameter<edm::ParameterSet>("genMatch");
137  sumPtToken_ = mayConsume<std::vector<double> >(genMatch_.getParameter<edm::InputTag>("sumPt"));
138  sumDRToken_ = mayConsume<std::vector<double> >(genMatch_.getParameter<edm::InputTag>("sumDR"));
139  }
140  // parameter subsets for kMvaDisc
141  if (cfg.exists("mvaDisc")) {
142  mvaDisc_ = cfg.getParameter<edm::ParameterSet>("mvaDisc");
143  methToken_ = mayConsume<std::string>(mvaDisc_.getParameter<edm::InputTag>("meth"));
144  discToken_ = mayConsume<std::vector<double> >(mvaDisc_.getParameter<edm::InputTag>("disc"));
145  }
146  // produces a TtEventEvent for:
147  // * TtSemiLeptonicEvent
148  // * TtFullLeptonicEvent
149  // * TtFullHadronicEvent
150  // from hypotheses and associated extra information
151  putToken_ = produces<C>();
152 }
153 
154 template <typename C>
156  C ttEvent;
157 
158  // set leptonic decay channels
160 
161  // set genEvent (if available)
163  if (!genEvt_.label().empty())
164  if (evt.getByToken(genEvtToken_, genEvt))
165  ttEvent.setGenEvent(genEvt);
166 
167  // add event hypotheses for all given
168  // hypothesis classes to the TtEvent
169  EventHypoIntToken hKey = hypKeyTokens_.begin();
170  EventHypoToken h = hypTokens_.begin();
171  for (; hKey != hypKeyTokens_.end(); ++hKey, ++h) {
172  const int key = evt.get(*hKey);
173 
174  const std::vector<TtEvent::HypoCombPair>& hypMatchVec = evt.get(*h);
175 
176  for (const auto& hm : hypMatchVec) {
177  ttEvent.addEventHypo(static_cast<TtEvent::HypoClassKey>(key), hm);
178  }
179  }
180 
181  // set kKinFit extras
182  if (ttEvent.isHypoAvailable(TtEvent::kKinFit)) {
183  ttEvent.setFitChi2(evt.get(fitChi2Token_));
184  ttEvent.setFitProb(evt.get(fitProbToken_));
185  }
186 
187  // set kHitFit extras
188  if (ttEvent.isHypoAvailable(TtEvent::kHitFit)) {
189  ttEvent.setHitFitChi2(evt.get(hitFitChi2Token_));
190  ttEvent.setHitFitProb(evt.get(hitFitProbToken_));
191  ttEvent.setHitFitMT(evt.get(hitFitMTToken_));
192  ttEvent.setHitFitSigMT(evt.get(hitFitSigMTToken_));
193  }
194 
195  // set kGenMatch extras
196  if (ttEvent.isHypoAvailable(TtEvent::kGenMatch)) {
197  ttEvent.setGenMatchSumPt(evt.get(sumPtToken_));
198  ttEvent.setGenMatchSumDR(evt.get(sumDRToken_));
199  }
200 
201  // set kMvaDisc extras
202  if (ttEvent.isHypoAvailable(TtEvent::kMVADisc)) {
203  ttEvent.setMvaMethod(evt.get(methToken_));
204  ttEvent.setMvaDiscriminators(evt.get(discToken_));
205  }
206 
207  // fill data members that are decay-channel specific
208  fillSpecific(ttEvent, evt);
209 
210  // print summary via MessageLogger if verbosity_>0
211  ttEvent.print(verbosity_);
212 
213  // write object into the edm::Event
214  evt.emplace(putToken_, std::move(ttEvent));
215 }
216 
217 template <>
219 }
220 
221 template <>
223  // set kKinSolution extras
224  if (ttEvent.isHypoAvailable(TtEvent::kKinSolution)) {
225  ttEvent.setSolWeight(evt.get(solWeightToken_));
226  ttEvent.setWrongCharge(evt.get(wrongChargeToken_));
227  }
228 }
229 
230 template <>
232  EventHypoIntToken hKey = hypKeyTokens_.begin();
233  EventHypoIntToken hNeutr = hypNeutrTokens_.begin();
234  EventHypoIntToken hJet = hypJetTokens_.begin();
235  for (; hKey != hypKeyTokens_.end(); ++hKey, ++hNeutr, ++hJet) {
236  const int key = evt.get(*hKey);
237 
238  // set number of real neutrino solutions for all hypotheses
239  const int numberOfRealNeutrinoSolutions = evt.get(*hNeutr);
240  ttEvent.setNumberOfRealNeutrinoSolutions(static_cast<TtEvent::HypoClassKey>(key), numberOfRealNeutrinoSolutions);
241 
242  // set number of considered jets for all hypotheses
243  const int numberOfConsideredJets = evt.get(*hJet);
244  ttEvent.setNumberOfConsideredJets(static_cast<TtEvent::HypoClassKey>(key), numberOfConsideredJets);
245  }
246 }
247 
254 
edm::EDGetTokenT< std::vector< double > > solWeightToken_
Definition: TtEvtBuilder.cc:76
void setWrongCharge(const bool &val)
set right or wrong charge combination of kKinSolution hypothesis
bool isHypoAvailable(const std::string &key, const unsigned &cmb=0) const
Definition: TtEvent.h:70
tuple cfg
Definition: looper.py:296
decayChnTop2_(cfg.getParameter< int >("decayChannel2"))
edm::EDGetTokenT< std::vector< double > > hitFitMTToken_
Definition: TtEvtBuilder.cc:71
edm::EDGetTokenT< std::string > methToken_
Definition: TtEvtBuilder.cc:86
TtEvtBuilder(const edm::ParameterSet &)
default constructor
Definition: TtEvtBuilder.cc:93
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
hypJetTokens_(edm::vector_transform(cfg.getParameter< std::vector< edm::InputTag > >("hypotheses"), [this](edm::InputTag const &tag){return consumes< int >(edm::InputTag(tag.label(),"NumberOfConsideredJets"));}))
edm::InputTag genEvt_
TtGenEvent.
Definition: TtEvtBuilder.cc:54
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
genEvtToken_(mayConsume< TtGenEvent >(genEvt_))
void setNumberOfRealNeutrinoSolutions(const HypoClassKey &key, const int &nr)
set number of real neutrino solutions for a given hypo class
Class derived from the TtEvent for the semileptonic decay channel.
std::vector< edm::EDGetTokenT< std::vector< TtEvent::HypoCombPair > > > hypTokens_
Definition: TtEvtBuilder.cc:48
hypTokens_(edm::vector_transform(cfg.getParameter< std::vector< edm::InputTag > >("hypotheses"), [this](edm::InputTag const &tag){return consumes< std::vector< TtEvent::HypoCombPair > >(tag);}))
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
std::vector< edm::EDGetTokenT< int > >::const_iterator EventHypoIntToken
Definition: TtEvtBuilder.cc:51
edm::ParameterSet genMatch_
Definition: TtEvtBuilder.cc:80
decayChnTop1_(cfg.getParameter< int >("decayChannel1"))
hypNeutrTokens_(edm::vector_transform(cfg.getParameter< std::vector< edm::InputTag > >("hypotheses"), [this](edm::InputTag const &tag){return consumes< int >(edm::InputTag(tag.label(),"NumberOfRealNeutrinoSolutions"));}))
std::vector< edm::EDGetTokenT< int > > hypNeutrTokens_
Definition: TtEvtBuilder.cc:49
edm::EDGetTokenT< std::vector< double > > sumDRToken_
Definition: TtEvtBuilder.cc:82
Template class to fill the TtEvent structure.
Definition: TtEvtBuilder.cc:31
std::vector< edm::EDGetTokenT< int > > hypKeyTokens_
vector of hypothesis class names
Definition: TtEvtBuilder.cc:47
edm::EDGetTokenT< std::vector< double > > discToken_
Definition: TtEvtBuilder.cc:87
def move
Definition: eostools.py:511
tuple key
prepare the HTCondor submission files and eventually submit them
std::vector< edm::EDGetTokenT< int > > hypJetTokens_
Definition: TtEvtBuilder.cc:50
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:346
int verbosity_
vebosity level
Definition: TtEvtBuilder.cc:45
edm::EDGetTokenT< std::vector< double > > hitFitSigMTToken_
Definition: TtEvtBuilder.cc:72
edm::ParameterSet mvaDisc_
Definition: TtEvtBuilder.cc:85
edm::EDGetTokenT< bool > wrongChargeToken_
Definition: TtEvtBuilder.cc:77
Class derived from the TtEvent for the full leptonic decay channel.
edm::EDGetTokenT< std::vector< double > > fitChi2Token_
Definition: TtEvtBuilder.cc:64
void setSolWeight(const std::vector< double > &val)
set weight of kKinSolution hypothesis
edm::ParameterSet kinFit_
Definition: TtEvtBuilder.cc:63
edm::ParameterSet kinSolution_
Definition: TtEvtBuilder.cc:75
Class derived from the TtEvent for the full hadronic decay channel.
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
OrphanHandle< PROD > emplace(EDPutTokenT< PROD > token, Args &&...args)
puts a new product
Definition: Event.h:433
edm::EDGetTokenT< std::vector< double > > hitFitChi2Token_
Definition: TtEvtBuilder.cc:69
std::vector< edm::EDGetTokenT< std::vector< TtEvent::HypoCombPair > > >::const_iterator EventHypoToken
Definition: TtEvtBuilder.cc:52
edm::ParameterSet hitFit_
Definition: TtEvtBuilder.cc:68
edm::EDGetTokenT< std::vector< double > > fitProbToken_
Definition: TtEvtBuilder.cc:65
void setNumberOfConsideredJets(const HypoClassKey &key, const unsigned int nJets)
set number of jets considered when building a given hypothesis
Definition: TtEvent.h:162
std::string const & label() const
Definition: InputTag.h:36
void fillSpecific(C &, const edm::Event &) const
fill data members that are decay-channel specific
genEvt_(cfg.getParameter< edm::InputTag >("genEvent"))
edm::EDPutTokenT< C > putToken_
Definition: TtEvtBuilder.cc:89
edm::EDGetTokenT< std::vector< double > > sumPtToken_
Definition: TtEvtBuilder.cc:81
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
edm::EDGetTokenT< TtGenEvent > genEvtToken_
Definition: TtEvtBuilder.cc:55
edm::EDGetTokenT< std::vector< double > > hitFitProbToken_
Definition: TtEvtBuilder.cc:70