CMS 3D CMS Logo

HGCGeometryValidation.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <iostream>
3 
6 
10 
23 
27 
29 
30 const double mmtocm = 0.1;
31 
33 public:
34  explicit HGCGeometryValidation(const edm::ParameterSet &);
35  ~HGCGeometryValidation() override = default;
36  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
37 
38 protected:
39  void dqmBeginRun(edm::Run const &, edm::EventSetup const &) override;
40  void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
41  void analyze(const edm::Event &, const edm::EventSetup &) override;
42 
43 private:
45  const std::vector<std::string> geometrySource_;
46  const int verbosity_;
47  const std::vector<edm::ESGetToken<HGCalDDDConstants, IdealGeometryRecord>> geomToken_;
48 
49  //HGCal geometry scheme
50  std::vector<const HGCalDDDConstants *> hgcGeometry_;
51 
52  //histogram related stuff
63 
68 };
69 
71  : g4Token_(consumes<PHGCalValidInfo>(cfg.getParameter<edm::InputTag>("g4Source"))),
72  geometrySource_(cfg.getUntrackedParameter<std::vector<std::string>>("geometrySource")),
73  verbosity_(cfg.getUntrackedParameter<int>("verbosity", 0)),
74  geomToken_{edm::vector_transform(geometrySource_, [this](const std::string &name) {
75  return esConsumes<HGCalDDDConstants, IdealGeometryRecord, edm::Transition::BeginRun>(edm::ESInputTag{"", name});
76  })} {}
77 
80  desc.setUnknown();
81  descriptions.addDefault(desc);
82 }
83 
85  for (size_t i = 0; i < geometrySource_.size(); i++) {
86  hgcGeometry_.emplace_back(&iSetup.getData(geomToken_[i]));
87  }
88 }
89 
91  iB.setCurrentFolder("HGCAL/HGCalSimHitsV/Geometry");
92 
93  //initiating histograms
94  heeTotEdepStep = iB.book1D("heeTotEdepStep", "", 100, 0, 100);
95  hefTotEdepStep = iB.book1D("hefTotEdepStep", "", 100, 0, 100);
96  hebTotEdepStep = iB.book1D("hebTotEdepStep", "", 100, 0, 100);
97 
98  hebLayerVsEnStep = iB.book2D("hebLayerVsEnStep", "", 25, 0, 25, 100, 0, 0.01);
99  hefLayerVsEnStep = iB.book2D("hefLayerVsEnStep", "", 36, 0, 36, 100, 0, 0.01);
100  heeLayerVsEnStep = iB.book2D("heeLayerVsEnStep", "", 84, 0, 84, 100, 0, 0.01);
101 
102  heeXG4VsId = iB.book2D("heeXG4VsId", "", 600, -300, 300, 600, -300, 300);
103  heeYG4VsId = iB.book2D("heeYG4VsId", "", 600, -300, 300, 600, -300, 300);
104  heeZG4VsId = iB.book2D("heeZG4VsId", "", 3000, 320, 350, 3000, 320, 350);
105 
106  hefXG4VsId = iB.book2D("hefXG4VsId", "", 600, -300, 300, 600, -300, 300);
107  hefYG4VsId = iB.book2D("hefYG4VsId", "", 600, -300, 300, 600, -300, 300);
108  hefZG4VsId = iB.book2D("hefZG4VsId", "", 6000, 350, 410, 6000, 350, 410);
109 
110  hebXG4VsId = iB.book2D("hebXG4VsId", "", 600, -300, 300, 600, -300, 300);
111  hebYG4VsId = iB.book2D("hebYG4VsId", "", 600, -300, 300, 600, -300, 300);
112  hebZG4VsId = iB.book2D("hebZG4VsId", "", 220, 400, 620, 220, 400, 620);
113 
114  heedzVsZ = iB.book2D("heedzVsZ", "", 600, 320, 350, 100, -1, 1);
115  heedyVsY = iB.book2D("heedyVsY", "", 400, -200, 200, 100, -1, 1);
116  heedxVsX = iB.book2D("heedxVsX", "", 400, -200, 200, 100, -1, 1);
117 
118  hefdzVsZ = iB.book2D("hefdzVsZ", "", 1200, 350, 410, 100, -1, 1);
119  hefdyVsY = iB.book2D("hefdyVsY", "", 400, -200, 200, 100, -1, 1);
120  hefdxVsX = iB.book2D("hefdxVsX", "", 400, -200, 200, 100, -1, 1);
121 
122  hebdzVsZ = iB.book2D("hebdzVsZ", "", 220, 400, 620, 100, -5, 5);
123  hebdyVsY = iB.book2D("hebdyVsY", "", 400, -200, 200, 100, -5, 5);
124  hebdxVsX = iB.book2D("hebdxVsX", "", 400, -200, 200, 100, -5, 5);
125 
126  heedzVsLayer = iB.book2D("heedzVsLayer", "", 100, 0, 100, 100, -1, 1);
127  hefdzVsLayer = iB.book2D("hefdzVsLayer", "", 40, 0, 40, 100, -1, 1);
128  hebdzVsLayer = iB.book2D("hebdzVsLayer", "", 50, 0, 25, 100, -5, 5);
129 
130  heedyVsLayer = iB.book2D("heedyVsLayer", "", 100, 0, 100, 100, -1, 1);
131  hefdyVsLayer = iB.book2D("hefdyVsLayer", "", 40, 0, 40, 100, -1, 1);
132  hebdyVsLayer = iB.book2D("hebdyVsLayer", "", 50, 0, 25, 100, -5, 5);
133 
134  heedxVsLayer = iB.book2D("heedxVsLayer", "", 100, 0, 100, 100, -1, 1);
135  hefdxVsLayer = iB.book2D("hefdxVsLayer", "", 40, 0, 40, 500, -1, 1);
136  hebdxVsLayer = iB.book2D("hebdxVsLayer", "", 50, 0, 25, 500, -5, 5.0);
137 
138  heedX = iB.book1D("heedX", "", 100, -1, 1);
139  heedY = iB.book1D("heedY", "", 100, -1, 1);
140  heedZ = iB.book1D("heedZ", "", 100, -1, 1);
141 
142  hefdX = iB.book1D("hefdX", "", 100, -1, 1);
143  hefdY = iB.book1D("hefdY", "", 100, -1, 1);
144  hefdZ = iB.book1D("hefdZ", "", 100, -1, 1);
145 
146  hebdX = iB.book1D("hebdX", "", 100, -1, 1);
147  hebdY = iB.book1D("hebdY", "", 100, -1, 1);
148  hebdZ = iB.book1D("hebdZ", "", 100, -1, 1);
149 }
150 
152  //Accessing G4 information
153  const edm::Handle<PHGCalValidInfo> &infoLayer = iEvent.getHandle(g4Token_);
154 
155  if (infoLayer.isValid()) {
156  //step vertex information
157  std::vector<float> hitVtxX = infoLayer->hitvtxX();
158  std::vector<float> hitVtxY = infoLayer->hitvtxY();
159  std::vector<float> hitVtxZ = infoLayer->hitvtxZ();
160  const std::vector<unsigned int> &hitDet = infoLayer->hitDets();
161  const std::vector<unsigned int> &hitIdx = infoLayer->hitIndex();
162 
163  //energy information
164  const std::vector<float> &edepLayerEE = infoLayer->eehgcEdep();
165  const std::vector<float> &edepLayerHE = infoLayer->hefhgcEdep();
166  const std::vector<float> &edepLayerHB = infoLayer->hebhgcEdep();
167 
168  unsigned int i;
169  for (i = 0; i < edepLayerEE.size(); i++) {
170  heeLayerVsEnStep->Fill(i, edepLayerEE[i]);
171  }
172 
173  for (i = 0; i < edepLayerHE.size(); i++) {
174  hefLayerVsEnStep->Fill(i, edepLayerHE[i]);
175  }
176 
177  for (i = 0; i < edepLayerHB.size(); i++) {
178  hebLayerVsEnStep->Fill(i, edepLayerHB[i]);
179  }
180 
181  //fill total energy deposited
182  heeTotEdepStep->Fill((double)infoLayer->eeTotEdep());
183  hefTotEdepStep->Fill((double)infoLayer->hefTotEdep());
184  hebTotEdepStep->Fill((double)infoLayer->hebTotEdep());
185 
186  //loop over all hits
187  for (unsigned int i = 0; i < hitVtxX.size(); i++) {
188  hitVtxX[i] *= mmtocm;
189  hitVtxY[i] *= mmtocm;
190  hitVtxZ[i] *= mmtocm;
191 
192  double xx, yy;
193  int dtype(0), layer(0), zside(1);
194  std::pair<float, float> xy;
195  if ((hitDet[i] == static_cast<unsigned int>(DetId::HGCalEE)) ||
196  (hitDet[i] == static_cast<unsigned int>(DetId::HGCalHSi))) {
197  HGCSiliconDetId id(hitIdx[i]);
198  dtype = (id.det() == DetId::HGCalEE) ? 0 : 1;
199  layer = id.layer();
200  zside = id.zside();
201  xy = hgcGeometry_[dtype]->locateCell(
202  zside, layer, id.waferU(), id.waferV(), id.cellU(), id.cellV(), true, true, false, false);
203  } else {
204  HGCScintillatorDetId id(hitIdx[i]);
205  dtype = 2;
206  layer = id.layer();
207  zside = id.zside();
208  xy = hgcGeometry_[dtype]->locateCellTrap(zside, layer, id.ietaAbs(), id.iphi(), true, false);
209  }
210  double zz = hgcGeometry_[dtype]->waferZ(layer, true); //cm
211  if (zside < 0)
212  zz = -zz;
213  xx = (zside < 0) ? -xy.first : xy.first;
214  yy = xy.second;
215 
216  if (dtype == 0) {
217  heedzVsZ->Fill(zz, (hitVtxZ[i] - zz));
218  heedyVsY->Fill(yy, (hitVtxY[i] - yy));
219  heedxVsX->Fill(xx, (hitVtxX[i] - xx));
220 
221  heeXG4VsId->Fill(hitVtxX[i], xx);
222  heeYG4VsId->Fill(hitVtxY[i], yy);
223  heeZG4VsId->Fill(hitVtxZ[i], zz);
224 
225  heedzVsLayer->Fill(layer, (hitVtxZ[i] - zz));
226  heedyVsLayer->Fill(layer, (hitVtxY[i] - yy));
227  heedxVsLayer->Fill(layer, (hitVtxX[i] - xx));
228 
229  heedX->Fill((hitVtxX[i] - xx));
230  heedZ->Fill((hitVtxZ[i] - zz));
231  heedY->Fill((hitVtxY[i] - yy));
232 
233  } else if (dtype == 1) {
234  hefdzVsZ->Fill(zz, (hitVtxZ[i] - zz));
235  hefdyVsY->Fill(yy, (hitVtxY[i] - yy));
236  hefdxVsX->Fill(xx, (hitVtxX[i] - xx));
237 
238  hefXG4VsId->Fill(hitVtxX[i], xx);
239  hefYG4VsId->Fill(hitVtxY[i], yy);
240  hefZG4VsId->Fill(hitVtxZ[i], zz);
241 
242  hefdzVsLayer->Fill(layer, (hitVtxZ[i] - zz));
243  hefdyVsLayer->Fill(layer, (hitVtxY[i] - yy));
244  hefdxVsLayer->Fill(layer, (hitVtxX[i] - xx));
245 
246  hefdX->Fill((hitVtxX[i] - xx));
247  hefdZ->Fill((hitVtxZ[i] - zz));
248  hefdY->Fill((hitVtxY[i] - yy));
249 
250  } else {
251  hebdzVsZ->Fill(zz, (hitVtxZ[i] - zz));
252  hebdyVsY->Fill(yy, (hitVtxY[i] - yy));
253  hebdxVsX->Fill(xx, (hitVtxX[i] - xx));
254 
255  hebXG4VsId->Fill(hitVtxX[i], xx);
256  hebYG4VsId->Fill(hitVtxY[i], yy);
257  hebZG4VsId->Fill(hitVtxZ[i], zz);
258 
259  hebdzVsLayer->Fill(layer, (hitVtxZ[i] - zz));
260  hebdyVsLayer->Fill(layer, (hitVtxY[i] - yy));
261  hebdxVsLayer->Fill(layer, (hitVtxX[i] - xx));
262 
263  hebdX->Fill((hitVtxX[i] - xx));
264  hebdZ->Fill((hitVtxZ[i] - zz));
265  hebdY->Fill((hitVtxY[i] - yy));
266  }
267  } //end G4 hits
268 
269  } else {
270  if (verbosity_ > 0)
271  edm::LogVerbatim("HGCalValid") << "HGCGeometryValidation::No PHGCalInfo";
272  }
273 }
274 
275 //define this as a plug-in
Log< level::Info, true > LogVerbatim
std::vector< unsigned int > hitIndex() const
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
float hefTotEdep() const
std::vector< float > hefhgcEdep() const
std::vector< float > hitvtxZ() const
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
int32_t waferU(const int32_t index)
const edm::EDGetTokenT< PHGCalValidInfo > g4Token_
int zside(DetId const &)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
std::vector< unsigned int > hitDets() const
void Fill(long long x)
ALPAKA_FN_ACC ALPAKA_FN_INLINE uint32_t ietaAbs(uint32_t id)
const std::vector< std::string > geometrySource_
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
MonitorElement * hefLayerVsEnStep
MonitorElement * heeLayerVsEnStep
std::vector< float > eehgcEdep() const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const double mmtocm
~HGCGeometryValidation() override=default
void analyze(const edm::Event &, const edm::EventSetup &) override
std::vector< float > hebhgcEdep() const
MonitorElement * hebTotEdepStep
std::vector< float > hitvtxY() 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:212
bool isValid() const
Definition: HandleBase.h:70
MonitorElement * hebLayerVsEnStep
const std::vector< edm::ESGetToken< HGCalDDDConstants, IdealGeometryRecord > > geomToken_
HLT enums.
std::vector< const HGCalDDDConstants * > hgcGeometry_
int32_t waferV(const int32_t index)
float hebTotEdep() const
HGCGeometryValidation(const edm::ParameterSet &)
float eeTotEdep() const
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
if(threadIdxLocalY==0 &&threadIdxLocalX==0)
Definition: Run.h:45
std::vector< float > hitvtxX() const