Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #ifndef LHEWEIGHTVALIDATION_H
0002 #define LHEWEIGHTVALIDATION_H
0003 
0004 // framework & common header files
0005 #include "FWCore/Framework/interface/Event.h"
0006 #include "FWCore/Framework/interface/Run.h"
0007 #include "FWCore/Framework/interface/Frameworkfwd.h"
0008 
0009 #include "DataFormats/Common/interface/Handle.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 #include "FWCore/Utilities/interface/InputTag.h"
0012 
0013 //DQM services
0014 #include "DQMServices/Core/interface/DQMStore.h"
0015 #include "FWCore/ServiceRegistry/interface/Service.h"
0016 #include "DQMServices/Core/interface/DQMOneEDAnalyzer.h"
0017 
0018 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
0019 #include "DataFormats/HepMCCandidate/interface/GenParticleFwd.h"
0020 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
0021 #include "DataFormats/JetReco/interface/GenJetCollection.h"
0022 
0023 #include "Validation/EventGenerator/interface/DQMHelper.h"
0024 
0025 #include "SimDataFormats/GeneratorProducts/interface/LHEEventProduct.h"
0026 #include "SimDataFormats/GeneratorProducts/interface/LHERunInfoProduct.h"
0027 
0028 class TH1F;  // forward declaration for ROOT
0029 
0030 class LheWeightValidation : public DQMOneEDAnalyzer<> {
0031 public:
0032   explicit LheWeightValidation(const edm::ParameterSet&);
0033   ~LheWeightValidation() override = default;
0034   void analyze(const edm::Event&, const edm::EventSetup&) override;
0035   void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
0036   void dqmBeginRun(const edm::Run&, const edm::EventSetup&) override;
0037   void dqmEndRun(const edm::Run&, const edm::EventSetup&) override;
0038 
0039 private:
0040   void bookTemplates(DQMHelper& aDqmHelper,
0041                      std::vector<std::unique_ptr<TH1F>>& scaleVar,
0042                      std::vector<std::unique_ptr<TH1F>>& pdfVar,
0043                      std::vector<MonitorElement*>& tmps,
0044                      const std::string& name,
0045                      const std::string& title,
0046                      int nbin,
0047                      float low,
0048                      float high,
0049                      const std::string& xtitle,
0050                      const std::string& ytitle);
0051 
0052   void fillTemplates(std::vector<std::unique_ptr<TH1F>>& scaleVar,
0053                      std::vector<std::unique_ptr<TH1F>>& pdfVar,
0054                      std::vector<MonitorElement*>& tmps,
0055                      float obs);
0056 
0057   void envelop(const std::vector<std::unique_ptr<TH1F>>& var, std::vector<MonitorElement*>& tmps);
0058   void pdfRMS(const std::vector<std::unique_ptr<TH1F>>& var, std::vector<MonitorElement*>& tmps);
0059 
0060   double weight_, orgWgt_;
0061   std::vector<LHEEventProduct::WGT> weights_;
0062 
0063   MonitorElement* nEvt_;
0064   MonitorElement* nlogWgt_;
0065   MonitorElement* wgtVal_;
0066   std::vector<MonitorElement*> leadLepPtTemp_;
0067   std::vector<MonitorElement*> leadLepEtaTemp_;
0068   std::vector<MonitorElement*> jetMultTemp_;
0069   std::vector<MonitorElement*> leadJetPtTemp_;
0070   std::vector<MonitorElement*> leadJetEtaTemp_;
0071 
0072   std::vector<std::unique_ptr<TH1F>> leadLepPtScaleVar_;
0073   std::vector<std::unique_ptr<TH1F>> leadLepPtPdfVar_;
0074   std::vector<std::unique_ptr<TH1F>> leadLepEtaScaleVar_;
0075   std::vector<std::unique_ptr<TH1F>> leadLepEtaPdfVar_;
0076   std::vector<std::unique_ptr<TH1F>> jetMultScaleVar_;
0077   std::vector<std::unique_ptr<TH1F>> jetMultPdfVar_;
0078   std::vector<std::unique_ptr<TH1F>> leadJetPtScaleVar_;
0079   std::vector<std::unique_ptr<TH1F>> leadJetPtPdfVar_;
0080   std::vector<std::unique_ptr<TH1F>> leadJetEtaScaleVar_;
0081   std::vector<std::unique_ptr<TH1F>> leadJetEtaPdfVar_;
0082 
0083   const edm::InputTag lheLabel_;
0084   const edm::EDGetTokenT<reco::GenParticleCollection> genParticleToken_;
0085   const edm::EDGetTokenT<LHEEventProduct> lheEvtToken_;
0086   const edm::EDGetTokenT<LHERunInfoProduct> lheRunToken_;
0087   const edm::EDGetTokenT<reco::GenJetCollection> genJetToken_;
0088 
0089   const bool dumpLHEheader_;
0090   const int leadLepPtNbin_, rapidityNbin_;
0091   const double leadLepPtRange_, leadLepPtCut_, lepEtaCut_, rapidityRange_;
0092   const int nJetsNbin_, jetPtNbin_;
0093   const double jetPtCut_, jetEtaCut_, jetPtRange_;
0094 
0095   const int nScaleVar_;  // including Nominal
0096   const int idxPdfStart_, idxPdfEnd_, nPdfVar_;
0097 };
0098 
0099 #endif