Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:24:12

0001 // -*- C++ -*-
0002 //
0003 // Package:    TagProbeFitTreeProducer
0004 // Class:      TagProbeFitTreeProducer
0005 //
0006 /**\class TagProbeFitTreeProducer TagProbeFitTreeProducer.cc
0007 
0008  Description: <one line class summary>
0009 
0010  Implementation:
0011      <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Sep 15 09:45
0015 //         Created:  Mon Sep 15 09:49:08 CEST 2008
0016 //
0017 //
0018 
0019 // system include files
0020 #include <memory>
0021 #include <cctype>
0022 
0023 // user include files
0024 #include "FWCore/Framework/interface/Frameworkfwd.h"
0025 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0026 
0027 #include "FWCore/Framework/interface/Event.h"
0028 #include "FWCore/Framework/interface/MakerMacros.h"
0029 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0030 #include "FWCore/Utilities/interface/InputTag.h"
0031 
0032 #include "DataFormats/Candidate/interface/CandidateFwd.h"
0033 #include "DataFormats/Candidate/interface/Candidate.h"
0034 
0035 #include "DataFormats/Common/interface/Association.h"
0036 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
0037 
0038 #include "PhysicsTools/TagAndProbe/interface/TPTreeFiller.h"
0039 #include "PhysicsTools/TagAndProbe/interface/TagProbePairMaker.h"
0040 
0041 #include <set>
0042 
0043 //
0044 // class decleration
0045 //
0046 
0047 class TagProbeFitTreeProducer : public edm::one::EDAnalyzer<> {
0048 public:
0049   explicit TagProbeFitTreeProducer(const edm::ParameterSet&);
0050   ~TagProbeFitTreeProducer() override;
0051 
0052 private:
0053   void analyze(const edm::Event&, const edm::EventSetup&) override;
0054   void endJob() override;
0055 
0056   //---- MC truth information
0057   /// Is this sample MC?
0058   bool isMC_;
0059   /// Token foran edm::Association<reco::GenParticle> from tags & probes to MC truth
0060   edm::EDGetTokenT<edm::Association<std::vector<reco::GenParticle> > > tagMatchesToken_, probeMatchesToken_;
0061   /// Possible pdgids for the mother. If empty, any truth-matched mu will be considered good
0062   std::set<int32_t> motherPdgId_;
0063   /// Return true if ref is not null and has an ancestor with pdgId inside 'motherPdgId_'
0064   bool checkMother(const reco::GenParticleRef& ref) const;
0065 
0066   //---- Unbiased MC truth information
0067   /// Do we have to compute this
0068   bool makeMCUnbiasTree_;
0069   /// Check mother pdgId in unbiased inefficiency measurement
0070   bool checkMotherInUnbiasEff_;
0071   /// InputTag to the collection of all probes
0072   edm::EDGetTokenT<reco::CandidateView> allProbesToken_;
0073 
0074   /// The object that produces pairs of tags and probes, making any arbitration needed
0075   tnp::TagProbePairMaker tagProbePairMaker_;
0076   /// The object that actually computes variables and fills the tree for T&P
0077   std::unique_ptr<tnp::TPTreeFiller> treeFiller_;
0078   /// The object that actually computes variables and fills the tree for unbiased MC
0079   std::unique_ptr<tnp::BaseTreeFiller> mcUnbiasFiller_;
0080   std::unique_ptr<tnp::BaseTreeFiller> oldTagFiller_;
0081   std::unique_ptr<tnp::BaseTreeFiller> tagFiller_;
0082   std::unique_ptr<tnp::BaseTreeFiller> pairFiller_;
0083   std::unique_ptr<tnp::BaseTreeFiller> mcFiller_;
0084 };
0085 
0086 //
0087 // constructors and destructor
0088 //
0089 TagProbeFitTreeProducer::TagProbeFitTreeProducer(const edm::ParameterSet& iConfig)
0090     : isMC_(iConfig.getParameter<bool>("isMC")),
0091       makeMCUnbiasTree_(isMC_ ? iConfig.getParameter<bool>("makeMCUnbiasTree") : false),
0092       checkMotherInUnbiasEff_(makeMCUnbiasTree_ ? iConfig.getParameter<bool>("checkMotherInUnbiasEff") : false),
0093       tagProbePairMaker_(iConfig, consumesCollector()),
0094       treeFiller_(new tnp::TPTreeFiller(iConfig, consumesCollector())),
0095       oldTagFiller_((iConfig.existsAs<bool>("fillTagTree") && iConfig.getParameter<bool>("fillTagTree"))
0096                         ? new tnp::BaseTreeFiller("tag_tree", iConfig, consumesCollector())
0097                         : nullptr) {
0098   if (isMC_) {
0099     // For mc efficiency we need the MC matches for tags & probes
0100     tagMatchesToken_ =
0101         consumes<edm::Association<std::vector<reco::GenParticle> > >(iConfig.getParameter<edm::InputTag>("tagMatches"));
0102     probeMatchesToken_ = consumes<edm::Association<std::vector<reco::GenParticle> > >(
0103         iConfig.getParameter<edm::InputTag>("probeMatches"));
0104     //.. and the pdgids of the possible mothers
0105     if (iConfig.existsAs<int32_t>("motherPdgId")) {
0106       motherPdgId_.insert(iConfig.getParameter<int32_t>("motherPdgId"));
0107     } else {
0108       std::vector<int32_t> motherIds = iConfig.getParameter<std::vector<int32_t> >("motherPdgId");
0109       motherPdgId_.insert(motherIds.begin(), motherIds.end());
0110     }
0111 
0112     // For unbiased efficiency we also need the collection of all probes
0113     if (makeMCUnbiasTree_) {
0114       allProbesToken_ = consumes<reco::CandidateView>(iConfig.getParameter<edm::InputTag>("allProbes"));
0115       mcUnbiasFiller_ = std::make_unique<tnp::BaseTreeFiller>("mcUnbias_tree", iConfig, consumesCollector());
0116     }
0117   }
0118 
0119   edm::ParameterSet tagPSet;
0120   if (iConfig.existsAs<edm::ParameterSet>("tagVariables"))
0121     tagPSet.addParameter<edm::ParameterSet>("variables", iConfig.getParameter<edm::ParameterSet>("tagVariables"));
0122   if (iConfig.existsAs<edm::ParameterSet>("tagFlags"))
0123     tagPSet.addParameter<edm::ParameterSet>("flags", iConfig.getParameter<edm::ParameterSet>("tagFlags"));
0124   if (!tagPSet.empty()) {
0125     tagFiller_ = std::make_unique<tnp::BaseTreeFiller>(*treeFiller_, tagPSet, consumesCollector(), "tag_");
0126   }
0127   edm::ParameterSet mcPSet;
0128   if (iConfig.existsAs<edm::ParameterSet>("mcVariables"))
0129     mcPSet.addParameter<edm::ParameterSet>("variables", iConfig.getParameter<edm::ParameterSet>("mcVariables"));
0130   if (iConfig.existsAs<edm::ParameterSet>("mcFlags"))
0131     mcPSet.addParameter<edm::ParameterSet>("flags", iConfig.getParameter<edm::ParameterSet>("mcFlags"));
0132   if (!mcPSet.empty()) {
0133     mcFiller_ = std::make_unique<tnp::BaseTreeFiller>(*treeFiller_, mcPSet, consumesCollector(), "mc_");
0134   }
0135   edm::ParameterSet pairPSet;
0136   if (iConfig.existsAs<edm::ParameterSet>("pairVariables"))
0137     pairPSet.addParameter<edm::ParameterSet>("variables", iConfig.getParameter<edm::ParameterSet>("pairVariables"));
0138   if (iConfig.existsAs<edm::ParameterSet>("pairFlags"))
0139     pairPSet.addParameter<edm::ParameterSet>("flags", iConfig.getParameter<edm::ParameterSet>("pairFlags"));
0140   if (!pairPSet.empty()) {
0141     pairFiller_ = std::make_unique<tnp::BaseTreeFiller>(*treeFiller_, pairPSet, consumesCollector(), "pair_");
0142   }
0143 }
0144 
0145 TagProbeFitTreeProducer::~TagProbeFitTreeProducer() {}
0146 
0147 //
0148 // member functions
0149 //
0150 
0151 // ------------ method called to for each event  ------------
0152 void TagProbeFitTreeProducer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0153   using namespace edm;
0154   using namespace std;
0155   Handle<reco::CandidateView> src, allProbes;
0156   Handle<Association<vector<reco::GenParticle> > > tagMatches, probeMatches;
0157   treeFiller_->init(iEvent);  // read out info from the event if needed (external vars, list of passing probes, ...)
0158   if (oldTagFiller_.get())
0159     oldTagFiller_->init(iEvent);
0160   if (tagFiller_.get())
0161     tagFiller_->init(iEvent);
0162   if (pairFiller_.get())
0163     pairFiller_->init(iEvent);
0164   if (mcFiller_.get())
0165     mcFiller_->init(iEvent);
0166 
0167   // on mc we want to load also the MC match info
0168   if (isMC_) {
0169     iEvent.getByToken(tagMatchesToken_, tagMatches);
0170     iEvent.getByToken(probeMatchesToken_, probeMatches);
0171   }
0172 
0173   // get the list of (tag+probe) pairs, performing arbitration
0174   tnp::TagProbePairs pairs = tagProbePairMaker_.run(iEvent);
0175   // loop on them to fill the tree
0176   for (tnp::TagProbePairs::const_iterator it = pairs.begin(), ed = pairs.end(); it != ed; ++it) {
0177     // on mc, fill mc info (on non-mc, let it to 'true', the treeFiller will ignore it anyway
0178     bool mcTrue = false;
0179     float mcMass = 0.f;
0180     if (isMC_) {
0181       reco::GenParticleRef mtag = (*tagMatches)[it->tag], mprobe = (*probeMatches)[it->probe];
0182       mcTrue = checkMother(mtag) && checkMother(mprobe);
0183       if (mcTrue) {
0184         mcMass = (mtag->p4() + mprobe->p4()).mass();
0185         if (mcFiller_.get())
0186           mcFiller_->fill(reco::CandidateBaseRef(mprobe));
0187       }
0188     }
0189     // fill in the variables for this t+p pair
0190     if (tagFiller_.get())
0191       tagFiller_->fill(it->tag);
0192     if (oldTagFiller_.get())
0193       oldTagFiller_->fill(it->tag);
0194     if (pairFiller_.get())
0195       pairFiller_->fill(it->pair);
0196     treeFiller_->fill(it->probe, it->mass, mcTrue, mcMass);
0197   }
0198 
0199   if (isMC_ && makeMCUnbiasTree_) {
0200     // read full collection of probes
0201     iEvent.getByToken(allProbesToken_, allProbes);
0202     // init the tree filler
0203     mcUnbiasFiller_->init(iEvent);
0204     // loop on probes
0205     for (size_t i = 0, n = allProbes->size(); i < n; ++i) {
0206       const reco::CandidateBaseRef& probe = allProbes->refAt(i);
0207       // check mc match, and possibly mother match
0208       reco::GenParticleRef probeMatch = (*probeMatches)[probe];
0209       bool probeOk = checkMotherInUnbiasEff_ ? checkMother(probeMatch) : probeMatch.isNonnull();
0210       // fill the tree only for good ones
0211       if (probeOk)
0212         mcUnbiasFiller_->fill(probe);
0213     }
0214   }
0215 }
0216 
0217 bool TagProbeFitTreeProducer::checkMother(const reco::GenParticleRef& ref) const {
0218   if (ref.isNull())
0219     return false;
0220   if (motherPdgId_.empty())
0221     return true;
0222   if (motherPdgId_.find(abs(ref->pdgId())) != motherPdgId_.end())
0223     return true;
0224   reco::GenParticle::mothers m = ref->motherRefVector();
0225   for (reco::GenParticle::mothers::const_iterator it = m.begin(), e = m.end(); it != e; ++it) {
0226     if (checkMother(*it))
0227       return true;
0228   }
0229   return false;
0230 }
0231 
0232 // ------------ method called once each job just after ending the event loop  ------------
0233 void TagProbeFitTreeProducer::endJob() {
0234   // ask to write the current PSet info into the TTree header
0235   treeFiller_->writeProvenance(edm::getProcessParameterSetContainingModule(moduleDescription()));
0236 }
0237 
0238 //define this as a plug-in
0239 DEFINE_FWK_MODULE(TagProbeFitTreeProducer);