CMS 3D CMS Logo

HGCalRecHitValidation.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HGCalRecHitValidation
4 // Class: HGCalRecHitValidation
5 //
11 //
12 // Original Author: Raman Khurana
14 // Created: Sunday, 17th Augst 2014 11:30:15 GMT
15 
16 // system include files
25 
34 
35 #include "CLHEP/Units/GlobalSystemOfUnits.h"
36 #include "TVector3.h"
37 #include <cmath>
38 
39 
42  nameDetector_(iConfig.getParameter<std::string>("DetectorName")),
43  verbosity_(iConfig.getUntrackedParameter<int>("Verbosity",0)) {
44  auto temp = iConfig.getParameter<edm::InputTag>("RecHitSource");
45  ifHCAL_ = iConfig.getParameter<bool>("ifHCAL");
46  if( nameDetector_ == "HGCalEESensitive" ||
47  nameDetector_ == "HGCalHESiliconSensitive" ||
48  nameDetector_ == "HGCalHEScintillatorSensitive" ) {
49  recHitSource_ = consumes<HGCRecHitCollection>(temp);
50  } else if ( nameDetector_ == "HCal" ) {
51  if (ifHCAL_) recHitSource_ = consumes<HBHERecHitCollection>(temp);
52  else recHitSource_ = consumes<HGChebRecHitCollection>(temp);
53  } else {
54  throw cms::Exception("BadHGCRecHitSource")
55  << "HGCal DetectorName given as " << nameDetector_ << " must be: "
56  << "\"HGCalHESiliconSensitive\", \"HGCalHESiliconSensitive\", "
57  << "\"HGCalHEScintillatorSensitive\", or \"HCal\"!";
58  }
59 }
60 
61 
63 
65  const edm::EventSetup& iSetup) {
66  OccupancyMap_plus.clear();
67  OccupancyMap_minus.clear();
68 
69  bool ok(true);
70  unsigned int ntot(0), nused(0);
71  if (nameDetector_ == "HCal") {
73  iSetup.get<CaloGeometryRecord>().get(geom);
74  if (!geom.isValid()) edm::LogWarning("HGCalValidation") << "Cannot get valid HGCalGeometry Object for " << nameDetector_;
75  const CaloGeometry* geom0 = geom.product();
76 
77  if (ifHCAL_) {
79  iEvent.getByToken(recHitSource_, hbhecoll);
80  if (hbhecoll.isValid()) {
81  if (verbosity_>0)
82  edm::LogInfo("HGCalValidation") << nameDetector_ << " with "
83  << hbhecoll->size() << " element(s)";
84  for (HBHERecHitCollection::const_iterator it=hbhecoll->begin();
85  it != hbhecoll->end(); ++it) {
86  DetId detId = it->id();
87  ntot++;
88  if (detId.subdetId() == HcalEndcap) {
89  nused++;
90  int layer = HcalDetId(detId).depth();
91  recHitValidation(detId, layer, geom0, it);
92  }
93  }
94  } else {
95  ok = false;
96  edm::LogWarning("HGCalValidation") << "HBHERecHitCollection Handle does not exist !!!";
97  }
98  } else {
100  iEvent.getByToken(recHitSource_, hbhecoll);
101  if (hbhecoll.isValid()) {
102  if (verbosity_>0)
103  edm::LogInfo("HGCalValidation") << nameDetector_ << " with "
104  << hbhecoll->size() << " element(s)";
105  for (HGChebRecHitCollection::const_iterator it=hbhecoll->begin();
106  it != hbhecoll->end(); ++it) {
107  DetId detId = it->id();
108  ntot++; nused++;
109  int layer = HcalDetId(detId).depth();
110  recHitValidation(detId, layer, geom0, it);
111  }
112  } else {
113  ok = false;
114  edm::LogWarning("HGCalValidation") << "HGChebRecHitCollection Handle does not exist !!!";
115  }
116  }
117  } else {
119  iSetup.get<IdealGeometryRecord>().get(nameDetector_, geom);
120  if (!geom.isValid()) edm::LogWarning("HGCalValidation") << "Cannot get valid HGCalGeometry Object for " << nameDetector_;
121  const HGCalGeometry* geom0 = geom.product();
122 
123  edm::Handle<HGCRecHitCollection> theRecHitContainers;
124  iEvent.getByToken(recHitSource_, theRecHitContainers);
125  if (theRecHitContainers.isValid()) {
126  if (verbosity_>0)
127  edm::LogInfo("HGCalValidation") << nameDetector_ << " with "
128  << theRecHitContainers->size()
129  << " element(s)";
130  for (HGCRecHitCollection::const_iterator it=theRecHitContainers->begin();
131  it !=theRecHitContainers->end(); ++it) {
132  ntot++; nused++;
133  DetId detId = it->id();
134  int layer = (detId.subdetId() == HGCEE) ? (HGCEEDetId(detId).layer()) : (HGCHEDetId(detId).layer());
135  recHitValidation(detId, layer, geom0, it);
136  }
137  } else {
138  ok = false;
139  edm::LogWarning("HGCalValidation") << "HGCRecHitCollection Handle does not exist !!!";
140  }
141  }
142  if (ok) fillHitsInfo();
143  edm::LogWarning("HGCalValidation") << "Event " << iEvent.id().event()
144  << " with " << ntot << " total and "
145  << nused << " used recHits";
146 }
147 
148 template<class T1, class T2>
150  const T1* geom, T2 it) {
151 
152  GlobalPoint global = geom->getPosition(detId);
153  double energy = it->energy();
154 
155  float globalx = global.x();
156  float globaly = global.y();
157  float globalz = global.z();
158 
159  HitsInfo hinfo;
160  hinfo.energy = energy;
161  hinfo.x = globalx;
162  hinfo.y = globaly;
163  hinfo.z = globalz;
164  hinfo.layer = layer;
165  hinfo.phi = global.phi();
166  hinfo.eta = global.eta();
167 
168  if (verbosity_>1)
169  edm::LogInfo("HGCalValidation") << " -------------------------- gx = "
170  << globalx << " gy = " << globaly
171  << " gz = " << globalz << " phi = "
172  << hinfo.phi << " eta = " << hinfo.eta;
173 
174  fillHitsInfo(hinfo);
175 
176  if (hinfo.eta > 0) fillOccupancyMap(OccupancyMap_plus, layer -1);
177  else fillOccupancyMap(OccupancyMap_minus, layer -1);
178 
179 }
180 
181 void HGCalRecHitValidation::fillOccupancyMap(std::map<int, int>& OccupancyMap, int layer){
182  if (OccupancyMap.find(layer) != OccupancyMap.end()) OccupancyMap[layer]++;
183  else OccupancyMap[layer] = 1;
184 }
185 
187 
188  for (auto itr = OccupancyMap_plus.begin() ; itr != OccupancyMap_plus.end(); ++itr) {
189  int layer = (*itr).first;
190  int occupancy = (*itr).second;
191  HitOccupancy_Plus_.at(layer)->Fill(occupancy);
192  }
193 
194  for (auto itr = OccupancyMap_minus.begin() ; itr != OccupancyMap_minus.end(); ++itr) {
195  int layer = (*itr).first;
196  int occupancy = (*itr).second;
197  HitOccupancy_Minus_.at(layer)->Fill(occupancy);
198  }
199 
200 }
201 
203 
204  unsigned int ilayer = hits.layer -1;
205  energy_.at(ilayer)->Fill(hits.energy);
206 
207  EtaPhi_Plus_.at(ilayer)->Fill(hits.eta , hits.phi);
208  EtaPhi_Minus_.at(ilayer)->Fill(hits.eta, hits.phi);
209 
210 }
211 
213  const edm::EventSetup& iSetup) {
214 
215  if (nameDetector_ == "HCal") {
217  iSetup.get<HcalRecNumberingRecord>().get( pHRNDC );
218  const HcalDDDRecConstants *hcons = &(*pHRNDC);
219  layers_ = hcons->getMaxDepth(1);
220  } else {
222  iSetup.get<IdealGeometryRecord>().get(nameDetector_, pHGDC);
223  const HGCalDDDConstants & hgcons_ = (*pHGDC);
224  layers_ = hgcons_.layers(true);
225  }
226 }
227 
229  edm::Run const&,
230  edm::EventSetup const&) {
231 
232  iB.setCurrentFolder("HGCalRecHitsV/"+nameDetector_);
233  std::ostringstream histoname;
234  for (unsigned int ilayer = 0; ilayer < layers_; ilayer++ ) {
235  histoname.str(""); histoname << "HitOccupancy_Plus_layer_" << ilayer;
236  HitOccupancy_Plus_.push_back(iB.book1D( histoname.str().c_str(), "RecHitOccupancy_Plus", 100, 0, 10000));
237  histoname.str(""); histoname << "HitOccupancy_Minus_layer_" << ilayer;
238  HitOccupancy_Minus_.push_back(iB.book1D( histoname.str().c_str(), "RecHitOccupancy_Minus", 100, 0, 10000));
239 
240  histoname.str(""); histoname << "EtaPhi_Plus_" << "layer_" << ilayer;
241  EtaPhi_Plus_.push_back(iB.book2D(histoname.str().c_str(), "Occupancy", 31, 1.45, 3.0, 72, -CLHEP::pi, CLHEP::pi));
242  histoname.str(""); histoname << "EtaPhi_Minus_" << "layer_" << ilayer;
243  EtaPhi_Minus_.push_back(iB.book2D(histoname.str().c_str(), "Occupancy", 31, -3.0, -1.45, 72, -CLHEP::pi, CLHEP::pi));
244 
245  histoname.str(""); histoname << "energy_layer_" << ilayer;
246  energy_.push_back(iB.book1D(histoname.str().c_str(),"energy_",100,0,0.002));
247  }//loop over layers ends here
248 
249  histoname.str(""); histoname << "SUMOfRecHitOccupancy_Plus";
250  MeanHitOccupancy_Plus_= iB.book1D( histoname.str().c_str(), "SUMOfRecHitOccupancy_Plus", layers_, -0.5, layers_-0.5);
251  histoname.str(""); histoname << "SUMOfRecHitOccupancy_Minus";
252  MeanHitOccupancy_Minus_ = iB.book1D( histoname.str().c_str(), "SUMOfRecHitOccupancy_Minus", layers_, -0.5,layers_-0.5);
253 }
254 
255 
256 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
259  desc.setUnknown();
260  descriptions.addDefault(desc);
261 }
262 
264 
265 //define this as a plug-in
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:41
HGCalRecHitValidation(const edm::ParameterSet &)
std::vector< MonitorElement * > EtaPhi_Minus_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
std::vector< HBHERecHit >::const_iterator const_iterator
T y() const
Definition: PV3DBase.h:63
std::vector< MonitorElement * > EtaPhi_Plus_
int layer() const
get the layer #
Definition: HGCHEDetId.h:50
std::vector< MonitorElement * > HitOccupancy_Minus_
const Double_t pi
MonitorElement * MeanHitOccupancy_Plus_
int getMaxDepth(const int type) const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int depth() const
get the tower depth
Definition: HcalDetId.cc:108
int iEvent
Definition: GenABIO.cc:230
std::map< int, int > OccupancyMap_minus
void addDefault(ParameterSetDescription const &psetDescription)
unsigned int layers(bool reco) const
T z() const
Definition: PV3DBase.h:64
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
MonitorElement * MeanHitOccupancy_Minus_
bool isValid() const
Definition: HandleBase.h:74
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
void dqmBeginRun(const edm::Run &, const edm::EventSetup &) override
const_iterator end() const
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: DetId.h:18
int layer() const
get the layer #
Definition: HGCEEDetId.h:49
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
void fillOccupancyMap(std::map< int, int > &OccupancyMap, int layer)
void analyze(const edm::Event &, const edm::EventSetup &) override
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
const T & get() const
Definition: EventSetup.h:56
T eta() const
Definition: PV3DBase.h:76
std::map< int, int > OccupancyMap_plus
edm::EventID id() const
Definition: EventBase.h:60
size_type size() const
std::vector< MonitorElement * > HitOccupancy_Plus_
bool isValid() const
Definition: ESHandle.h:47
T x() const
Definition: PV3DBase.h:62
T const * product() const
Definition: ESHandle.h:86
void recHitValidation(DetId &detId, int layer, const T1 *geom, T2 it)
std::vector< MonitorElement * > energy_
const_iterator begin() const
Definition: Run.h:42