Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0002 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/Framework/interface/EventSetup.h"
0005 #include "FWCore/Framework/interface/MakerMacros.h"
0006 
0007 #include "DataFormats/Common/interface/Handle.h"
0008 #include "DataFormats/Common/interface/AssociationVector.h"
0009 #include "DataFormats/Common/interface/RefToBaseProd.h"
0010 
0011 #include "DataFormats/TrackReco/interface/Track.h"
0012 #include "DataFormats/MuonReco/interface/Muon.h"
0013 #include "DataFormats/RecoCandidate/interface/IsoDeposit.h"
0014 #include "DataFormats/RecoCandidate/interface/IsoDepositFwd.h"
0015 #include "DataFormats/RecoCandidate/interface/IsoDepositVetos.h"
0016 
0017 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0018 #include <string>
0019 
0020 #include "DataFormats/Candidate/interface/CandAssociation.h"
0021 #include "FWCore/ServiceRegistry/interface/Service.h"
0022 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0023 #include <TH1.h>
0024 #include <TProfile.h>
0025 
0026 class CandIsoComparer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0027 public:
0028   CandIsoComparer(const edm::ParameterSet&);
0029 
0030   virtual ~CandIsoComparer();
0031 
0032   virtual void analyze(const edm::Event&, const edm::EventSetup&);
0033   virtual void endJob();
0034 
0035 private:
0036   std::string label_;
0037   edm::EDGetTokenT<reco::CandViewDoubleAssociations> src1Token_;
0038   edm::EDGetTokenT<reco::CandViewDoubleAssociations> src2Token_;
0039   TH1 *h1_[2], *h2_[2], *hd_[2];
0040   TProfile *p1_[2], *p2_[2], *pd_[2];
0041 };
0042 
0043 /// constructor with config
0044 CandIsoComparer::CandIsoComparer(const edm::ParameterSet& par)
0045     : src1Token_(consumes<reco::CandViewDoubleAssociations>(par.getParameter<edm::InputTag>("src1"))),
0046       src2Token_(consumes<reco::CandViewDoubleAssociations>(par.getParameter<edm::InputTag>("src2"))) {
0047   label_ = par.getParameter<std::string>("@module_label");
0048 
0049   usesResource(TFileService::kSharedResource);
0050   edm::Service<TFileService> fs;
0051   double max = par.getParameter<double>("max");
0052   double rmax = par.getParameter<double>("rmax");
0053   double dmax = par.getParameter<double>("dmax");
0054   int32_t bins = par.getParameter<int32_t>("bins");
0055 
0056   h1_[0] = fs->make<TH1F>("first", "first", bins, 0, max);
0057   h2_[0] = fs->make<TH1F>("second", "second", bins, 0, max);
0058   hd_[0] = fs->make<TH1F>("diff", "diff", bins, -dmax, dmax);
0059   h1_[1] = fs->make<TH1F>("firstRel", "firstRel", bins, 0, rmax);
0060   h2_[1] = fs->make<TH1F>("secondRel", "secondRel", bins, 0, rmax);
0061   hd_[1] = fs->make<TH1F>("diffRel", "diffRel", bins, -dmax * rmax / max, dmax * rmax / max);
0062 
0063   p1_[0] = fs->make<TProfile>("firstEta", "firstEta", bins, -3.0, 3.0);
0064   p2_[0] = fs->make<TProfile>("secondEta", "secondEta", bins, -3.0, 3.0);
0065   pd_[0] = fs->make<TProfile>("diffEta", "diffEta", bins, -3.0, 3.0);
0066   p1_[1] = fs->make<TProfile>("firstPt", "firstPt", bins * 2, 0, 120.0);
0067   p2_[1] = fs->make<TProfile>("secondPt", "secondPt", bins * 2, 0, 120.0);
0068   pd_[1] = fs->make<TProfile>("diffPt", "diffPt", bins * 2, 0, 120.0);
0069 }
0070 
0071 /// destructor
0072 CandIsoComparer::~CandIsoComparer() {}
0073 
0074 void CandIsoComparer::endJob() { std::cout << "MODULE " << label_ << " RMS DIFF = " << hd_[0]->GetRMS() << std::endl; }
0075 
0076 /// build deposits
0077 void CandIsoComparer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0078   using namespace reco::isodeposit;
0079   edm::Handle<reco::CandViewDoubleAssociations> hDeps1, hDeps2;
0080   iEvent.getByToken(src1Token_, hDeps1);
0081   iEvent.getByToken(src2Token_, hDeps2);
0082   for (size_t dep = 0; dep < hDeps1->size(); ++dep) {
0083     const reco::CandidateBaseRef& cand = hDeps1->key(dep);
0084     double iso1 = (*hDeps1)[cand];
0085     double iso2 = (*hDeps2)[cand];
0086     h1_[0]->Fill(iso1);
0087     h2_[0]->Fill(iso2);
0088     hd_[0]->Fill(iso1 - iso2);
0089     h1_[1]->Fill(iso1 / cand->pt());
0090     h2_[1]->Fill(iso2 / cand->pt());
0091     hd_[1]->Fill((iso1 - iso2) / cand->pt());
0092     p1_[0]->Fill(cand->eta(), iso1);
0093     p2_[0]->Fill(cand->eta(), iso2);
0094     pd_[0]->Fill(cand->eta(), iso1 - iso2);
0095     p1_[1]->Fill(cand->pt(), iso1);
0096     p2_[1]->Fill(cand->pt(), iso2);
0097     pd_[1]->Fill(cand->pt(), iso1 - iso2);
0098   }
0099 }
0100 
0101 DEFINE_FWK_MODULE(CandIsoComparer);