Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "DataFormats/Common/interface/Handle.h"
0002 #include "DataFormats/MuonReco/interface/Muon.h"
0003 #include "PhysicsTools/PatExamples/interface/PatMuonAnalyzer.h"
0004 
0005 /// default constructor
0006 PatMuonAnalyzer::PatMuonAnalyzer(const edm::ParameterSet& cfg, TFileDirectory& fs)
0007     : edm::BasicAnalyzer::BasicAnalyzer(cfg, fs), muons_(cfg.getParameter<edm::InputTag>("muons")) {
0008   hists_["muonPt"] = fs.make<TH1F>("muonPt", "pt", 100, 0., 300.);
0009   hists_["muonEta"] = fs.make<TH1F>("muonEta", "eta", 100, -3., 3.);
0010   hists_["muonPhi"] = fs.make<TH1F>("muonPhi", "phi", 100, -5., 5.);
0011 }
0012 PatMuonAnalyzer::PatMuonAnalyzer(const edm::ParameterSet& cfg, TFileDirectory& fs, edm::ConsumesCollector&& iC)
0013     : edm::BasicAnalyzer::BasicAnalyzer(cfg, fs),
0014       muons_(cfg.getParameter<edm::InputTag>("muons")),
0015       muonsToken_(iC.consumes<std::vector<pat::Muon> >(muons_)) {
0016   hists_["muonPt"] = fs.make<TH1F>("muonPt", "pt", 100, 0., 300.);
0017   hists_["muonEta"] = fs.make<TH1F>("muonEta", "eta", 100, -3., 3.);
0018   hists_["muonPhi"] = fs.make<TH1F>("muonPhi", "phi", 100, -5., 5.);
0019 }
0020 
0021 /// everything that needs to be done during the event loop
0022 void PatMuonAnalyzer::analyze(const edm::EventBase& event) {
0023   // define what muon you are using; this is necessary as FWLite is not
0024   // capable of reading edm::Views
0025   using pat::Muon;
0026 
0027   // Handle to the muon collection
0028   edm::Handle<std::vector<Muon> > muons;
0029   event.getByLabel(muons_, muons);
0030 
0031   // loop muon collection and fill histograms
0032   for (std::vector<Muon>::const_iterator mu1 = muons->begin(); mu1 != muons->end(); ++mu1) {
0033     hists_["muonPt"]->Fill(mu1->pt());
0034     hists_["muonEta"]->Fill(mu1->eta());
0035     hists_["muonPhi"]->Fill(mu1->phi());
0036   }
0037 }