CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HGCGeometryValidation.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <iostream>
3 
6 
11 
23 
29 
31 
32 #include "CLHEP/Geometry/Point3D.h"
33 #include "CLHEP/Geometry/Vector3D.h"
34 #include "CLHEP/Units/GlobalSystemOfUnits.h"
35 #include "CLHEP/Units/GlobalPhysicalConstants.h"
36 
37 const double mmtocm = 0.1;
38 
40 public:
41  explicit HGCGeometryValidation(const edm::ParameterSet &);
42  ~HGCGeometryValidation() override;
43  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
44 
45 protected:
46  void dqmBeginRun(edm::Run const &, edm::EventSetup const &) override;
47  void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
48  void analyze(const edm::Event &, const edm::EventSetup &) override;
49 
50 private:
52  std::vector<std::string> geometrySource_;
54  std::vector<edm::ESGetToken<HGCalDDDConstants, IdealGeometryRecord>> geomToken_;
55 
56  //HGCal geometry scheme
57  std::vector<const HGCalDDDConstants *> hgcGeometry_;
58 
59  //histogram related stuff
70 
75 };
76 
78  g4Token_ = consumes<PHGCalValidInfo>(cfg.getParameter<edm::InputTag>("g4Source"));
79  geometrySource_ = cfg.getUntrackedParameter<std::vector<std::string>>("geometrySource");
80  verbosity_ = cfg.getUntrackedParameter<int>("verbosity", 0);
81  for (const auto &name : geometrySource_)
82  geomToken_.emplace_back(
83  esConsumes<HGCalDDDConstants, IdealGeometryRecord, edm::Transition::BeginRun>(edm::ESInputTag{"", name}));
84 }
85 
87 
90  desc.setUnknown();
91  descriptions.addDefault(desc);
92 }
93 
95  //initiating hgcnumbering
96  for (size_t i = 0; i < geometrySource_.size(); i++) {
97  hgcGeometry_.emplace_back(&iSetup.getData(geomToken_[i]));
98  }
99 }
100 
102  iB.setCurrentFolder("HGCAL/HGCalSimHitsV/Geometry");
103 
104  //initiating histograms
105  heeTotEdepStep = iB.book1D("heeTotEdepStep", "", 100, 0, 100);
106  hefTotEdepStep = iB.book1D("hefTotEdepStep", "", 100, 0, 100);
107  hebTotEdepStep = iB.book1D("hebTotEdepStep", "", 100, 0, 100);
108 
109  hebLayerVsEnStep = iB.book2D("hebLayerVsEnStep", "", 25, 0, 25, 100, 0, 0.01);
110  hefLayerVsEnStep = iB.book2D("hefLayerVsEnStep", "", 36, 0, 36, 100, 0, 0.01);
111  heeLayerVsEnStep = iB.book2D("heeLayerVsEnStep", "", 84, 0, 84, 100, 0, 0.01);
112 
113  heeXG4VsId = iB.book2D("heeXG4VsId", "", 600, -300, 300, 600, -300, 300);
114  heeYG4VsId = iB.book2D("heeYG4VsId", "", 600, -300, 300, 600, -300, 300);
115  heeZG4VsId = iB.book2D("heeZG4VsId", "", 3000, 320, 350, 3000, 320, 350);
116 
117  hefXG4VsId = iB.book2D("hefXG4VsId", "", 600, -300, 300, 600, -300, 300);
118  hefYG4VsId = iB.book2D("hefYG4VsId", "", 600, -300, 300, 600, -300, 300);
119  hefZG4VsId = iB.book2D("hefZG4VsId", "", 6000, 350, 410, 6000, 350, 410);
120 
121  hebXG4VsId = iB.book2D("hebXG4VsId", "", 600, -300, 300, 600, -300, 300);
122  hebYG4VsId = iB.book2D("hebYG4VsId", "", 600, -300, 300, 600, -300, 300);
123  hebZG4VsId = iB.book2D("hebZG4VsId", "", 220, 400, 620, 220, 400, 620);
124 
125  heedzVsZ = iB.book2D("heedzVsZ", "", 600, 320, 350, 100, -1, 1);
126  heedyVsY = iB.book2D("heedyVsY", "", 400, -200, 200, 100, -1, 1);
127  heedxVsX = iB.book2D("heedxVsX", "", 400, -200, 200, 100, -1, 1);
128 
129  hefdzVsZ = iB.book2D("hefdzVsZ", "", 1200, 350, 410, 100, -1, 1);
130  hefdyVsY = iB.book2D("hefdyVsY", "", 400, -200, 200, 100, -1, 1);
131  hefdxVsX = iB.book2D("hefdxVsX", "", 400, -200, 200, 100, -1, 1);
132 
133  hebdzVsZ = iB.book2D("hebdzVsZ", "", 220, 400, 620, 100, -5, 5);
134  hebdyVsY = iB.book2D("hebdyVsY", "", 400, -200, 200, 100, -5, 5);
135  hebdxVsX = iB.book2D("hebdxVsX", "", 400, -200, 200, 100, -5, 5);
136 
137  heedzVsLayer = iB.book2D("heedzVsLayer", "", 100, 0, 100, 100, -1, 1);
138  hefdzVsLayer = iB.book2D("hefdzVsLayer", "", 40, 0, 40, 100, -1, 1);
139  hebdzVsLayer = iB.book2D("hebdzVsLayer", "", 50, 0, 25, 100, -5, 5);
140 
141  heedyVsLayer = iB.book2D("heedyVsLayer", "", 100, 0, 100, 100, -1, 1);
142  hefdyVsLayer = iB.book2D("hefdyVsLayer", "", 40, 0, 40, 100, -1, 1);
143  hebdyVsLayer = iB.book2D("hebdyVsLayer", "", 50, 0, 25, 100, -5, 5);
144 
145  heedxVsLayer = iB.book2D("heedxVsLayer", "", 100, 0, 100, 100, -1, 1);
146  hefdxVsLayer = iB.book2D("hefdxVsLayer", "", 40, 0, 40, 500, -1, 1);
147  hebdxVsLayer = iB.book2D("hebdxVsLayer", "", 50, 0, 25, 500, -5, 5.0);
148 
149  heedX = iB.book1D("heedX", "", 100, -1, 1);
150  heedY = iB.book1D("heedY", "", 100, -1, 1);
151  heedZ = iB.book1D("heedZ", "", 100, -1, 1);
152 
153  hefdX = iB.book1D("hefdX", "", 100, -1, 1);
154  hefdY = iB.book1D("hefdY", "", 100, -1, 1);
155  hefdZ = iB.book1D("hefdZ", "", 100, -1, 1);
156 
157  hebdX = iB.book1D("hebdX", "", 100, -1, 1);
158  hebdY = iB.book1D("hebdY", "", 100, -1, 1);
159  hebdZ = iB.book1D("hebdZ", "", 100, -1, 1);
160 }
161 
163  //Accessing G4 information
165  iEvent.getByToken(g4Token_, infoLayer);
166 
167  if (infoLayer.isValid()) {
168  //step vertex information
169  std::vector<float> hitVtxX = infoLayer->hitvtxX();
170  std::vector<float> hitVtxY = infoLayer->hitvtxY();
171  std::vector<float> hitVtxZ = infoLayer->hitvtxZ();
172  std::vector<unsigned int> hitDet = infoLayer->hitDets();
173  std::vector<unsigned int> hitIdx = infoLayer->hitIndex();
174 
175  //energy information
176  std::vector<float> edepLayerEE = infoLayer->eehgcEdep();
177  std::vector<float> edepLayerHE = infoLayer->hefhgcEdep();
178  std::vector<float> edepLayerHB = infoLayer->hebhgcEdep();
179 
180  unsigned int i;
181  for (i = 0; i < edepLayerEE.size(); i++) {
182  heeLayerVsEnStep->Fill(i, edepLayerEE.at(i));
183  }
184 
185  for (i = 0; i < edepLayerHE.size(); i++) {
186  hefLayerVsEnStep->Fill(i, edepLayerHE.at(i));
187  }
188 
189  for (i = 0; i < edepLayerHB.size(); i++) {
190  hebLayerVsEnStep->Fill(i, edepLayerHB.at(i));
191  }
192 
193  //fill total energy deposited
194  heeTotEdepStep->Fill((double)infoLayer->eeTotEdep());
195  hefTotEdepStep->Fill((double)infoLayer->hefTotEdep());
196  hebTotEdepStep->Fill((double)infoLayer->hebTotEdep());
197 
198  //loop over all hits
199  for (unsigned int i = 0; i < hitVtxX.size(); i++) {
200  hitVtxX.at(i) = mmtocm * hitVtxX.at(i);
201  hitVtxY.at(i) = mmtocm * hitVtxY.at(i);
202  hitVtxZ.at(i) = mmtocm * hitVtxZ.at(i);
203 
204  if ((hitDet.at(i) == (unsigned int)(DetId::Forward)) || (hitDet.at(i) == (unsigned int)(DetId::HGCalEE)) ||
205  (hitDet.at(i) == (unsigned int)(DetId::HGCalHSi)) || (hitDet.at(i) == (unsigned int)(DetId::HGCalHSc))) {
206  double xx, yy;
207  int dtype(0), layer(0), zside(1);
208  std::pair<float, float> xy;
209  if (hitDet.at(i) == (unsigned int)(DetId::Forward)) {
210  int subdet, wafer, celltype, cell;
211  HGCalTestNumbering::unpackHexagonIndex(hitIdx.at(i), subdet, zside, layer, wafer, celltype, cell);
212  dtype = (subdet == (int)(HGCEE)) ? 0 : 1;
213  xy = hgcGeometry_[dtype]->locateCell(cell, layer, wafer, true); //cm
214  } else if ((hitDet.at(i) == (unsigned int)(DetId::HGCalEE)) ||
215  (hitDet.at(i) == (unsigned int)(DetId::HGCalHSi))) {
216  HGCSiliconDetId id(hitIdx.at(i));
217  dtype = (id.det() == DetId::HGCalEE) ? 0 : 1;
218  layer = id.layer();
219  zside = id.zside();
220  xy = hgcGeometry_[dtype]->locateCell(layer, id.waferU(), id.waferV(), id.cellU(), id.cellV(), true, true);
221  } else {
222  HGCScintillatorDetId id(hitIdx.at(i));
223  dtype = 2;
224  layer = id.layer();
225  zside = id.zside();
226  xy = hgcGeometry_[dtype]->locateCellTrap(layer, id.ietaAbs(), id.iphi(), true);
227  }
228  double zz = hgcGeometry_[dtype]->waferZ(layer, true); //cm
229  if (zside < 0)
230  zz = -zz;
231  xx = (zside < 0) ? -xy.first : xy.first;
232  yy = xy.second;
233 
234  if (dtype == 0) {
235  heedzVsZ->Fill(zz, (hitVtxZ.at(i) - zz));
236  heedyVsY->Fill(yy, (hitVtxY.at(i) - yy));
237  heedxVsX->Fill(xx, (hitVtxX.at(i) - xx));
238 
239  heeXG4VsId->Fill(hitVtxX.at(i), xx);
240  heeYG4VsId->Fill(hitVtxY.at(i), yy);
241  heeZG4VsId->Fill(hitVtxZ.at(i), zz);
242 
243  heedzVsLayer->Fill(layer, (hitVtxZ.at(i) - zz));
244  heedyVsLayer->Fill(layer, (hitVtxY.at(i) - yy));
245  heedxVsLayer->Fill(layer, (hitVtxX.at(i) - xx));
246 
247  heedX->Fill((hitVtxX.at(i) - xx));
248  heedZ->Fill((hitVtxZ.at(i) - zz));
249  heedY->Fill((hitVtxY.at(i) - yy));
250 
251  } else if (dtype == 1) {
252  hefdzVsZ->Fill(zz, (hitVtxZ.at(i) - zz));
253  hefdyVsY->Fill(yy, (hitVtxY.at(i) - yy));
254  hefdxVsX->Fill(xx, (hitVtxX.at(i) - xx));
255 
256  hefXG4VsId->Fill(hitVtxX.at(i), xx);
257  hefYG4VsId->Fill(hitVtxY.at(i), yy);
258  hefZG4VsId->Fill(hitVtxZ.at(i), zz);
259 
260  hefdzVsLayer->Fill(layer, (hitVtxZ.at(i) - zz));
261  hefdyVsLayer->Fill(layer, (hitVtxY.at(i) - yy));
262  hefdxVsLayer->Fill(layer, (hitVtxX.at(i) - xx));
263 
264  hefdX->Fill((hitVtxX.at(i) - xx));
265  hefdZ->Fill((hitVtxZ.at(i) - zz));
266  hefdY->Fill((hitVtxY.at(i) - yy));
267 
268  } else {
269  hebdzVsZ->Fill(zz, (hitVtxZ.at(i) - zz));
270  hebdyVsY->Fill(yy, (hitVtxY.at(i) - yy));
271  hebdxVsX->Fill(xx, (hitVtxX.at(i) - xx));
272 
273  hebXG4VsId->Fill(hitVtxX.at(i), xx);
274  hebYG4VsId->Fill(hitVtxY.at(i), yy);
275  hebZG4VsId->Fill(hitVtxZ.at(i), zz);
276 
277  hebdzVsLayer->Fill(layer, (hitVtxZ.at(i) - zz));
278  hebdyVsLayer->Fill(layer, (hitVtxY.at(i) - yy));
279  hebdxVsLayer->Fill(layer, (hitVtxX.at(i) - xx));
280 
281  hebdX->Fill((hitVtxX.at(i) - xx));
282  hebdZ->Fill((hitVtxZ.at(i) - zz));
283  hebdY->Fill((hitVtxY.at(i) - yy));
284  }
285  }
286  } //end G4 hits
287 
288  } else {
289  if (verbosity_ > 0)
290  edm::LogVerbatim("HGCalValid") << "HGCGeometryValidation::No PHGCalInfo";
291  }
292 }
293 
294 //define this as a plug-in
Log< level::Info, true > LogVerbatim
T getUntrackedParameter(std::string const &, T const &) const
std::vector< edm::ESGetToken< HGCalDDDConstants, IdealGeometryRecord > > geomToken_
tuple cfg
Definition: looper.py:296
edm::EDGetTokenT< PHGCalValidInfo > g4Token_
uint16_t *__restrict__ id
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
int32_t waferU(const int32_t index)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
int zside(DetId const &)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
constexpr std::array< uint8_t, layerIndexSize > layer
void Fill(long long x)
bool getData(T &iHolder) const
Definition: EventSetup.h:122
int iEvent
Definition: GenABIO.cc:224
void addDefault(ParameterSetDescription const &psetDescription)
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
void dqmBeginRun(edm::Run const &, edm::EventSetup const &) override
MonitorElement * hefTotEdepStep
if(conf_.getParameter< bool >("UseStripCablingDB"))
MonitorElement * hefLayerVsEnStep
MonitorElement * heeLayerVsEnStep
const double mmtocm
bool isValid() const
Definition: HandleBase.h:70
void analyze(const edm::Event &, const edm::EventSetup &) override
MonitorElement * hebTotEdepStep
Basic2DVector< T > xy() const
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:177
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
MonitorElement * hebLayerVsEnStep
std::vector< const HGCalDDDConstants * > hgcGeometry_
int32_t waferV(const int32_t index)
HGCGeometryValidation(const edm::ParameterSet &)
std::vector< std::string > geometrySource_
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
MonitorElement * heeTotEdepStep
static void unpackHexagonIndex(const uint32_t &idx, int &subdet, int &z, int &lay, int &wafer, int &celltyp, int &cell)
Definition: Run.h:45