CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ElectronGeneralAnalyzer.cc
Go to the documentation of this file.
1 
3 
5 
14 
21 
22 //#include "CLHEP/Units/GlobalPhysicalConstants.h"
23 //#include "TMath.h"
24 
25 #include <iostream>
26 
27 using namespace reco ;
28 
31  {
32  // collection input tags
33  electronCollection_ = conf.getParameter<edm::InputTag>("ElectronCollection");
34  matchingObjectCollection_ = conf.getParameter<edm::InputTag>("MatchingObjectCollection");
35  trackCollection_ = conf.getParameter<edm::InputTag>("TrackCollection");
36  vertexCollection_ = conf.getParameter<edm::InputTag>("VertexCollection");
37  gsftrackCollection_ = conf.getParameter<edm::InputTag>("GsfTrackCollection");
38  beamSpotTag_ = conf.getParameter<edm::InputTag>("BeamSpot");
39  triggerResults_ = conf.getParameter<edm::InputTag>("TriggerResults");
40 
41 // // for trigger
42 // HLTPathsByName_= conf.getParameter<std::vector<std::string > >("HltPaths");
43 // HLTPathsByIndex_.resize(HLTPathsByName_.size());
44  }
45 
47  {}
48 
50  {
51  h2_ele_beamSpotXvsY = bookH2("beamSpotXvsY","beam spot x vs y",100,-0.2,0.2,100,-0.2,0.2,"x (cm)","y (cm)") ;
52  py_ele_nElectronsVsLs = bookP1("nElectronsVsLs","# gsf electrons vs LS",150,0.,1500.,0.,20.,"LS","<N_{ele}>") ;
53  py_ele_nClustersVsLs = bookP1("nClustersVsLs","# clusters vs LS",150,0.,1500.,0.,100.,"LS","<N_{SC}>") ;
54  py_ele_nGsfTracksVsLs = bookP1("nGsfTracksVsLs","# gsf tracks vs LS",150,0.,1500.,0.,20.,"LS","<N_{GSF tk}>") ;
55  py_ele_nTracksVsLs = bookP1("nTracksVsLs","# tracks vs LS",150,0.,1500.,0.,100.,"LS","<N_{gen tk}>") ;
56  py_ele_nVerticesVsLs = bookP1("nVerticesVsLs","# vertices vs LS",150,0.,1500.,0.,10.,"LS","<N_{vert}>") ;
57  h1_ele_triggers = bookH1("triggers","hlt triggers",256,0.,256.,"HLT bit") ;
58  }
59 
61  {
63  iEvent.getByLabel(electronCollection_,gsfElectrons) ;
65  iEvent.getByLabel(matchingObjectCollection_,recoClusters) ;
67  iEvent.getByLabel(trackCollection_,tracks);
69  iEvent.getByLabel(gsftrackCollection_,gsfTracks);
71  iEvent.getByLabel(vertexCollection_,vertices);
72  edm::Handle<reco::BeamSpot> recoBeamSpotHandle ;
73  iEvent.getByLabel(beamSpotTag_,recoBeamSpotHandle) ;
74  const BeamSpot bs = *recoBeamSpotHandle ;
75 
76  int ievt = iEvent.id().event();
77  int irun = iEvent.id().run();
78  int ils = iEvent.luminosityBlock();
79 
80  edm::LogInfo("ElectronGeneralAnalyzer::analyze")
81  <<"Treating "<<gsfElectrons.product()->size()<<" electrons"
82  <<" from event "<<ievt<<" in run "<<irun<<" and lumiblock "<<ils ;
83 
84  h2_ele_beamSpotXvsY->Fill(bs.position().x(),bs.position().y());
85  py_ele_nElectronsVsLs->Fill(float(ils),(*gsfElectrons).size());
86  py_ele_nClustersVsLs->Fill(float(ils),(*recoClusters).size());
87  py_ele_nGsfTracksVsLs->Fill(float(ils),(*gsfTracks).size());
88  py_ele_nTracksVsLs->Fill(float(ils),(*tracks).size());
89  py_ele_nVerticesVsLs->Fill(float(ils),(*vertices).size());
90 
91  // trigger
93  iEvent.getByLabel(triggerResults_,triggerResults) ;
94  if (triggerResults.isValid())
95  {
96  unsigned int i, n = triggerResults->size() ;
97  for ( i=0 ; i!=n ; ++i )
98  {
99  if (triggerResults->accept(i))
100  { h1_ele_triggers->Fill(float(i)) ; }
101  }
102  }
103  }
104 
RunNumber_t run() const
Definition: EventID.h:42
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:44
MonitorElement * py_ele_nGsfTracksVsLs
int i
Definition: DBlmapReader.cc:9
MonitorElement * bookP1(const std::string &name, const std::string &title, int nchX, double lowX, double highX, double lowY, double highY, const std::string &titleX="", const std::string &titleY="", Option_t *option="E1 P")
MonitorElement * h2_ele_beamSpotXvsY
MonitorElement * bookH1(const std::string &name, const std::string &title, int nchX, double lowX, double highX, const std::string &titleX="", const std::string &titleY="Events", Option_t *option="E1 P")
MonitorElement * py_ele_nVerticesVsLs
MonitorElement * py_ele_nElectronsVsLs
ElectronGeneralAnalyzer(const edm::ParameterSet &conf)
edm::LuminosityBlockNumber_t luminosityBlock() const
Definition: EventBase.h:59
MonitorElement * py_ele_nClustersVsLs
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:243
virtual void analyze(const edm::Event &e, const edm::EventSetup &c)
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
tuple conf
Definition: dbtoconf.py:185
tuple tracks
Definition: testEve_cfg.py:39
T const * product() const
Definition: Handle.h:74
MonitorElement * bookH2(const std::string &name, const std::string &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, const std::string &titleX="", const std::string &titleY="", Option_t *option="COLZ")
edm::EventID id() const
Definition: EventBase.h:56
const Point & position() const
position
Definition: BeamSpot.h:63