Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:31:14

0001 #include <vector>
0002 
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/Framework/interface/global/EDProducer.h"
0005 #include "FWCore/Framework/interface/Frameworkfwd.h"
0006 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0007 #include "FWCore/Utilities/interface/transform.h"
0008 #include "AnalysisDataFormats/TopObjects/interface/TtEvent.h"
0009 #include "AnalysisDataFormats/TopObjects/interface/TtGenEvent.h"
0010 #include "AnalysisDataFormats/TopObjects/interface/TtFullHadronicEvent.h"
0011 #include "AnalysisDataFormats/TopObjects/interface/TtSemiLeptonicEvent.h"
0012 #include "AnalysisDataFormats/TopObjects/interface/TtFullLeptonicEvent.h"
0013 
0014 /**
0015    \class   TtEvtBuilder TtEvtBuilder.h "TopQuarkAnalysis/TopEventProducers/interface/TtEvtBuilder.h"
0016 
0017    \brief    Template class to fill the TtEvent structure
0018 
0019    Template class to fill the TtEvent structure for:
0020 
0021    * TtSemiLeptonicEvent
0022    * TtFullLeptonicEvent
0023    * TtFullHadronicEvent
0024 
0025    event hypothesis, genEvent and extra information (if
0026    available) are read from the event and contracted into
0027    the TtEvent
0028 */
0029 
0030 template <typename C>
0031 class TtEvtBuilder : public edm::global::EDProducer<> {
0032 public:
0033   /// default constructor
0034   explicit TtEvtBuilder(const edm::ParameterSet&);
0035 
0036 private:
0037   /// produce function (this one is not even accessible for
0038   /// derived classes)
0039   void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0040   /// fill data members that are decay-channel specific
0041   void fillSpecific(C&, const edm::Event&) const;
0042 
0043 private:
0044   /// vebosity level
0045   int verbosity_;
0046   /// vector of hypothesis class names
0047   std::vector<edm::EDGetTokenT<int> > hypKeyTokens_;
0048   std::vector<edm::EDGetTokenT<std::vector<TtEvent::HypoCombPair> > > hypTokens_;
0049   std::vector<edm::EDGetTokenT<int> > hypNeutrTokens_;
0050   std::vector<edm::EDGetTokenT<int> > hypJetTokens_;
0051   typedef std::vector<edm::EDGetTokenT<int> >::const_iterator EventHypoIntToken;
0052   typedef std::vector<edm::EDGetTokenT<std::vector<TtEvent::HypoCombPair> > >::const_iterator EventHypoToken;
0053   /// TtGenEvent
0054   edm::InputTag genEvt_;
0055   edm::EDGetTokenT<TtGenEvent> genEvtToken_;
0056   /// decay channels of the two top decay branches; to be
0057   /// filled according to WDecay::LeptonTypes in TtGenEvent
0058   int decayChnTop1_;
0059   int decayChnTop2_;
0060 
0061   /// input parameters for the kKinFit
0062   /// hypothesis class extras
0063   edm::ParameterSet kinFit_;
0064   edm::EDGetTokenT<std::vector<double> > fitChi2Token_;
0065   edm::EDGetTokenT<std::vector<double> > fitProbToken_;
0066   /// input parameters for the kHitFit
0067   /// hypothesis class extras
0068   edm::ParameterSet hitFit_;
0069   edm::EDGetTokenT<std::vector<double> > hitFitChi2Token_;
0070   edm::EDGetTokenT<std::vector<double> > hitFitProbToken_;
0071   edm::EDGetTokenT<std::vector<double> > hitFitMTToken_;
0072   edm::EDGetTokenT<std::vector<double> > hitFitSigMTToken_;
0073   /// input parameters for the kKinSolution
0074   /// hypothesis class extras
0075   edm::ParameterSet kinSolution_;
0076   edm::EDGetTokenT<std::vector<double> > solWeightToken_;
0077   edm::EDGetTokenT<bool> wrongChargeToken_;
0078   /// input parameters for the kGenMatch
0079   /// hypothesis class extras
0080   edm::ParameterSet genMatch_;
0081   edm::EDGetTokenT<std::vector<double> > sumPtToken_;
0082   edm::EDGetTokenT<std::vector<double> > sumDRToken_;
0083   /// input parameters for the kMVADisc
0084   /// hypothesis class extras
0085   edm::ParameterSet mvaDisc_;
0086   edm::EDGetTokenT<std::string> methToken_;
0087   edm::EDGetTokenT<std::vector<double> > discToken_;
0088 
0089   edm::EDPutTokenT<C> putToken_;
0090 };
0091 
0092 template <typename C>
0093 TtEvtBuilder<C>::TtEvtBuilder(const edm::ParameterSet& cfg)
0094     : verbosity_(cfg.getParameter<int>("verbosity")),
0095       hypKeyTokens_(edm::vector_transform(
0096           cfg.getParameter<std::vector<edm::InputTag> >("hypotheses"),
0097           [this](edm::InputTag const& tag) { return consumes<int>(edm::InputTag(tag.label(), "Key")); })),
0098       hypTokens_(edm::vector_transform(
0099           cfg.getParameter<std::vector<edm::InputTag> >("hypotheses"),
0100           [this](edm::InputTag const& tag) { return consumes<std::vector<TtEvent::HypoCombPair> >(tag); })),
0101       hypNeutrTokens_(edm::vector_transform(cfg.getParameter<std::vector<edm::InputTag> >("hypotheses"),
0102                                             [this](edm::InputTag const& tag) {
0103                                               return consumes<int>(
0104                                                   edm::InputTag(tag.label(), "NumberOfRealNeutrinoSolutions"));
0105                                             })),
0106       hypJetTokens_(edm::vector_transform(cfg.getParameter<std::vector<edm::InputTag> >("hypotheses"),
0107                                           [this](edm::InputTag const& tag) {
0108                                             return consumes<int>(edm::InputTag(tag.label(), "NumberOfConsideredJets"));
0109                                           })),
0110       genEvt_(cfg.getParameter<edm::InputTag>("genEvent")),
0111       genEvtToken_(mayConsume<TtGenEvent>(genEvt_)),
0112       decayChnTop1_(cfg.getParameter<int>("decayChannel1")),
0113       decayChnTop2_(cfg.getParameter<int>("decayChannel2")) {
0114   // parameter subsets for kKinFit
0115   if (cfg.exists("kinFit")) {
0116     kinFit_ = cfg.getParameter<edm::ParameterSet>("kinFit");
0117     fitChi2Token_ = mayConsume<std::vector<double> >(kinFit_.getParameter<edm::InputTag>("chi2"));
0118     fitProbToken_ = mayConsume<std::vector<double> >(kinFit_.getParameter<edm::InputTag>("prob"));
0119   }
0120   // parameter subsets for kHitFit
0121   if (cfg.exists("hitFit")) {
0122     hitFit_ = cfg.getParameter<edm::ParameterSet>("hitFit");
0123     hitFitChi2Token_ = mayConsume<std::vector<double> >(hitFit_.getParameter<edm::InputTag>("chi2"));
0124     hitFitProbToken_ = mayConsume<std::vector<double> >(hitFit_.getParameter<edm::InputTag>("prob"));
0125     hitFitMTToken_ = mayConsume<std::vector<double> >(hitFit_.getParameter<edm::InputTag>("mt"));
0126     hitFitSigMTToken_ = mayConsume<std::vector<double> >(hitFit_.getParameter<edm::InputTag>("sigmt"));
0127   }
0128   // parameter subsets for kKinSolution
0129   if (cfg.exists("kinSolution")) {
0130     kinSolution_ = cfg.getParameter<edm::ParameterSet>("kinSolution");
0131     solWeightToken_ = mayConsume<std::vector<double> >(kinSolution_.getParameter<edm::InputTag>("solWeight"));
0132     wrongChargeToken_ = mayConsume<bool>(kinSolution_.getParameter<edm::InputTag>("wrongCharge"));
0133   }
0134   // parameter subsets for kGenMatch
0135   if (cfg.exists("genMatch")) {
0136     genMatch_ = cfg.getParameter<edm::ParameterSet>("genMatch");
0137     sumPtToken_ = mayConsume<std::vector<double> >(genMatch_.getParameter<edm::InputTag>("sumPt"));
0138     sumDRToken_ = mayConsume<std::vector<double> >(genMatch_.getParameter<edm::InputTag>("sumDR"));
0139   }
0140   // parameter subsets for kMvaDisc
0141   if (cfg.exists("mvaDisc")) {
0142     mvaDisc_ = cfg.getParameter<edm::ParameterSet>("mvaDisc");
0143     methToken_ = mayConsume<std::string>(mvaDisc_.getParameter<edm::InputTag>("meth"));
0144     discToken_ = mayConsume<std::vector<double> >(mvaDisc_.getParameter<edm::InputTag>("disc"));
0145   }
0146   // produces a TtEventEvent for:
0147   //  * TtSemiLeptonicEvent
0148   //  * TtFullLeptonicEvent
0149   //  * TtFullHadronicEvent
0150   // from hypotheses and associated extra information
0151   putToken_ = produces<C>();
0152 }
0153 
0154 template <typename C>
0155 void TtEvtBuilder<C>::produce(edm::StreamID, edm::Event& evt, const edm::EventSetup& setup) const {
0156   C ttEvent;
0157 
0158   // set leptonic decay channels
0159   ttEvent.setLepDecays(WDecay::LeptonType(decayChnTop1_), WDecay::LeptonType(decayChnTop2_));
0160 
0161   // set genEvent (if available)
0162   edm::Handle<TtGenEvent> genEvt;
0163   if (!genEvt_.label().empty())
0164     if (evt.getByToken(genEvtToken_, genEvt))
0165       ttEvent.setGenEvent(genEvt);
0166 
0167   // add event hypotheses for all given
0168   // hypothesis classes to the TtEvent
0169   EventHypoIntToken hKey = hypKeyTokens_.begin();
0170   EventHypoToken h = hypTokens_.begin();
0171   for (; hKey != hypKeyTokens_.end(); ++hKey, ++h) {
0172     const int key = evt.get(*hKey);
0173 
0174     const std::vector<TtEvent::HypoCombPair>& hypMatchVec = evt.get(*h);
0175 
0176     for (const auto& hm : hypMatchVec) {
0177       ttEvent.addEventHypo(static_cast<TtEvent::HypoClassKey>(key), hm);
0178     }
0179   }
0180 
0181   // set kKinFit extras
0182   if (ttEvent.isHypoAvailable(TtEvent::kKinFit)) {
0183     ttEvent.setFitChi2(evt.get(fitChi2Token_));
0184     ttEvent.setFitProb(evt.get(fitProbToken_));
0185   }
0186 
0187   // set kHitFit extras
0188   if (ttEvent.isHypoAvailable(TtEvent::kHitFit)) {
0189     ttEvent.setHitFitChi2(evt.get(hitFitChi2Token_));
0190     ttEvent.setHitFitProb(evt.get(hitFitProbToken_));
0191     ttEvent.setHitFitMT(evt.get(hitFitMTToken_));
0192     ttEvent.setHitFitSigMT(evt.get(hitFitSigMTToken_));
0193   }
0194 
0195   // set kGenMatch extras
0196   if (ttEvent.isHypoAvailable(TtEvent::kGenMatch)) {
0197     ttEvent.setGenMatchSumPt(evt.get(sumPtToken_));
0198     ttEvent.setGenMatchSumDR(evt.get(sumDRToken_));
0199   }
0200 
0201   // set kMvaDisc extras
0202   if (ttEvent.isHypoAvailable(TtEvent::kMVADisc)) {
0203     ttEvent.setMvaMethod(evt.get(methToken_));
0204     ttEvent.setMvaDiscriminators(evt.get(discToken_));
0205   }
0206 
0207   // fill data members that are decay-channel specific
0208   fillSpecific(ttEvent, evt);
0209 
0210   // print summary via MessageLogger if verbosity_>0
0211   ttEvent.print(verbosity_);
0212 
0213   // write object into the edm::Event
0214   evt.emplace(putToken_, std::move(ttEvent));
0215 }
0216 
0217 template <>
0218 inline void TtEvtBuilder<TtFullHadronicEvent>::fillSpecific(TtFullHadronicEvent& ttEvent, const edm::Event& evt) const {
0219 }
0220 
0221 template <>
0222 inline void TtEvtBuilder<TtFullLeptonicEvent>::fillSpecific(TtFullLeptonicEvent& ttEvent, const edm::Event& evt) const {
0223   // set kKinSolution extras
0224   if (ttEvent.isHypoAvailable(TtEvent::kKinSolution)) {
0225     ttEvent.setSolWeight(evt.get(solWeightToken_));
0226     ttEvent.setWrongCharge(evt.get(wrongChargeToken_));
0227   }
0228 }
0229 
0230 template <>
0231 inline void TtEvtBuilder<TtSemiLeptonicEvent>::fillSpecific(TtSemiLeptonicEvent& ttEvent, const edm::Event& evt) const {
0232   EventHypoIntToken hKey = hypKeyTokens_.begin();
0233   EventHypoIntToken hNeutr = hypNeutrTokens_.begin();
0234   EventHypoIntToken hJet = hypJetTokens_.begin();
0235   for (; hKey != hypKeyTokens_.end(); ++hKey, ++hNeutr, ++hJet) {
0236     const int key = evt.get(*hKey);
0237 
0238     // set number of real neutrino solutions for all hypotheses
0239     const int numberOfRealNeutrinoSolutions = evt.get(*hNeutr);
0240     ttEvent.setNumberOfRealNeutrinoSolutions(static_cast<TtEvent::HypoClassKey>(key), numberOfRealNeutrinoSolutions);
0241 
0242     // set number of considered jets for all hypotheses
0243     const int numberOfConsideredJets = evt.get(*hJet);
0244     ttEvent.setNumberOfConsideredJets(static_cast<TtEvent::HypoClassKey>(key), numberOfConsideredJets);
0245   }
0246 }
0247 
0248 #include "AnalysisDataFormats/TopObjects/interface/TtFullHadronicEvent.h"
0249 #include "AnalysisDataFormats/TopObjects/interface/TtFullLeptonicEvent.h"
0250 #include "AnalysisDataFormats/TopObjects/interface/TtSemiLeptonicEvent.h"
0251 using TtFullHadEvtBuilder = TtEvtBuilder<TtFullHadronicEvent>;
0252 using TtFullLepEvtBuilder = TtEvtBuilder<TtFullLeptonicEvent>;
0253 using TtSemiLepEvtBuilder = TtEvtBuilder<TtSemiLeptonicEvent>;
0254 
0255 #include "FWCore/Framework/interface/MakerMacros.h"
0256 DEFINE_FWK_MODULE(TtFullHadEvtBuilder);
0257 DEFINE_FWK_MODULE(TtFullLepEvtBuilder);
0258 DEFINE_FWK_MODULE(TtSemiLepEvtBuilder);