Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:47

0001 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0002 #include "DataFormats/PatCandidates/interface/Jet.h"
0003 #include "DQMOffline/RecoB/interface/JetTagPlotter.h"
0004 #include "FWCore/Framework/interface/Frameworkfwd.h"
0005 #include "FWCore/Framework/interface/MakerMacros.h"
0006 #include "FWCore/Framework/interface/Event.h"
0007 
0008 /** \class MiniAODTaggerAnalyzer
0009  *
0010  *  Tagger analyzer to run on MiniAOD
0011  *
0012  */
0013 
0014 class MiniAODTaggerAnalyzer : public DQMEDAnalyzer {
0015 public:
0016   explicit MiniAODTaggerAnalyzer(const edm::ParameterSet& pSet);
0017   ~MiniAODTaggerAnalyzer() override = default;
0018 
0019   void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
0020 
0021 private:
0022   void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
0023   typedef std::vector<std::string> vstring;
0024 
0025   // using JetTagPlotter object for all the hard work ;)
0026   std::unique_ptr<JetTagPlotter> jetTagPlotter_;
0027 
0028   const edm::EDGetTokenT<std::vector<pat::Jet> > jetToken_;
0029   const edm::ParameterSet discrParameters_;
0030 
0031   const std::string folder_;
0032   const vstring discrNumerator_;
0033   const vstring discrDenominator_;
0034 
0035   const int mclevel_;
0036   const bool doCTagPlots_;
0037   const bool dodifferentialPlots_;
0038   const double discrCut_;
0039 
0040   const bool etaActive_;
0041   const double etaMin_;
0042   const double etaMax_;
0043   const bool ptActive_;
0044   const double ptMin_;
0045   const double ptMax_;
0046 };
0047 
0048 MiniAODTaggerAnalyzer::MiniAODTaggerAnalyzer(const edm::ParameterSet& pSet)
0049     : jetToken_(consumes<std::vector<pat::Jet> >(pSet.getParameter<edm::InputTag>("JetTag"))),
0050       discrParameters_(pSet.getParameter<edm::ParameterSet>("parameters")),
0051       folder_(pSet.getParameter<std::string>("folder")),
0052       discrNumerator_(pSet.getParameter<vstring>("numerator")),
0053       discrDenominator_(pSet.getParameter<vstring>("denominator")),
0054       mclevel_(pSet.getParameter<int>("MClevel")),
0055       doCTagPlots_(pSet.getParameter<bool>("CTagPlots")),
0056       dodifferentialPlots_(pSet.getParameter<bool>("differentialPlots")),
0057       discrCut_(pSet.getParameter<double>("discrCut")),
0058       etaActive_(pSet.getParameter<bool>("etaActive")),
0059       etaMin_(pSet.getParameter<double>("etaMin")),
0060       etaMax_(pSet.getParameter<double>("etaMax")),
0061       ptActive_(pSet.getParameter<bool>("ptActive")),
0062       ptMin_(pSet.getParameter<double>("ptMin")),
0063       ptMax_(pSet.getParameter<double>("ptMax")) {}
0064 
0065 void MiniAODTaggerAnalyzer::bookHistograms(DQMStore::IBooker& ibook, edm::Run const& run, edm::EventSetup const& es) {
0066   jetTagPlotter_ = std::make_unique<JetTagPlotter>(folder_,
0067                                                    EtaPtBin(etaActive_, etaMin_, etaMax_, ptActive_, ptMin_, ptMax_),
0068                                                    discrParameters_,
0069                                                    mclevel_,
0070                                                    false,
0071                                                    ibook,
0072                                                    doCTagPlots_,
0073                                                    dodifferentialPlots_,
0074                                                    discrCut_);
0075 }
0076 
0077 void MiniAODTaggerAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0078   edm::Handle<std::vector<pat::Jet> > jetCollection;
0079   iEvent.getByToken(jetToken_, jetCollection);
0080 
0081   // Loop over the pat::Jets
0082   for (std::vector<pat::Jet>::const_iterator jet = jetCollection->begin(); jet != jetCollection->end(); ++jet) {
0083     // apply basic jet cuts
0084     if (jet->pt() > ptMin_ && std::abs(jet->eta()) < etaMax_) {
0085       // fill numerator
0086       float numerator = 0;
0087       for (const auto& discrLabel : discrNumerator_) {
0088         numerator += jet->bDiscriminator(discrLabel);
0089       }
0090 
0091       // fill denominator
0092       float denominator;
0093       if (discrDenominator_.empty()) {
0094         denominator = 1;  // no division performed
0095       } else {
0096         denominator = 0;
0097 
0098         for (const auto& discrLabel : discrDenominator_) {
0099           denominator += jet->bDiscriminator(discrLabel);
0100         }
0101       }
0102 
0103       const float jec = 1.;  // JEC not implemented!
0104 
0105       // only add to histograms when discriminator values are valid
0106       if (numerator >= 0 && denominator > 0) {
0107         const reco::Jet& recoJet = *jet;
0108         if (jetTagPlotter_->etaPtBin().inBin(recoJet, jec)) {
0109           jetTagPlotter_->analyzeTag(recoJet, jec, numerator / denominator, jet->partonFlavour());
0110         }
0111       }
0112     }
0113   }
0114 
0115   // fill JetMultiplicity (once per event)
0116   if (mclevel_ > 0) {
0117     jetTagPlotter_->analyzeTag(1.);
0118   } else {
0119     jetTagPlotter_->analyzeTag();
0120   }
0121 }
0122 
0123 //define this as a plug-in
0124 DEFINE_FWK_MODULE(MiniAODTaggerAnalyzer);