CMS 3D CMS Logo

HGCalRecHitValidation.cc
Go to the documentation of this file.
1 // system include files
2 #include <cmath>
3 #include <fstream>
4 #include <iostream>
5 #include <map>
6 #include <string>
7 #include <vector>
8 
15 
17 
20 
29 
37 
38 #include "CLHEP/Units/GlobalSystemOfUnits.h"
39 #include "TVector3.h"
40 
42 
43 public:
44  struct energysum{
46  double e15, e25, e50, e100, e250, e1000;
47  };
48 
49  struct HitsInfo{
51  x=y=z=time=energy=phi=eta=0.0;
52  layer=0;
53  }
54  float x, y, z, time, energy, phi, eta ;
55  float layer;
56  };
57 
58 
60  ~HGCalRecHitValidation() override {}
61 
62  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
63  void dqmBeginRun(const edm::Run&, const edm::EventSetup&) override;
64  void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
65  void analyze(const edm::Event&, const edm::EventSetup&) override;
66 
67 private:
68  template<class T1, class T2>
69  void recHitValidation(DetId & detId, int layer, const T1* geom, T2 it);
70  void fillHitsInfo();
71  void fillHitsInfo(HitsInfo& hits);
72  void fillOccupancyMap(std::map<int, int>& OccupancyMap, int layer);
73 
74  // ----------member data ---------------------------
77  bool ifHCAL_;
79  unsigned int layers_;
80  std::map<int, int> OccupancyMap_plus;
81  std::map<int, int> OccupancyMap_minus;
82 
83  std::vector<MonitorElement*> EtaPhi_Plus_;
84  std::vector<MonitorElement*> EtaPhi_Minus_;
85  std::vector<MonitorElement*> energy_;
86  std::vector<MonitorElement*> HitOccupancy_Plus_;
87  std::vector<MonitorElement*> HitOccupancy_Minus_;
90 };
91 
93  nameDetector_(iConfig.getParameter<std::string>("DetectorName")),
94  ifHCAL_(iConfig.getParameter<bool>("ifHCAL")),
95  verbosity_(iConfig.getUntrackedParameter<int>("Verbosity",0)) {
96 
97  auto temp = iConfig.getParameter<edm::InputTag>("RecHitSource");
98  if (nameDetector_ == "HGCalEESensitive" ||
99  nameDetector_ == "HGCalHESiliconSensitive" ||
100  nameDetector_ == "HGCalHEScintillatorSensitive") {
101  recHitSource_ = consumes<HGCRecHitCollection>(temp);
102  } else if (nameDetector_ == "HCal") {
103  if (ifHCAL_) recHitSource_ = consumes<HBHERecHitCollection>(temp);
104  else recHitSource_ = consumes<HGChebRecHitCollection>(temp);
105  } else {
106  throw cms::Exception("BadHGCRecHitSource")
107  << "HGCal DetectorName given as " << nameDetector_ << " must be: "
108  << "\"HGCalHESiliconSensitive\", \"HGCalHESiliconSensitive\", "
109  << "\"HGCalHEScintillatorSensitive\", or \"HCal\"!";
110  }
111 }
112 
113 
116  desc.add<std::string>("DetectorName","HGCalEESensitive");
117  desc.add<edm::InputTag>("RecHitSource",edm::InputTag("HGCalRecHit","HGCEERecHits"));
118  desc.add<bool>("ifHCAL",false);
119  desc.addUntracked<int>("Verbosity",0);
120  descriptions.add("hgcalRecHitValidationEE",desc);
121 }
122 
124  const edm::EventSetup& iSetup) {
125  OccupancyMap_plus.clear();
126  OccupancyMap_minus.clear();
127 
128  bool ok(true);
129  unsigned int ntot(0), nused(0);
130  if (nameDetector_ == "HCal") {
132  iSetup.get<CaloGeometryRecord>().get(geom);
133  if (!geom.isValid())
134  edm::LogWarning("HGCalValidation") << "Cannot get valid HGCalGeometry "
135  << "Object for " << nameDetector_;
136  const CaloGeometry* geom0 = geom.product();
137 
138  if (ifHCAL_) {
140  iEvent.getByToken(recHitSource_, hbhecoll);
141  if (hbhecoll.isValid()) {
142  if (verbosity_>0)
143  edm::LogVerbatim("HGCalValidation") << nameDetector_ << " with "
144  << hbhecoll->size()
145  << " element(s)";
146  for (const auto & it : *(hbhecoll.product())) {
147  DetId detId = it.id();
148  ntot++;
149  if (detId.subdetId() == HcalEndcap) {
150  nused++;
151  int layer = HcalDetId(detId).depth();
152  recHitValidation(detId, layer, geom0, &it);
153  }
154  }
155  } else {
156  ok = false;
157  edm::LogWarning("HGCalValidation") << "HBHERecHitCollection Handle does not exist !!!";
158  }
159  } else {
161  iEvent.getByToken(recHitSource_, hbhecoll);
162  if (hbhecoll.isValid()) {
163  if (verbosity_>0)
164  edm::LogVerbatim("HGCalValidation") << nameDetector_ << " with "
165  << hbhecoll->size()
166  << " element(s)";
167  for (const auto & it : *(hbhecoll.product())) {
168  DetId detId = it.id();
169  ntot++; nused++;
170  int layer = HcalDetId(detId).depth();
171  recHitValidation(detId, layer, geom0, &it);
172  }
173  } else {
174  ok = false;
175  edm::LogWarning("HGCalValidation") << "HGChebRecHitCollection Handle does not exist !!!";
176  }
177  }
178  } else {
180  iSetup.get<IdealGeometryRecord>().get(nameDetector_, geom);
181  if (!geom.isValid()) edm::LogWarning("HGCalValidation") << "Cannot get valid HGCalGeometry Object for " << nameDetector_;
182  const HGCalGeometry* geom0 = geom.product();
183 
184  edm::Handle<HGCRecHitCollection> theRecHitContainers;
185  iEvent.getByToken(recHitSource_, theRecHitContainers);
186  if (theRecHitContainers.isValid()) {
187  if (verbosity_>0)
188  edm::LogVerbatim("HGCalValidation") << nameDetector_ << " with "
189  << theRecHitContainers->size()
190  << " element(s)";
191  for (const auto & it : *(theRecHitContainers.product())) {
192  ntot++; nused++;
193  DetId detId = it.id();
194  int layer = HGCalDetId(detId).layer();
195  recHitValidation(detId, layer, geom0, &it);
196  }
197  } else {
198  ok = false;
199  edm::LogWarning("HGCalValidation") << "HGCRecHitCollection Handle does not exist !!!";
200  }
201  }
202  if (ok) fillHitsInfo();
203  edm::LogVerbatim("HGCalValidation") << "Event " << iEvent.id().event()
204  << " with " << ntot << " total and "
205  << nused << " used recHits";
206 }
207 
208 template<class T1, class T2>
210  const T1* geom, T2 it) {
211 
212  const GlobalPoint& global = geom->getPosition(detId);
213  double energy = it->energy();
214 
215  float globalx = global.x();
216  float globaly = global.y();
217  float globalz = global.z();
218 
219  HitsInfo hinfo;
220  hinfo.energy = energy;
221  hinfo.x = globalx;
222  hinfo.y = globaly;
223  hinfo.z = globalz;
224  hinfo.layer = layer;
225  hinfo.phi = global.phi();
226  hinfo.eta = global.eta();
227 
228  if (verbosity_>1)
229  edm::LogVerbatim("HGCalValidation") << "-------------------------- gx = "
230  << globalx << " gy = " << globaly
231  << " gz = " << globalz << " phi = "
232  << hinfo.phi << " eta = " << hinfo.eta;
233 
234  fillHitsInfo(hinfo);
235 
236  if (hinfo.eta > 0) fillOccupancyMap(OccupancyMap_plus, layer -1);
237  else fillOccupancyMap(OccupancyMap_minus, layer -1);
238 
239 }
240 
241 void HGCalRecHitValidation::fillOccupancyMap(std::map<int, int>& OccupancyMap, int layer){
242  if (OccupancyMap.find(layer) != OccupancyMap.end()) OccupancyMap[layer]++;
243  else OccupancyMap[layer] = 1;
244 }
245 
247 
248  for (auto itr = OccupancyMap_plus.begin() ; itr != OccupancyMap_plus.end(); ++itr) {
249  int layer = (*itr).first;
250  int occupancy = (*itr).second;
251  HitOccupancy_Plus_.at(layer)->Fill(occupancy);
252  }
253 
254  for (auto itr = OccupancyMap_minus.begin() ; itr != OccupancyMap_minus.end(); ++itr) {
255  int layer = (*itr).first;
256  int occupancy = (*itr).second;
257  HitOccupancy_Minus_.at(layer)->Fill(occupancy);
258  }
259 
260 }
261 
263 
264  unsigned int ilayer = hits.layer -1;
265  energy_.at(ilayer)->Fill(hits.energy);
266 
267  EtaPhi_Plus_.at(ilayer)->Fill(hits.eta , hits.phi);
268  EtaPhi_Minus_.at(ilayer)->Fill(hits.eta, hits.phi);
269 
270 }
271 
273  const edm::EventSetup& iSetup) {
274 
275  if (nameDetector_ == "HCal") {
277  iSetup.get<HcalRecNumberingRecord>().get( pHRNDC );
278  const HcalDDDRecConstants *hcons = &(*pHRNDC);
279  layers_ = hcons->getMaxDepth(1);
280  } else {
282  iSetup.get<IdealGeometryRecord>().get(nameDetector_, pHGDC);
283  const HGCalDDDConstants & hgcons_ = (*pHGDC);
284  layers_ = hgcons_.layers(true);
285  }
286 }
287 
289  edm::Run const&,
290  edm::EventSetup const&) {
291 
292  iB.setCurrentFolder("HGCalRecHitsV/"+nameDetector_);
293  std::ostringstream histoname;
294  for (unsigned int ilayer = 0; ilayer < layers_; ilayer++ ) {
295  histoname.str(""); histoname << "HitOccupancy_Plus_layer_" << ilayer;
296  HitOccupancy_Plus_.push_back(iB.book1D( histoname.str().c_str(), "RecHitOccupancy_Plus", 100, 0, 10000));
297  histoname.str(""); histoname << "HitOccupancy_Minus_layer_" << ilayer;
298  HitOccupancy_Minus_.push_back(iB.book1D( histoname.str().c_str(), "RecHitOccupancy_Minus", 100, 0, 10000));
299 
300  histoname.str(""); histoname << "EtaPhi_Plus_" << "layer_" << ilayer;
301  EtaPhi_Plus_.push_back(iB.book2D(histoname.str().c_str(), "Occupancy", 31, 1.45, 3.0, 72, -CLHEP::pi, CLHEP::pi));
302  histoname.str(""); histoname << "EtaPhi_Minus_" << "layer_" << ilayer;
303  EtaPhi_Minus_.push_back(iB.book2D(histoname.str().c_str(), "Occupancy", 31, -3.0, -1.45, 72, -CLHEP::pi, CLHEP::pi));
304 
305  histoname.str(""); histoname << "energy_layer_" << ilayer;
306  energy_.push_back(iB.book1D(histoname.str().c_str(),"energy_",100,0,0.002));
307  }//loop over layers ends here
308 
309  histoname.str(""); histoname << "SUMOfRecHitOccupancy_Plus";
310  MeanHitOccupancy_Plus_= iB.book1D( histoname.str().c_str(), "SUMOfRecHitOccupancy_Plus", layers_, -0.5, layers_-0.5);
311  histoname.str(""); histoname << "SUMOfRecHitOccupancy_Minus";
312  MeanHitOccupancy_Minus_ = iB.book1D( histoname.str().c_str(), "SUMOfRecHitOccupancy_Minus", layers_, -0.5,layers_-0.5);
313 }
314 
315 
316 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
317 
319 
320 //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_
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
std::vector< MonitorElement * > EtaPhi_Plus_
std::vector< MonitorElement * > HitOccupancy_Minus_
const Double_t pi
MonitorElement * MeanHitOccupancy_Plus_
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
unsigned int layers(bool reco) const
T z() const
Definition: PV3DBase.h:64
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:118
MonitorElement * MeanHitOccupancy_Minus_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
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:38
void dqmBeginRun(const edm::Run &, const edm::EventSetup &) override
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: DetId.h:18
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:279
void fillOccupancyMap(std::map< int, int > &OccupancyMap, int layer)
void analyze(const edm::Event &, const edm::EventSetup &) override
T const * product() const
Definition: Handle.h:81
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:136
int getMaxDepth(const int &type) const
const T & get() const
Definition: EventSetup.h:59
void add(std::string const &label, ParameterSetDescription const &psetDescription)
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_
int layer() const
get the layer #
Definition: HGCalDetId.h:48
Definition: Run.h:43