Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "FWCore/Framework/interface/Event.h"
0002 #include "FWCore/Framework/interface/stream/EDProducer.h"
0003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0005 
0006 #include "PhysicsTools/JetMCUtils/interface/combination.h"
0007 #include "AnalysisDataFormats/TopObjects/interface/TtSemiLepEvtPartons.h"
0008 #include "DataFormats/PatCandidates/interface/Lepton.h"
0009 #include "AnalysisDataFormats/TopObjects/interface/TtSemiEvtSolution.h"
0010 
0011 #include "TopQuarkAnalysis/TopHitFit/interface/RunHitFit.h"
0012 #include "TopQuarkAnalysis/TopHitFit/interface/Top_Decaykin.h"
0013 #include "TopQuarkAnalysis/TopHitFit/interface/LeptonTranslatorBase.h"
0014 #include "TopQuarkAnalysis/TopHitFit/interface/JetTranslatorBase.h"
0015 #include "TopQuarkAnalysis/TopHitFit/interface/METTranslatorBase.h"
0016 
0017 template <typename LeptonCollection>
0018 class TtSemiLepHitFitProducer : public edm::stream::EDProducer<> {
0019 public:
0020   explicit TtSemiLepHitFitProducer(const edm::ParameterSet&);
0021 
0022 private:
0023   // produce
0024   void produce(edm::Event&, const edm::EventSetup&) override;
0025 
0026   edm::EDGetTokenT<std::vector<pat::Jet> > jetsToken_;
0027   edm::EDGetTokenT<LeptonCollection> lepsToken_;
0028   edm::EDGetTokenT<std::vector<pat::MET> > metsToken_;
0029 
0030   /// maximal number of jets (-1 possible to indicate 'all')
0031   int maxNJets_;
0032   /// maximal number of combinations to be written to the event
0033   int maxNComb_;
0034 
0035   /// maximum eta value for muons, needed to limited range in which resolutions are provided
0036   double maxEtaMu_;
0037   /// maximum eta value for electrons, needed to limited range in which resolutions are provided
0038   double maxEtaEle_;
0039   /// maximum eta value for jets, needed to limited range in which resolutions are provided
0040   double maxEtaJet_;
0041 
0042   /// input tag for b-tagging algorithm
0043   std::string bTagAlgo_;
0044   /// min value of bTag for a b-jet
0045   double minBTagValueBJet_;
0046   /// max value of bTag for a non-b-jet
0047   double maxBTagValueNonBJet_;
0048   /// switch to tell whether to use b-tagging or not
0049   bool useBTag_;
0050 
0051   /// constraints
0052   double mW_;
0053   double mTop_;
0054 
0055   /// jet correction level
0056   std::string jetCorrectionLevel_;
0057 
0058   /// jet energy scale
0059   double jes_;
0060   double jesB_;
0061 
0062   struct FitResult {
0063     int Status;
0064     double Chi2;
0065     double Prob;
0066     double MT;
0067     double SigMT;
0068     pat::Particle HadB;
0069     pat::Particle HadP;
0070     pat::Particle HadQ;
0071     pat::Particle LepB;
0072     pat::Particle LepL;
0073     pat::Particle LepN;
0074     std::vector<int> JetCombi;
0075     bool operator<(const FitResult& rhs) { return Chi2 < rhs.Chi2; };
0076   };
0077 
0078   typedef hitfit::RunHitFit<pat::Electron, pat::Muon, pat::Jet, pat::MET> PatHitFit;
0079 
0080   edm::FileInPath hitfitDefault_;
0081   edm::FileInPath hitfitElectronResolution_;
0082   edm::FileInPath hitfitMuonResolution_;
0083   edm::FileInPath hitfitUdscJetResolution_;
0084   edm::FileInPath hitfitBJetResolution_;
0085   edm::FileInPath hitfitMETResolution_;
0086 
0087   hitfit::LeptonTranslatorBase<pat::Electron> electronTranslator_;
0088   hitfit::LeptonTranslatorBase<pat::Muon> muonTranslator_;
0089   hitfit::JetTranslatorBase<pat::Jet> jetTranslator_;
0090   hitfit::METTranslatorBase<pat::MET> metTranslator_;
0091 
0092   std::unique_ptr<PatHitFit> HitFit;
0093 };
0094 
0095 template <typename LeptonCollection>
0096 TtSemiLepHitFitProducer<LeptonCollection>::TtSemiLepHitFitProducer(const edm::ParameterSet& cfg)
0097     : jetsToken_(consumes<std::vector<pat::Jet> >(cfg.getParameter<edm::InputTag>("jets"))),
0098       lepsToken_(consumes<LeptonCollection>(cfg.getParameter<edm::InputTag>("leps"))),
0099       metsToken_(consumes<std::vector<pat::MET> >(cfg.getParameter<edm::InputTag>("mets"))),
0100       maxNJets_(cfg.getParameter<int>("maxNJets")),
0101       maxNComb_(cfg.getParameter<int>("maxNComb")),
0102       bTagAlgo_(cfg.getParameter<std::string>("bTagAlgo")),
0103       minBTagValueBJet_(cfg.getParameter<double>("minBDiscBJets")),
0104       maxBTagValueNonBJet_(cfg.getParameter<double>("maxBDiscLightJets")),
0105       useBTag_(cfg.getParameter<bool>("useBTagging")),
0106       mW_(cfg.getParameter<double>("mW")),
0107       mTop_(cfg.getParameter<double>("mTop")),
0108       jetCorrectionLevel_(cfg.getParameter<std::string>("jetCorrectionLevel")),
0109       jes_(cfg.getParameter<double>("jes")),
0110       jesB_(cfg.getParameter<double>("jesB")),
0111 
0112       // The following five initializers read the config parameters for the
0113       // ASCII text files which contains the physics object resolutions.
0114       hitfitDefault_(cfg.getUntrackedParameter<edm::FileInPath>(
0115           std::string("hitfitDefault"),
0116           edm::FileInPath(std::string("TopQuarkAnalysis/TopHitFit/data/setting/RunHitFitConfiguration.txt")))),
0117       hitfitElectronResolution_(cfg.getUntrackedParameter<edm::FileInPath>(
0118           std::string("hitfitElectronResolution"),
0119           edm::FileInPath(std::string("TopQuarkAnalysis/TopHitFit/data/resolution/tqafElectronResolution.txt")))),
0120       hitfitMuonResolution_(cfg.getUntrackedParameter<edm::FileInPath>(
0121           std::string("hitfitMuonResolution"),
0122           edm::FileInPath(std::string("TopQuarkAnalysis/TopHitFit/data/resolution/tqafMuonResolution.txt")))),
0123       hitfitUdscJetResolution_(cfg.getUntrackedParameter<edm::FileInPath>(
0124           std::string("hitfitUdscJetResolution"),
0125           edm::FileInPath(std::string("TopQuarkAnalysis/TopHitFit/data/resolution/tqafUdscJetResolution.txt")))),
0126       hitfitBJetResolution_(cfg.getUntrackedParameter<edm::FileInPath>(
0127           std::string("hitfitBJetResolution"),
0128           edm::FileInPath(std::string("TopQuarkAnalysis/TopHitFit/data/resolution/tqafBJetResolution.txt")))),
0129       hitfitMETResolution_(cfg.getUntrackedParameter<edm::FileInPath>(
0130           std::string("hitfitMETResolution"),
0131           edm::FileInPath(std::string("TopQuarkAnalysis/TopHitFit/data/resolution/tqafKtResolution.txt")))),
0132 
0133       // The following four initializers instantiate the translator between PAT objects
0134       // and HitFit objects using the ASCII text files which contains the resolutions.
0135       electronTranslator_(hitfitElectronResolution_.fullPath()),
0136       muonTranslator_(hitfitMuonResolution_.fullPath()),
0137       jetTranslator_(
0138           hitfitUdscJetResolution_.fullPath(), hitfitBJetResolution_.fullPath(), jetCorrectionLevel_, jes_, jesB_),
0139       metTranslator_(hitfitMETResolution_.fullPath())
0140 
0141 {
0142   // Create an instance of RunHitFit and initialize it.
0143   HitFit = std::make_unique<PatHitFit>(
0144       electronTranslator_, muonTranslator_, jetTranslator_, metTranslator_, hitfitDefault_.fullPath(), mW_, mW_, mTop_);
0145 
0146   maxEtaMu_ = 2.4;
0147   maxEtaEle_ = 2.5;
0148   maxEtaJet_ = 2.5;
0149 
0150   edm::LogVerbatim("TopHitFit") << "\n"
0151                                 << "+++++++++++ TtSemiLepHitFitProducer ++++++++++++ \n"
0152                                 << " Due to the eta ranges for which resolutions     \n"
0153                                 << " are provided in                                 \n"
0154                                 << " TopQuarkAnalysis/TopHitFit/data/resolution/     \n"
0155                                 << " so far, the following cuts are currently        \n"
0156                                 << " implemented in the TtSemiLepHitFitProducer:     \n"
0157                                 << " |eta(muons    )| <= " << maxEtaMu_ << " \n"
0158                                 << " |eta(electrons)| <= " << maxEtaEle_ << " \n"
0159                                 << " |eta(jets     )| <= " << maxEtaJet_ << " \n"
0160                                 << "++++++++++++++++++++++++++++++++++++++++++++++++ \n";
0161 
0162   produces<std::vector<pat::Particle> >("PartonsHadP");
0163   produces<std::vector<pat::Particle> >("PartonsHadQ");
0164   produces<std::vector<pat::Particle> >("PartonsHadB");
0165   produces<std::vector<pat::Particle> >("PartonsLepB");
0166   produces<std::vector<pat::Particle> >("Leptons");
0167   produces<std::vector<pat::Particle> >("Neutrinos");
0168 
0169   produces<std::vector<std::vector<int> > >();
0170   produces<std::vector<double> >("Chi2");
0171   produces<std::vector<double> >("Prob");
0172   produces<std::vector<double> >("MT");
0173   produces<std::vector<double> >("SigMT");
0174   produces<std::vector<int> >("Status");
0175   produces<int>("NumberOfConsideredJets");
0176 }
0177 
0178 template <typename LeptonCollection>
0179 void TtSemiLepHitFitProducer<LeptonCollection>::produce(edm::Event& evt, const edm::EventSetup& setup) {
0180   std::unique_ptr<std::vector<pat::Particle> > pPartonsHadP(new std::vector<pat::Particle>);
0181   std::unique_ptr<std::vector<pat::Particle> > pPartonsHadQ(new std::vector<pat::Particle>);
0182   std::unique_ptr<std::vector<pat::Particle> > pPartonsHadB(new std::vector<pat::Particle>);
0183   std::unique_ptr<std::vector<pat::Particle> > pPartonsLepB(new std::vector<pat::Particle>);
0184   std::unique_ptr<std::vector<pat::Particle> > pLeptons(new std::vector<pat::Particle>);
0185   std::unique_ptr<std::vector<pat::Particle> > pNeutrinos(new std::vector<pat::Particle>);
0186 
0187   std::unique_ptr<std::vector<std::vector<int> > > pCombi(new std::vector<std::vector<int> >);
0188   std::unique_ptr<std::vector<double> > pChi2(new std::vector<double>);
0189   std::unique_ptr<std::vector<double> > pProb(new std::vector<double>);
0190   std::unique_ptr<std::vector<double> > pMT(new std::vector<double>);
0191   std::unique_ptr<std::vector<double> > pSigMT(new std::vector<double>);
0192   std::unique_ptr<std::vector<int> > pStatus(new std::vector<int>);
0193   std::unique_ptr<int> pJetsConsidered(new int);
0194 
0195   edm::Handle<std::vector<pat::Jet> > jets;
0196   evt.getByToken(jetsToken_, jets);
0197 
0198   edm::Handle<std::vector<pat::MET> > mets;
0199   evt.getByToken(metsToken_, mets);
0200 
0201   edm::Handle<LeptonCollection> leps;
0202   evt.getByToken(lepsToken_, leps);
0203 
0204   // -----------------------------------------------------
0205   // skip events with no appropriate lepton candidate in
0206   // or empty MET or less jets than partons
0207   // -----------------------------------------------------
0208 
0209   const unsigned int nPartons = 4;
0210 
0211   // Clear the internal state
0212   HitFit->clear();
0213 
0214   // Add lepton into HitFit
0215   bool foundLepton = false;
0216   if (!leps->empty()) {
0217     double maxEtaLep = maxEtaMu_;
0218     if (!dynamic_cast<const reco::Muon*>(&((*leps)[0])))  // assume electron if it is not a muon
0219       maxEtaLep = maxEtaEle_;
0220     for (unsigned iLep = 0; iLep < (*leps).size() && !foundLepton; ++iLep) {
0221       if (std::abs((*leps)[iLep].eta()) <= maxEtaLep) {
0222         HitFit->AddLepton((*leps)[iLep]);
0223         foundLepton = true;
0224       }
0225     }
0226   }
0227 
0228   // Add jets into HitFit
0229   unsigned int nJetsFound = 0;
0230   for (unsigned iJet = 0; iJet < (*jets).size() && (int)nJetsFound != maxNJets_; ++iJet) {
0231     double jet_eta = (*jets)[iJet].eta();
0232     if ((*jets)[iJet].isCaloJet()) {
0233       jet_eta = ((reco::CaloJet*)((*jets)[iJet]).originalObject())->detectorP4().eta();
0234     }
0235     if (std::abs(jet_eta) <= maxEtaJet_) {
0236       HitFit->AddJet((*jets)[iJet]);
0237       nJetsFound++;
0238     }
0239   }
0240   *pJetsConsidered = nJetsFound;
0241 
0242   // Add missing transverse energy into HitFit
0243   if (!mets->empty())
0244     HitFit->SetMet((*mets)[0]);
0245 
0246   if (!foundLepton || mets->empty() || nJetsFound < nPartons) {
0247     // the kinFit getters return empty objects here
0248     pPartonsHadP->push_back(pat::Particle());
0249     pPartonsHadQ->push_back(pat::Particle());
0250     pPartonsHadB->push_back(pat::Particle());
0251     pPartonsLepB->push_back(pat::Particle());
0252     pLeptons->push_back(pat::Particle());
0253     pNeutrinos->push_back(pat::Particle());
0254     // indices referring to the jet combination
0255     std::vector<int> invalidCombi;
0256     for (unsigned int i = 0; i < nPartons; ++i)
0257       invalidCombi.push_back(-1);
0258     pCombi->push_back(invalidCombi);
0259     // chi2
0260     pChi2->push_back(-1.);
0261     // chi2 probability
0262     pProb->push_back(-1.);
0263     // fitted top mass
0264     pMT->push_back(-1.);
0265     pSigMT->push_back(-1.);
0266     // status of the fitter
0267     pStatus->push_back(-1);
0268     // feed out all products
0269     evt.put(std::move(pCombi));
0270     evt.put(std::move(pPartonsHadP), "PartonsHadP");
0271     evt.put(std::move(pPartonsHadQ), "PartonsHadQ");
0272     evt.put(std::move(pPartonsHadB), "PartonsHadB");
0273     evt.put(std::move(pPartonsLepB), "PartonsLepB");
0274     evt.put(std::move(pLeptons), "Leptons");
0275     evt.put(std::move(pNeutrinos), "Neutrinos");
0276     evt.put(std::move(pChi2), "Chi2");
0277     evt.put(std::move(pProb), "Prob");
0278     evt.put(std::move(pMT), "MT");
0279     evt.put(std::move(pSigMT), "SigMT");
0280     evt.put(std::move(pStatus), "Status");
0281     evt.put(std::move(pJetsConsidered), "NumberOfConsideredJets");
0282     return;
0283   }
0284 
0285   std::list<FitResult> FitResultList;
0286 
0287   //
0288   // BEGIN DECLARATION OF VARIABLES FROM KINEMATIC FIT
0289   //
0290 
0291   // In this part are variables from the
0292   // kinematic fit procedure
0293 
0294   // Number of all permutations of the event
0295   size_t nHitFit = 0;
0296 
0297   // Number of jets in the event
0298   size_t nHitFitJet = 0;
0299 
0300   // Results of the fit for all jet permutations of the event
0301   std::vector<hitfit::Fit_Result> hitFitResult;
0302 
0303   //
0304   // R U N   H I T F I T
0305   //
0306   // Run the kinematic fit and get how many permutations are possible
0307   // in the fit
0308 
0309   nHitFit = HitFit->FitAllPermutation();
0310 
0311   //
0312   // BEGIN PART WHICH EXTRACTS INFORMATION FROM HITFIT
0313   //
0314 
0315   // Get the number of jets
0316   nHitFitJet = HitFit->GetUnfittedEvent()[0].njets();
0317 
0318   // Get the fit results for all permutations
0319   hitFitResult = HitFit->GetFitAllPermutation();
0320 
0321   // Loop over all permutations and extract the information
0322   for (size_t fit = 0; fit != nHitFit; ++fit) {
0323     // Get the event after the fit
0324     hitfit::Lepjets_Event fittedEvent = hitFitResult[fit].ev();
0325 
0326     /*
0327         Get jet permutation according to TQAF convention
0328         11 : leptonic b
0329         12 : hadronic b
0330         13 : hadronic W
0331         14 : hadronic W
0332       */
0333     std::vector<int> hitCombi(4);
0334     for (size_t jet = 0; jet != nHitFitJet; ++jet) {
0335       int jet_type = fittedEvent.jet(jet).type();
0336 
0337       switch (jet_type) {
0338         case 11:
0339           hitCombi[TtSemiLepEvtPartons::LepB] = jet;
0340           break;
0341         case 12:
0342           hitCombi[TtSemiLepEvtPartons::HadB] = jet;
0343           break;
0344         case 13:
0345           hitCombi[TtSemiLepEvtPartons::LightQ] = jet;
0346           break;
0347         case 14:
0348           hitCombi[TtSemiLepEvtPartons::LightQBar] = jet;
0349           break;
0350       }
0351     }
0352 
0353     // Store the kinematic quantities in the corresponding containers.
0354 
0355     hitfit::Lepjets_Event_Jet hadP_ = fittedEvent.jet(hitCombi[TtSemiLepEvtPartons::LightQ]);
0356     hitfit::Lepjets_Event_Jet hadQ_ = fittedEvent.jet(hitCombi[TtSemiLepEvtPartons::LightQBar]);
0357     hitfit::Lepjets_Event_Jet hadB_ = fittedEvent.jet(hitCombi[TtSemiLepEvtPartons::HadB]);
0358     hitfit::Lepjets_Event_Jet lepB_ = fittedEvent.jet(hitCombi[TtSemiLepEvtPartons::LepB]);
0359     hitfit::Lepjets_Event_Lep lepL_ = fittedEvent.lep(0);
0360 
0361     /*
0362   /// input tag for b-tagging algorithm
0363   std::string bTagAlgo_;
0364   /// min value of bTag for a b-jet
0365   double minBTagValueBJet_;
0366   /// max value of bTag for a non-b-jet
0367   double maxBTagValueNonBJet_;
0368   /// switch to tell whether to use b-tagging or not
0369   bool useBTag_;
0370       */
0371 
0372     if (hitFitResult[fit].chisq() > 0  // only take into account converged fits
0373         && (!useBTag_ ||
0374             (useBTag_  // use btag information if chosen
0375              && jets->at(hitCombi[TtSemiLepEvtPartons::LightQ]).bDiscriminator(bTagAlgo_) < maxBTagValueNonBJet_ &&
0376              jets->at(hitCombi[TtSemiLepEvtPartons::LightQBar]).bDiscriminator(bTagAlgo_) < maxBTagValueNonBJet_ &&
0377              jets->at(hitCombi[TtSemiLepEvtPartons::HadB]).bDiscriminator(bTagAlgo_) > minBTagValueBJet_ &&
0378              jets->at(hitCombi[TtSemiLepEvtPartons::LepB]).bDiscriminator(bTagAlgo_) > minBTagValueBJet_))) {
0379       FitResult result;
0380       result.Status = 0;
0381       result.Chi2 = hitFitResult[fit].chisq();
0382       result.Prob = exp(-1.0 * (hitFitResult[fit].chisq()) / 2.0);
0383       result.MT = hitFitResult[fit].mt();
0384       result.SigMT = hitFitResult[fit].sigmt();
0385       result.HadB = pat::Particle(reco::LeafCandidate(
0386           0, math::XYZTLorentzVector(hadB_.p().x(), hadB_.p().y(), hadB_.p().z(), hadB_.p().t()), math::XYZPoint()));
0387       result.HadP = pat::Particle(reco::LeafCandidate(
0388           0, math::XYZTLorentzVector(hadP_.p().x(), hadP_.p().y(), hadP_.p().z(), hadP_.p().t()), math::XYZPoint()));
0389       result.HadQ = pat::Particle(reco::LeafCandidate(
0390           0, math::XYZTLorentzVector(hadQ_.p().x(), hadQ_.p().y(), hadQ_.p().z(), hadQ_.p().t()), math::XYZPoint()));
0391       result.LepB = pat::Particle(reco::LeafCandidate(
0392           0, math::XYZTLorentzVector(lepB_.p().x(), lepB_.p().y(), lepB_.p().z(), lepB_.p().t()), math::XYZPoint()));
0393       result.LepL = pat::Particle(reco::LeafCandidate(
0394           0, math::XYZTLorentzVector(lepL_.p().x(), lepL_.p().y(), lepL_.p().z(), lepL_.p().t()), math::XYZPoint()));
0395       result.LepN = pat::Particle(reco::LeafCandidate(
0396           0,
0397           math::XYZTLorentzVector(
0398               fittedEvent.met().x(), fittedEvent.met().y(), fittedEvent.met().z(), fittedEvent.met().t()),
0399           math::XYZPoint()));
0400       result.JetCombi = hitCombi;
0401 
0402       FitResultList.push_back(result);
0403     }
0404   }
0405 
0406   // sort results w.r.t. chi2 values
0407   FitResultList.sort();
0408 
0409   // -----------------------------------------------------
0410   // feed out result
0411   // starting with the JetComb having the smallest chi2
0412   // -----------------------------------------------------
0413 
0414   if (((unsigned)FitResultList.size()) < 1) {  // in case no fit results were stored in the list (all fits aborted)
0415     pPartonsHadP->push_back(pat::Particle());
0416     pPartonsHadQ->push_back(pat::Particle());
0417     pPartonsHadB->push_back(pat::Particle());
0418     pPartonsLepB->push_back(pat::Particle());
0419     pLeptons->push_back(pat::Particle());
0420     pNeutrinos->push_back(pat::Particle());
0421     // indices referring to the jet combination
0422     std::vector<int> invalidCombi;
0423     for (unsigned int i = 0; i < nPartons; ++i)
0424       invalidCombi.push_back(-1);
0425     pCombi->push_back(invalidCombi);
0426     // chi2
0427     pChi2->push_back(-1.);
0428     // chi2 probability
0429     pProb->push_back(-1.);
0430     // fitted top mass
0431     pMT->push_back(-1.);
0432     pSigMT->push_back(-1.);
0433     // status of the fitter
0434     pStatus->push_back(-1);
0435   } else {
0436     unsigned int iComb = 0;
0437     for (typename std::list<FitResult>::const_iterator result = FitResultList.begin(); result != FitResultList.end();
0438          ++result) {
0439       if (maxNComb_ >= 1 && iComb == (unsigned int)maxNComb_)
0440         break;
0441       iComb++;
0442       // partons
0443       pPartonsHadP->push_back(result->HadP);
0444       pPartonsHadQ->push_back(result->HadQ);
0445       pPartonsHadB->push_back(result->HadB);
0446       pPartonsLepB->push_back(result->LepB);
0447       // lepton
0448       pLeptons->push_back(result->LepL);
0449       // neutrino
0450       pNeutrinos->push_back(result->LepN);
0451       // indices referring to the jet combination
0452       pCombi->push_back(result->JetCombi);
0453       // chi2
0454       pChi2->push_back(result->Chi2);
0455       // chi2 probability
0456       pProb->push_back(result->Prob);
0457       // fitted top mass
0458       pMT->push_back(result->MT);
0459       pSigMT->push_back(result->SigMT);
0460       // status of the fitter
0461       pStatus->push_back(result->Status);
0462     }
0463   }
0464   evt.put(std::move(pCombi));
0465   evt.put(std::move(pPartonsHadP), "PartonsHadP");
0466   evt.put(std::move(pPartonsHadQ), "PartonsHadQ");
0467   evt.put(std::move(pPartonsHadB), "PartonsHadB");
0468   evt.put(std::move(pPartonsLepB), "PartonsLepB");
0469   evt.put(std::move(pLeptons), "Leptons");
0470   evt.put(std::move(pNeutrinos), "Neutrinos");
0471   evt.put(std::move(pChi2), "Chi2");
0472   evt.put(std::move(pProb), "Prob");
0473   evt.put(std::move(pMT), "MT");
0474   evt.put(std::move(pSigMT), "SigMT");
0475   evt.put(std::move(pStatus), "Status");
0476   evt.put(std::move(pJetsConsidered), "NumberOfConsideredJets");
0477 }
0478 
0479 #include "DataFormats/PatCandidates/interface/Muon.h"
0480 #include "DataFormats/PatCandidates/interface/Electron.h"
0481 typedef TtSemiLepHitFitProducer<std::vector<pat::Muon> > TtSemiLepHitFitProducerMuon;
0482 typedef TtSemiLepHitFitProducer<std::vector<pat::Electron> > TtSemiLepHitFitProducerElectron;
0483 
0484 #include "FWCore/Framework/interface/MakerMacros.h"
0485 DEFINE_FWK_MODULE(TtSemiLepHitFitProducerMuon);
0486 DEFINE_FWK_MODULE(TtSemiLepHitFitProducerElectron);