CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EBRecoSummary.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: EBRecoSummary
4 // Class: EBRecoSummary
5 // Original Author: Martina Malberti
6 //
7 // system include files
8 #include <memory>
9 
10 // user include files
13 
18 
22 
23 
32 
46 
49 
52 
57 
58 #include "TVector3.h"
59 
60 #include <iostream>
61 #include <cmath>
62 #include <fstream>
63 #include <string>
64 
65 //
66 // constructors and destructor
67 //
69 {
70 
71  prefixME_ = ps.getUntrackedParameter<std::string>("prefixME", "");
72 
73  //now do what ever initialization is needed
74  recHitCollection_EB_ = ps.getParameter<edm::InputTag>("recHitCollection_EB");
75  redRecHitCollection_EB_ = ps.getParameter<edm::InputTag>("redRecHitCollection_EB");
76  basicClusterCollection_EB_ = ps.getParameter<edm::InputTag>("basicClusterCollection_EB");
77  superClusterCollection_EB_ = ps.getParameter<edm::InputTag>("superClusterCollection_EB");
78 
79  ethrEB_ = ps.getParameter<double>("ethrEB");
80 
81  scEtThrEB_ = ps.getParameter<double>("scEtThrEB");
82 
83  // DQM Store -------------------
85 
86  // Monitor Elements (ex THXD)
87  dqmStore_->setCurrentFolder(prefixME_ + "/EBRecoSummary"); // to organise the histos in folders
88 
89  // ReducedRecHits ----------------------------------------------
90  // ... barrel
91  h_redRecHits_EB_recoFlag = dqmStore_->book1D("redRecHits_EB_recoFlag","redRecHits_EB_recoFlag",16,-0.5,15.5);
92 
93  // RecHits ----------------------------------------------
94  // ... barrel
95  h_recHits_EB_energyMax = dqmStore_->book1D("recHits_EB_energyMax","recHitsEB_energyMax",110,-10,100);
96  h_recHits_EB_Chi2 = dqmStore_->book1D("recHits_EB_Chi2","recHits_EB_Chi2",200,0,100);
97  h_recHits_EB_time = dqmStore_->book1D("recHits_EB_time","recHits_EB_time",200,-50,50);
98  h_recHits_EB_E1oE4 = dqmStore_->book1D("recHits_EB_E1oE4","recHitsEB_E1oE4",150, 0, 1.5);
99  h_recHits_EB_recoFlag = dqmStore_->book1D("recHits_EB_recoFlag","recHits_EB_recoFlag",16,-0.5,15.5);
100 
101  // Basic Clusters ----------------------------------------------
102  // ... associated barrel rec hits
103  h_basicClusters_recHits_EB_recoFlag = dqmStore_->book1D("basicClusters_recHits_EB_recoFlag","basicClusters_recHits_EB_recoFlag",16,-0.5,15.5);
104 
105  // Super Clusters ----------------------------------------------
106  // ... barrel
107  h_superClusters_EB_nBC = dqmStore_->book1D("superClusters_EB_nBC","superClusters_EB_nBC",100,0.,100.);
108  h_superClusters_EB_E1oE4 = dqmStore_->book1D("superClusters_EB_E1oE4","superClusters_EB_E1oE4",150,0,1.5);
109 
110  h_superClusters_eta = dqmStore_->book1D("superClusters_eta","superClusters_eta",150,-3.,3.);
111  h_superClusters_EB_phi = dqmStore_->book1D("superClusters_EB_phi","superClusters_EB_phi",360,-3.1415927,3.1415927);
112 
113 }
114 
115 
116 
118 {
119  // do anything here that needs to be done at desctruction time
120  // (e.g. close files, deallocate resources etc.)
121 }
122 
123 
124 //
125 // member functions
126 //
127 
128 // ------------ method called to for each event ------------
129 void EBRecoSummary::analyze(const edm::Event& ev, const edm::EventSetup& iSetup)
130 {
131 
132  //Get the magnetic field
133  edm::ESHandle<MagneticField> theMagField;
134  iSetup.get<IdealMagneticFieldRecord>().get(theMagField);
135 
136  // calo topology
137  edm::ESHandle<CaloTopology> pTopology;
138  iSetup.get<CaloTopologyRecord>().get(pTopology);
139  const CaloTopology *topology = pTopology.product();
140 
141  // --- REDUCED REC HITS -------------------------------------------------------------------------------------
143  ev.getByLabel( redRecHitCollection_EB_, redRecHitsEB );
144  const EcalRecHitCollection* theBarrelEcalredRecHits = redRecHitsEB.product () ;
145  if ( ! redRecHitsEB.isValid() ) {
146  edm::LogWarning("EBRecoSummary") << "redRecHitsEB not found";
147  }
148 
149  for ( EcalRecHitCollection::const_iterator itr = theBarrelEcalredRecHits->begin () ;
150  itr != theBarrelEcalredRecHits->end () ;++itr)
151  {
152 
153  h_redRecHits_EB_recoFlag->Fill( itr -> recoFlag() );
154 
155  }
156 
157  // --- REC HITS -------------------------------------------------------------------------------------
158 
159  // ... barrel
161  ev.getByLabel( recHitCollection_EB_, recHitsEB );
162  const EcalRecHitCollection* theBarrelEcalRecHits = recHitsEB.product () ;
163  if ( ! recHitsEB.isValid() ) {
164  edm::LogWarning("EBRecoSummary") << "recHitsEB not found";
165  }
166 
167  float maxRecHitEnergyEB = -999.;
168 
169  EBDetId ebid_MrecHitEB;
170 
171  for ( EcalRecHitCollection::const_iterator itr = theBarrelEcalRecHits->begin () ;
172  itr != theBarrelEcalRecHits->end () ;++itr)
173  {
174 
175  EBDetId ebid( itr -> id() );
176 
177  h_recHits_EB_recoFlag -> Fill( itr -> recoFlag() );
178 
179  // max E rec hit
180  if (itr -> energy() > maxRecHitEnergyEB ){
181  maxRecHitEnergyEB = itr -> energy() ;
182  }
183 
184  if ( itr -> energy() > ethrEB_ ){
185  h_recHits_EB_Chi2 -> Fill( itr -> chi2() );
186  h_recHits_EB_time -> Fill( itr -> time() );
187  }
188 
189  float R4 = EcalTools::swissCross( ebid, *theBarrelEcalRecHits, 0. );
190 
191  if ( itr -> energy() > 3. && abs(ebid.ieta())!=85 ) h_recHits_EB_E1oE4-> Fill( R4 );
192 
193  }
194 
195  h_recHits_EB_energyMax -> Fill( maxRecHitEnergyEB );
196 
197  //--- BASIC CLUSTERS --------------------------------------------------------------
198 
199  // ... barrel
201  ev.getByLabel( basicClusterCollection_EB_, basicClusters_EB_h );
202  const reco::BasicClusterCollection* theBarrelBasicClusters = basicClusters_EB_h.product () ;
203  if ( ! basicClusters_EB_h.isValid() ) {
204  edm::LogWarning("EBRecoSummary") << "basicClusters_EB_h not found";
205  }
206 
207  for (reco::BasicClusterCollection::const_iterator itBC = theBarrelBasicClusters->begin();
208  itBC != theBarrelBasicClusters->end(); ++itBC ) {
209 
210  //Get the associated RecHits
211  const std::vector<std::pair<DetId,float> > & hits= itBC->hitsAndFractions();
212  for (std::vector<std::pair<DetId,float> > ::const_iterator rh = hits.begin(); rh!=hits.end(); ++rh){
213 
214  EBRecHitCollection::const_iterator itrechit = theBarrelEcalRecHits->find((*rh).first);
215  if (itrechit==theBarrelEcalRecHits->end()) continue;
216  h_basicClusters_recHits_EB_recoFlag -> Fill ( itrechit -> recoFlag() );
217 
218  }
219 
220  }
221 
222  // Super Clusters
223  // ... barrel
225  ev.getByLabel( superClusterCollection_EB_, superClusters_EB_h );
226  const reco::SuperClusterCollection* theBarrelSuperClusters = superClusters_EB_h.product () ;
227  if ( ! superClusters_EB_h.isValid() ) {
228  edm::LogWarning("EBRecoSummary") << "superClusters_EB_h not found";
229  }
230 
231  for (reco::SuperClusterCollection::const_iterator itSC = theBarrelSuperClusters->begin();
232  itSC != theBarrelSuperClusters->end(); ++itSC ) {
233 
234  double scEt = itSC -> energy() * sin(2.*atan( exp(- itSC->position().eta() )));
235 
236  if (scEt < scEtThrEB_ ) continue;
237 
238  h_superClusters_EB_nBC -> Fill( itSC -> clustersSize());
239  h_superClusters_eta -> Fill( itSC -> eta() );
240  h_superClusters_EB_phi -> Fill( itSC -> phi() );
241 
242  float E1 = EcalClusterTools::eMax ( *itSC, theBarrelEcalRecHits);
243  float E4 = EcalClusterTools::eTop ( *itSC, theBarrelEcalRecHits, topology )+
244  EcalClusterTools::eRight ( *itSC, theBarrelEcalRecHits, topology )+
245  EcalClusterTools::eBottom( *itSC, theBarrelEcalRecHits, topology )+
246  EcalClusterTools::eLeft ( *itSC, theBarrelEcalRecHits, topology );
247 
248  if ( E1 > 3. ) h_superClusters_EB_E1oE4 -> Fill( 1.- E4/E1);
249 
250  }
251 
252 }
253 
254 
255 // ------------ method called once each job just before starting event loop ------------
256  void
258 {
259 }
260 
261 // ------------ method called once each job just after ending the event loop ------------
262 void
264 {}
265 
266 // ----------additional functions-------------------
267 
268 void EBRecoSummary::convxtalid(Int_t &nphi,Int_t &neta)
269 {
270  // Barrel only
271  // Output nphi 0...359; neta 0...84; nside=+1 (for eta>0), or 0 (for eta<0).
272  // neta will be [-85,-1] , or [0,84], the minus sign indicates the z<0 side.
273 
274  if(neta > 0) neta -= 1;
275  if(nphi > 359) nphi=nphi-360;
276 
277 } //end of convxtalid
278 
279 int EBRecoSummary::diff_neta_s(Int_t neta1, Int_t neta2){
280  Int_t mdiff;
281  mdiff=(neta1-neta2);
282  return mdiff;
283 }
284 
285 // Calculate the distance in xtals taking into account the periodicity of the Barrel
286 int EBRecoSummary::diff_nphi_s(Int_t nphi1,Int_t nphi2) {
287  Int_t mdiff;
288  if(abs(nphi1-nphi2) < (360-abs(nphi1-nphi2))) {
289  mdiff=nphi1-nphi2;
290  }
291  else {
292  mdiff=360-abs(nphi1-nphi2);
293  if(nphi1>nphi2) mdiff=-mdiff;
294  }
295  return mdiff;
296 }
297 
298 //define this as a plug-in
MonitorElement * h_recHits_EB_time
Definition: EBRecoSummary.h:75
static float eBottom(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
T getParameter(std::string const &) const
MonitorElement * h_recHits_EB_energyMax
Definition: EBRecoSummary.h:73
std::string prefixME_
Definition: EBRecoSummary.h:63
T getUntrackedParameter(std::string const &, T const &) const
MonitorElement * h_recHits_EB_Chi2
Definition: EBRecoSummary.h:74
MonitorElement * h_superClusters_EB_nBC
Definition: EBRecoSummary.h:84
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
MonitorElement * h_superClusters_EB_E1oE4
Definition: EBRecoSummary.h:85
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
std::vector< EcalRecHit >::const_iterator const_iterator
#define abs(x)
Definition: mlp_lapack.h:159
DQMStore * dqmStore_
Definition: EBRecoSummary.h:61
MonitorElement * h_recHits_EB_E1oE4
Definition: EBRecoSummary.h:76
T eta() const
edm::InputTag superClusterCollection_EB_
Definition: EBRecoSummary.h:96
MonitorElement * h_recHits_EB_recoFlag
Definition: EBRecoSummary.h:77
MonitorElement * h_superClusters_eta
Definition: EBRecoSummary.h:87
int diff_neta_s(int, int)
void Fill(long long x)
virtual void beginJob()
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
EBRecoSummary(const edm::ParameterSet &)
virtual void analyze(const edm::Event &, const edm::EventSetup &)
virtual void endJob()
int ieta() const
get the crystal ieta
Definition: EBDetId.h:44
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
const_iterator end() const
static float eRight(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
static float eTop(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
static float eMax(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits)
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
edm::InputTag recHitCollection_EB_
Definition: EBRecoSummary.h:93
std::vector< BasicCluster > BasicClusterCollection
collection of BasicCluster objects
T const * product() const
Definition: Handle.h:74
void convxtalid(int &, int &)
static float swissCross(const DetId &id, const EcalRecHitCollection &recHits, float recHitThreshold, bool avoidIeta85=true)
the good old 1-e4/e1. Ignore hits below recHitThreshold
Definition: EcalTools.cc:12
edm::InputTag redRecHitCollection_EB_
Definition: EBRecoSummary.h:94
iterator find(key_type k)
MonitorElement * h_superClusters_EB_phi
Definition: EBRecoSummary.h:88
static float eLeft(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
MonitorElement * h_redRecHits_EB_recoFlag
Definition: EBRecoSummary.h:69
edm::InputTag basicClusterCollection_EB_
Definition: EBRecoSummary.h:95
MonitorElement * h_basicClusters_recHits_EB_recoFlag
Definition: EBRecoSummary.h:80
int diff_nphi_s(int, int)
const_iterator begin() const
Definition: DDAxes.h:10