CMS 3D CMS Logo

EgammaBasicClusters.cc
Go to the documentation of this file.
2 
5 
8 
10 
12 
14  barrelBasicClusterCollection_(consumes<reco::BasicClusterCollection>(ps.getParameter<edm::InputTag>("barrelBasicClusterCollection"))),
15  endcapBasicClusterCollection_(consumes<reco::BasicClusterCollection>(ps.getParameter<edm::InputTag>("endcapBasicClusterCollection"))),
16  hsSize_(ps, "Size"),
17  hsNumRecHits_(ps, "NumRecHits"),
18  hsET_(ps, "ET"),
19  hsEta_(ps, "Eta"),
20  hsPhi_(ps, "Phi"),
21  hsR_(ps, "R"),
22  hist_EB_BC_Size_(0),
23  hist_EE_BC_Size_(0),
24  hist_EB_BC_NumRecHits_(0),
25  hist_EE_BC_NumRecHits_(0),
26  hist_EB_BC_ET_(0),
27  hist_EE_BC_ET_(0),
28  hist_EB_BC_Eta_(0),
29  hist_EE_BC_Eta_(0),
30  hist_EB_BC_Phi_(0),
31  hist_EE_BC_Phi_(0),
32  hist_EB_BC_ET_vs_Eta_(0),
33  hist_EB_BC_ET_vs_Phi_(0),
34  hist_EE_BC_ET_vs_Eta_(0),
35  hist_EE_BC_ET_vs_Phi_(0),
36  hist_EE_BC_ET_vs_R_(0)
37 {
38 }
39 
41 {
42 }
43 
44 void
46 {
47  _ibooker.setCurrentFolder("EcalClusterV/EcalBasicClusters/");
48 
50  = _ibooker.book1D("hist_EB_BC_Size_","# Basic Clusters in Barrel",
53  = _ibooker.book1D("hist_EE_BC_Size_","# Basic Clusters in Endcap",
55 
57  = _ibooker.book1D("hist_EB_BC_NumRecHits_","# of RecHits in Basic Clusters in Barrel",
60  = _ibooker.book1D("hist_EE_BC_NumRecHits_","# of RecHits in Basic Clusters in Endcap",
62 
64  = _ibooker.book1D("hist_EB_BC_ET_","ET of Basic Clusters in Barrel",
67  = _ibooker.book1D("hist_EE_BC_ET_","ET of Basic Clusters in Endcap",
69 
71  = _ibooker.book1D("hist_EB_BC_Eta_","Eta of Basic Clusters in Barrel",
74  = _ibooker.book1D("hist_EE_BC_Eta_","Eta of Basic Clusters in Endcap",
76 
78  = _ibooker.book1D("hist_EB_BC_Phi_","Phi of Basic Clusters in Barrel",
81  = _ibooker.book1D("hist_EE_BC_Phi_","Phi of Basic Clusters in Endcap",
83 
84 
86  = _ibooker.book2D( "hist_EB_BC_ET_vs_Eta_", "Basic Cluster ET versus Eta in Barrel",
89 
91  = _ibooker.book2D( "hist_EB_BC_ET_vs_Phi_", "Basic Cluster ET versus Phi in Barrel",
94 
96  = _ibooker.book2D( "hist_EE_BC_ET_vs_Eta_", "Basic Cluster ET versus Eta in Endcap",
99 
101  = _ibooker.book2D( "hist_EE_BC_ET_vs_Phi_", "Basic Cluster ET versus Phi in Endcap",
104 
106  = _ibooker.book2D( "hist_EE_BC_ET_vs_R_", "Basic Cluster ET versus Radius in Endcap",
108  hsR_.bins, hsR_.min, hsR_.max );
109 }
110 
111 void
113 {
114  edm::Handle<reco::BasicClusterCollection> pBarrelBasicClusters;
115  evt.getByToken(barrelBasicClusterCollection_, pBarrelBasicClusters);
116  if (!pBarrelBasicClusters.isValid()) {
117 
118  Labels l;
120  edm::LogError("EgammaBasicClusters") << "Error! can't get collection with label "
121  << l.module;
122  }
123 
124  const reco::BasicClusterCollection* barrelBasicClusters = pBarrelBasicClusters.product();
125  hist_EB_BC_Size_->Fill(barrelBasicClusters->size());
126 
127  for(reco::BasicClusterCollection::const_iterator aClus = barrelBasicClusters->begin();
128  aClus != barrelBasicClusters->end(); aClus++)
129  {
130  hist_EB_BC_NumRecHits_->Fill(aClus->size());
131  hist_EB_BC_ET_->Fill(aClus->energy()/std::cosh(aClus->position().eta()));
132  hist_EB_BC_Eta_->Fill(aClus->position().eta());
133  hist_EB_BC_Phi_->Fill(aClus->position().phi());
134 
135  hist_EB_BC_ET_vs_Eta_->Fill( aClus->energy()/std::cosh(aClus->position().eta()), aClus->eta() );
136  hist_EB_BC_ET_vs_Phi_->Fill( aClus->energy()/std::cosh(aClus->position().eta()), aClus->phi() );
137  }
138 
139  edm::Handle<reco::BasicClusterCollection> pEndcapBasicClusters;
140 
141  evt.getByToken(endcapBasicClusterCollection_, pEndcapBasicClusters);
142  if (!pEndcapBasicClusters.isValid()) {
143 
144  Labels l;
146  edm::LogError("EgammaBasicClusters") << "Error! can't get collection with label "
147  << l.module;
148  }
149 
150  const reco::BasicClusterCollection* endcapBasicClusters = pEndcapBasicClusters.product();
151  hist_EE_BC_Size_->Fill(endcapBasicClusters->size());
152 
153  for(reco::BasicClusterCollection::const_iterator aClus = endcapBasicClusters->begin();
154  aClus != endcapBasicClusters->end(); aClus++)
155  {
156  hist_EE_BC_NumRecHits_->Fill(aClus->size());
157  hist_EE_BC_ET_->Fill(aClus->energy()/std::cosh(aClus->position().eta()));
158  hist_EE_BC_Eta_->Fill(aClus->position().eta());
159  hist_EE_BC_Phi_->Fill(aClus->position().phi());
160 
161  hist_EE_BC_ET_vs_Eta_->Fill( aClus->energy()/std::cosh(aClus->position().eta()), aClus->eta() );
162  hist_EE_BC_ET_vs_Phi_->Fill( aClus->energy()/std::cosh(aClus->position().eta()), aClus->phi() );
163  hist_EE_BC_ET_vs_R_->Fill( aClus->energy()/std::cosh(aClus->position().eta()),
164  std::sqrt( std::pow(aClus->x(),2) + std::pow(aClus->y(),2) ) );
165 
166  }
167 }
MonitorElement * hist_EE_BC_NumRecHits_
MonitorElement * hist_EB_BC_Eta_
MonitorElement * hist_EE_BC_ET_vs_R_
MonitorElement * hist_EE_BC_Phi_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
int bins
Definition: HistSpec.h:9
MonitorElement * hist_EB_BC_NumRecHits_
void analyze(const edm::Event &, const edm::EventSetup &) override
MonitorElement * hist_EE_BC_ET_vs_Phi_
MonitorElement * hist_EB_BC_ET_
double max
Definition: HistSpec.h:8
MonitorElement * hist_EB_BC_ET_vs_Eta_
MonitorElement * hist_EE_BC_ET_vs_Eta_
void Fill(long long x)
double min
Definition: HistSpec.h:7
MonitorElement * hist_EB_BC_ET_vs_Phi_
MonitorElement * hist_EE_BC_Eta_
T sqrt(T t)
Definition: SSEVec.h:18
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
bool isValid() const
Definition: HandleBase.h:74
char const * module
Definition: ProductLabels.h:5
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
EgammaBasicClusters(const edm::ParameterSet &)
MonitorElement * hist_EB_BC_Phi_
edm::EDGetTokenT< reco::BasicClusterCollection > barrelBasicClusterCollection_
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
T const * product() const
Definition: Handle.h:81
MonitorElement * hist_EE_BC_ET_
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
std::vector< BasicCluster > BasicClusterCollection
collection of BasicCluster objects
void labelsForToken(EDGetToken iToken, Labels &oLabels) const
MonitorElement * hist_EB_BC_Size_
fixed size matrix
HLT enums.
MonitorElement * hist_EE_BC_Size_
edm::EDGetTokenT< reco::BasicClusterCollection > endcapBasicClusterCollection_
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
Definition: Run.h:42