CMS 3D CMS Logo

HGCalSimHitValidation.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HGCalSimHitValidation
4 // Class: HGCalSimHitValidation
5 //
11 //
12 // Original Author: Raman Khurana
14 // Created: Fri, 31 Jan 2014 18:35:18 GMT
15 // $Id$
16 
17 // system include files
18 #include <cmath>
33 #include "CLHEP/Geometry/Point3D.h"
34 #include "CLHEP/Geometry/Vector3D.h"
35 #include "CLHEP/Units/GlobalSystemOfUnits.h"
36 #include "CLHEP/Units/GlobalPhysicalConstants.h"
38 
40  symmDet_(true) {
41 
42  nameDetector_ = iConfig.getParameter<std::string>("DetectorName");
43  caloHitSource_ = iConfig.getParameter<std::string>("CaloHitSource");
44  times_ = iConfig.getParameter<std::vector<double> >("TimeSlices");
45  verbosity_ = iConfig.getUntrackedParameter<int>("Verbosity",0);
46  testNumber_ = iConfig.getUntrackedParameter<bool>("TestNumber", true);
47  heRebuild_ = (nameDetector_ == "HCal") ? true : false;
48  tok_hepMC_ = consumes<edm::HepMCProduct>(edm::InputTag("generator"));
49  tok_hits_ = consumes<edm::PCaloHitContainer>(edm::InputTag("g4SimHits",caloHitSource_));
50  nTimes_ = (times_.size() > 6) ? 6 : times_.size();
51 }
52 
54 
56  const edm::EventSetup& iSetup) {
57 
58  //Generator input
60  iEvent.getByToken(tok_hepMC_,evtMC);
61  if (!evtMC.isValid()) {
62  edm::LogWarning("HGCalValidation") << "no HepMCProduct found\n";
63  } else {
64  const HepMC::GenEvent * myGenEvent = evtMC->GetEvent();
65  unsigned int k(0);
66  for (HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
67  p != myGenEvent->particles_end(); ++p, ++k) {
68  edm::LogInfo("HGCalValidation") << "Particle[" << k << "] with pt "
69  << (*p)->momentum().perp() << " eta "
70  << (*p)->momentum().eta() << " phi "
71  << (*p)->momentum().phi() << std::endl;
72  }
73  }
74 
75  //Now the hits
76  edm::Handle<edm::PCaloHitContainer> theCaloHitContainers;
77  iEvent.getByToken(tok_hits_, theCaloHitContainers);
78  if (theCaloHitContainers.isValid()) {
79  if (verbosity_>0)
80  edm::LogInfo("HGCalValidation") << " PcalohitItr = "
81  << theCaloHitContainers->size() << "\n";
82  std::vector<PCaloHit> caloHits;
83  caloHits.insert(caloHits.end(), theCaloHitContainers->begin(),
84  theCaloHitContainers->end());
85  if (heRebuild_ && testNumber_) {
86  for (unsigned int i=0; i<caloHits.size(); ++i) {
87  unsigned int id_ = caloHits[i].id();
89  if(hid.subdet()!=int(HcalEndcap)) hid = HcalDetId(HcalEmpty,hid.ieta(),hid.iphi(),hid.depth());
90  caloHits[i].setID(hid.rawId());
91  if (verbosity_>0)
92  edm::LogInfo("HGCalValidation") << "Hit[" << i << "] " << hid <<"\n";
93  }
94  }
95  analyzeHits(caloHits);
96  } else if (verbosity_>0) {
97  edm::LogInfo("HGCalValidation") << "PCaloHitContainer does not exist !!!\n";
98  }
99 }
100 
101 void HGCalSimHitValidation::analyzeHits (std::vector<PCaloHit>& hits) {
102 
103  std::map<int, int> OccupancyMap_plus, OccupancyMap_minus;
104  OccupancyMap_plus.clear(); OccupancyMap_minus.clear();
105 
106  std::map<uint32_t,std::pair<hitsinfo,energysum> > map_hits;
107  map_hits.clear();
108 
109  if (verbosity_ > 0)
110  edm::LogInfo("HGCalValidation") << nameDetector_ << " with " << hits.size()
111  << " PcaloHit elements\n";
112  unsigned int nused(0);
113  for (unsigned int i=0; i<hits.size(); i++) {
114  double energy = hits[i].energy();
115  double time = hits[i].time();
116  uint32_t id_ = hits[i].id();
117  int cell, sector, subsector, layer, zside;
118  int subdet(0);
119  if (heRebuild_) {
120  HcalDetId detId = HcalDetId(id_);
121  subdet = detId.subdet();
122  if (subdet != static_cast<int>(HcalEndcap)) continue;
123  cell = detId.ietaAbs();
124  sector = detId.iphi();
125  subsector = 1;
126  layer = detId.depth();
127  zside = detId.zside();
128  } else {
130  HGCalTestNumbering::unpackSquareIndex(id_, zside, layer, sector, subsector, cell);
131  } else {
132  HGCalTestNumbering::unpackHexagonIndex(id_, subdet, zside, layer, sector, subsector, cell);
133  }
134  }
135  nused++;
136  if (verbosity_>1)
137  edm::LogInfo("HGCalValidation") << "Detector " << nameDetector_
138  << " zside = " << zside
139  << " sector|wafer = " << sector
140  << " subsector|type = " << subsector
141  << " layer = " << layer
142  << " cell = " << cell
143  << " energy = " << energy
144  << " energyem = " << hits[i].energyEM()
145  << " energyhad = " << hits[i].energyHad()
146  << " time = " << time << "\n";
147 
148  HepGeom::Point3D<float> gcoord;
149  if (heRebuild_) {
150  std::pair<double,double> etaphi = hcons_->getEtaPhi(subdet,zside*cell,sector);
151  double rz = hcons_->getRZ(subdet,zside*cell,layer);
152 // std::cout << "i/p " << subdet << ":" << zside << ":" << cell << ":" << sector << ":" << layer << " o/p " << etaphi.first << ":" << etaphi.second << ":" << rz << std::endl;
153  gcoord = HepGeom::Point3D<float>(rz*cos(etaphi.second)/cosh(etaphi.first),
154  rz*sin(etaphi.second)/cosh(etaphi.first),
155  rz*tanh(etaphi.first));
156  } else {
158  std::pair<float,float> xy = hgcons_->locateCell(cell,layer,subsector,false);
159  const HepGeom::Point3D<float> lcoord(xy.first,xy.second,0);
160  int subs = (symmDet_ ? 0 : subsector);
161  id_ = HGCalTestNumbering::packSquareIndex(zside,layer,sector,subs,0);
162  gcoord = (transMap_[id_]*lcoord);
163  } else {
164  std::pair<float,float> xy = hgcons_->locateCell(cell,layer,sector,false);
165  double zp = hgcons_->waferZ(layer,false);
166  if (zside < 0) zp = -zp;
167  float xp = (zp < 0) ? -xy.first : xy.first;
168  gcoord = HepGeom::Point3D<float>(xp,xy.second,zp);
169  }
170  }
171  double tof = (gcoord.mag()*CLHEP::mm)/CLHEP::c_light;
172  if (verbosity_>1)
173  edm::LogInfo("HGCalValidation") << std::hex << id_ << std::dec
174  << " global coordinate " << gcoord
175  << " time " << time << ":" << tof <<"\n";
176  time -= tof;
177 
178  energysum esum;
179  hitsinfo hinfo;
180  if (map_hits.count(id_) != 0) {
181  hinfo = map_hits[id_].first;
182  esum = map_hits[id_].second;
183  } else {
184  hinfo.x = gcoord.x();
185  hinfo.y = gcoord.y();
186  hinfo.z = gcoord.z();
187  hinfo.sector = sector;
188  hinfo.cell = cell;
189  hinfo.layer = layer;
190  hinfo.phi = gcoord.getPhi();
191  hinfo.eta = gcoord.getEta();
192  }
193  esum.etotal += energy;
194  for (unsigned int k=0; k<nTimes_; ++k) {
195  if (time > 0 && time < times_[k]) esum.eTime[k] += energy;
196  }
197  if (verbosity_>1)
198  edm::LogInfo("HGCalValidation") << " -------------------------- gx = "
199  << hinfo.x << " gy = " << hinfo.y
200  << " gz = " << hinfo.z << " phi = "
201  << hinfo.phi << " eta = " << hinfo.eta
202  << std::endl;
203  map_hits[id_] = std::pair<hitsinfo,energysum>(hinfo,esum);
204  }
205  if (verbosity_>0)
206  edm::LogInfo("HGCalValidation") << nameDetector_ << " with "
207  << map_hits.size()
208  << " detector elements being hit\n";
209 
210  std::map<uint32_t,std::pair<hitsinfo,energysum> >::iterator itr;
211  for (itr = map_hits.begin() ; itr != map_hits.end(); ++itr) {
212  hitsinfo hinfo = (*itr).second.first;
213  energysum esum = (*itr).second.second;
214  int layer = hinfo.layer;
215 
216  for (unsigned int itimeslice = 0; itimeslice < nTimes_; itimeslice++ ) {
217  fillHitsInfo((*itr).second, itimeslice, esum.eTime[itimeslice]);
218  }
219 
220  double eta = hinfo.eta;
221 
222  if (eta > 0.0) fillOccupancyMap(OccupancyMap_plus, layer-1);
223  else fillOccupancyMap(OccupancyMap_minus, layer-1);
224  }
225  edm::LogInfo("HGCalValidation") << "With map:used:total " << hits.size()
226  << "|" << nused << "|" << map_hits.size()
227  << " hits\n";
228 
229  for (auto itr = OccupancyMap_plus.begin() ; itr != OccupancyMap_plus.end(); ++itr) {
230  int layer = (*itr).first;
231  int occupancy = (*itr).second;
232  HitOccupancy_Plus_.at(layer)->Fill(occupancy);
233  }
234  for (auto itr = OccupancyMap_minus.begin() ; itr != OccupancyMap_minus.end(); ++itr) {
235  int layer = (*itr).first;
236  int occupancy = (*itr).second;
237  HitOccupancy_Minus_.at(layer)->Fill(occupancy);
238  }
239 }
240 
241 void HGCalSimHitValidation::fillOccupancyMap(std::map<int, int>& OccupancyMap, int layer){
242  if (OccupancyMap.find(layer) != OccupancyMap.end()) {
243  OccupancyMap[layer] ++;
244  } else {
245  OccupancyMap[layer] = 1;
246  }
247 }
248 
249 void HGCalSimHitValidation::fillHitsInfo(std::pair<hitsinfo,energysum> hits,
250  unsigned int itimeslice, double esum){
251 
252  unsigned int ilayer = hits.first.layer - 1;
253  if (ilayer < layers_) {
254  energy_[itimeslice].at(ilayer)->Fill(esum);
255  if (itimeslice==0) {
256  EtaPhi_Plus_.at(ilayer)->Fill(hits.first.eta , hits.first.phi);
257  EtaPhi_Minus_.at(ilayer)->Fill(hits.first.eta , hits.first.phi);
258  }
259  } else {
260  if (verbosity_>0)
261  edm::LogInfo("HGCalValidation") << "Problematic Hit for "
262  << nameDetector_ << " at sector "
263  << hits.first.sector << " layer "
264  << hits.first.layer << " cell "
265  << hits.first.cell << " energy "
266  << hits.second.etotal << std::endl;
267  }
268 }
269 
271  if (verbosity_>0)
272  edm::LogInfo("HGCalValidation") << "Initialize HGCalDDDConstants for "
273  << nameDetector_ << " : " << hgcons_ <<"\n";
274 
275  if (hgcons_->geomMode() == HGCalGeometryMode::Square) {
276  const DDCompactView & cview = *ddViewH;
277  std::string attribute = "Volume";
279 
280  DDSpecificsMatchesValueFilter filter{DDValue(attribute, value, 0)};
281  DDFilteredView fv(cview,filter);
282  bool dodet = fv.firstChild();
283 
284  while (dodet) {
285  const DDSolid & sol = fv.logicalPart().solid();
286  std::string name = sol.name();
287  int isd = (name.find(nameDetector_) == std::string::npos) ? -1 : 1;
288  if (isd > 0) {
289  std::vector<int> copy = fv.copyNumbers();
290  int nsiz = (int)(copy.size());
291  int lay = (nsiz > 0) ? copy[nsiz-1] : -1;
292  int sec = (nsiz > 1) ? copy[nsiz-2] : -1;
293  int zp = (nsiz > 3) ? copy[nsiz-4] : -1;
294  if (zp !=1 ) zp = -1;
295  const DDTrap & trp = static_cast<DDTrap>(sol);
296  int subs = (trp.alpha1()>0 ? 1 : 0);
297  symmDet_ = (trp.alpha1()==0 ? true : false);
298  uint32_t id = HGCalTestNumbering::packSquareIndex(zp,lay,sec,subs,0);
299  DD3Vector x, y, z;
300  fv.rotation().GetComponents( x, y, z ) ;
301  const CLHEP::HepRep3x3 rotation ( x.X(), y.X(), z.X(),
302  x.Y(), y.Y(), z.Y(),
303  x.Z(), y.Z(), z.Z() );
304  const CLHEP::HepRotation hr ( rotation );
305  const CLHEP::Hep3Vector h3v ( fv.translation().X(),
306  fv.translation().Y(),
307  fv.translation().Z() ) ;
308  const HepGeom::Transform3D ht3d (hr, h3v);
309  transMap_.insert(std::make_pair(id,ht3d));
310  if (verbosity_>2)
311  edm::LogInfo("HGCalValidation") << HGCalDetId(id)
312  << " Transform using " << h3v
313  << " and " << hr << std::endl;
314  }
315  dodet = fv.next();
316  }
317  if (verbosity_>0)
318  edm::LogInfo("HGCalValidation") << "Finds " << transMap_.size()
319  << " elements and SymmDet_ = "
320  << symmDet_ << std::endl;
321  }
322  return true;
323 }
324 
325 // ------------ method called when starting to processes a run ------------
327  const edm::EventSetup& iSetup) {
328  if (heRebuild_) {
330  iSetup.get<HcalRecNumberingRecord>().get( pHRNDC );
331  hcons_ = &(*pHRNDC);
332  layers_ = hcons_->getMaxDepth(1);
333  } else {
335  iSetup.get<IdealGeometryRecord>().get(nameDetector_, pHGDC);
336  hgcons_ = &(*pHGDC);
337  layers_ = hgcons_->layers(false);
339  iSetup.get<IdealGeometryRecord>().get( pDD );
340  defineGeometry(pDD);
341  }
342  if (verbosity_>0)
343  edm::LogInfo("HGCalValidation") << nameDetector_ << " defined with "
344  << layers_ << " Layers\n";
345 }
346 
348  edm::Run const&,
349  edm::EventSetup const&) {
350 
351  iB.setCurrentFolder("HGCalSimHitsV/"+nameDetector_);
352 
353  std::ostringstream histoname;
354  for (unsigned int ilayer = 0; ilayer < layers_; ilayer++ ) {
355  histoname.str(""); histoname << "HitOccupancy_Plus_layer_" << ilayer;
356  HitOccupancy_Plus_.push_back(iB.book1D(histoname.str().c_str(), "HitOccupancy_Plus", 501, -0.5, 500.5));
357  histoname.str(""); histoname << "HitOccupancy_Minus_layer_" << ilayer;
358  HitOccupancy_Minus_.push_back(iB.book1D(histoname.str().c_str(), "HitOccupancy_Minus", 501, -0.5, 500.5));
359 
360  histoname.str(""); histoname << "EtaPhi_Plus_" << "layer_" << ilayer;
361  EtaPhi_Plus_.push_back(iB.book2D(histoname.str().c_str(), "Occupancy", 31, 1.45, 3.0, 72, -CLHEP::pi, CLHEP::pi));
362  histoname.str(""); histoname << "EtaPhi_Minus_" << "layer_" << ilayer;
363  EtaPhi_Minus_.push_back(iB.book2D(histoname.str().c_str(), "Occupancy", 31, -3.0, -1.45, 72, -CLHEP::pi, CLHEP::pi));
364 
365  for (unsigned int itimeslice = 0; itimeslice < nTimes_ ; itimeslice++ ) {
366  histoname.str(""); histoname << "energy_time_"<< itimeslice << "_layer_" << ilayer;
367  energy_[itimeslice].push_back(iB.book1D(histoname.str().c_str(),"energy_",100,0,0.1));
368  }
369  }
370 
371  MeanHitOccupancy_Plus_ = iB.book1D("MeanHitOccupancy_Plus", "MeanHitOccupancy_Plus", layers_, 0.5, layers_ + 0.5);
372  MeanHitOccupancy_Minus_ = iB.book1D("MeanHitOccupancy_Minus", "MeanHitOccupancy_Minus", layers_, 0.5, layers_ + 0.5);
373 }
374 
375 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
377  //The following says we do not know what parameters are allowed so do no validation
378  // Please change this to state exactly what you do use, even if it is no parameters
380  desc.setUnknown();
381  descriptions.addDefault(desc);
382 }
383 
385 //define this as a plug-in
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
MonitorElement * MeanHitOccupancy_Plus_
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::vector< MonitorElement * > energy_[6]
const DDLogicalPart & logicalPart() const
The logical-part of the current node in the filtered-view.
const N & name() const
Definition: DDBase.h:78
void fillOccupancyMap(std::map< int, int > &OccupancyMap, int layer)
std::pair< double, double > getEtaPhi(const int &subdet, const int &ieta, const int &iphi) const
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:49
std::vector< MonitorElement * > HitOccupancy_Plus_
void analyze(const edm::Event &, const edm::EventSetup &) override
ProductID id() const
Definition: HandleBase.cc:15
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
nav_type copyNumbers() const
return the stack of copy numbers
int zside() const
get the z-side of the cell (1/-1)
Definition: HcalDetId.cc:114
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::EDGetTokenT< edm::HepMCProduct > tok_hepMC_
const DDRotationMatrix & rotation() const
The absolute rotation of the current node.
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
void analyzeHits(std::vector< PCaloHit > &hits)
std::vector< MonitorElement * > HitOccupancy_Minus_
const DDSolid & solid(void) const
Returns a reference object of the solid being the shape of this LogicalPart.
std::pair< float, float > locateCell(int cell, int lay, int type, bool reco) const
type of data representation of DDCompactView
Definition: DDCompactView.h:90
A DDSolid represents the shape of a part.
Definition: DDSolid.h:38
MonitorElement * MeanHitOccupancy_Minus_
const Double_t pi
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
bool defineGeometry(edm::ESTransientHandle< DDCompactView > &ddViewH)
int depth() const
get the tower depth
Definition: HcalDetId.cc:129
int iEvent
Definition: GenABIO.cc:230
void addDefault(ParameterSetDescription const &psetDescription)
unsigned int layers(bool reco) const
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > DD3Vector
A DD Translation is currently implemented with Root Vector3D.
Definition: DDTranslation.h:6
bool next()
set current node to the next node in the filtered tree
const HcalDDDRecConstants * hcons_
susybsm::HSCParticleRef hr
Definition: classes.h:26
double getRZ(const int &subdet, const int &ieta, const int &depth) const
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
void dqmBeginRun(const edm::Run &, const edm::EventSetup &) override
int ieta() const
get the cell ieta
Definition: HcalDetId.h:56
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
Interface to a Trapezoid.
Definition: DDSolid.h:79
std::pair< T, T > etaphi(T x, T y, T z)
Definition: FastMath.h:128
HGCalGeometryMode::GeometryMode geomMode() const
Definition: value.py:1
std::vector< MonitorElement * > EtaPhi_Plus_
bool isValid() const
Definition: HandleBase.h:74
HGCalSimHitValidation(const edm::ParameterSet &)
double waferZ(int layer, bool reco) const
int k[5][pyjets_maxn]
int ietaAbs() const
get the absolute value of the cell ieta
Definition: HcalDetId.cc:119
int iphi() const
get the cell iphi
Definition: HcalDetId.cc:124
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
int getMaxDepth(const int &type) const
static void unpackSquareIndex(const uint32_t &idx, int &z, int &lay, int &sec, int &subsec, int &cell)
const T & get() const
Definition: EventSetup.h:55
double alpha1(void) const
Angle with respect to the y axis from the centre of the side at y=-pDy1 to the centre at y=+pDy1 of t...
Definition: DDSolid.cc:183
edm::EDGetTokenT< edm::PCaloHitContainer > tok_hits_
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
std::map< uint32_t, HepGeom::Transform3D > transMap_
bool firstChild()
set the current node to the first child ...
const HGCalDDDConstants * hgcons_
std::vector< MonitorElement * > EtaPhi_Minus_
const DDTranslation & translation() const
The absolute translation of the current node.
DetId relabel(const uint32_t testId) const
static uint32_t packSquareIndex(int z, int lay, int sec, int subsec, int cell)
void fillHitsInfo(std::pair< hitsinfo, energysum > hit_, unsigned int itimeslice, double esum)
std::vector< double > times_
static void unpackHexagonIndex(const uint32_t &idx, int &subdet, int &z, int &lay, int &wafer, int &celltyp, int &cell)
Definition: Run.h:43