CMS 3D CMS Logo

SiPixelCompareRecHits.cc
Go to the documentation of this file.
1 // TODO: change file name to SiPixelCompareRecHitsSoA.cc when CUDA code is removed
2 
21 
22 // TODO: change class name to SiPixelCompareRecHitsSoA when CUDA code is removed
23 template <typename T>
25 public:
27 
29  ~SiPixelCompareRecHits() override = default;
30  void dqmBeginRun(const edm::Run&, const edm::EventSetup&) override;
31  void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& iRun, edm::EventSetup const& iSetup) override;
32  void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
33  // analyzeSeparate is templated to accept distinct types of SoAs
34  // The default use case is to use rechits from Alpaka reconstructed on CPU and GPU;
35  template <typename U, typename V>
36  void analyzeSeparate(U tokenRef, V tokenTar, const edm::Event& iEvent);
37  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
38 
39 private:
42  // these two are both on Host but originally they have been produced on Host or on Device
46  const float mind2cut_;
48  static constexpr float micron_ = 10000.;
49  const TrackerGeometry* tkGeom_ = nullptr;
50  const TrackerTopology* tTopo_ = nullptr;
52  MonitorElement* hBchargeL_[4]; // max 4 barrel hits
57  MonitorElement* hFchargeD_[2][12]; // max 12 endcap disks
62  //differences
73 };
74 
75 //
76 // constructors
77 //
78 template <typename T>
81  topoToken_(esConsumes<TrackerTopology, TrackerTopologyRcd, edm::Transition::BeginRun>()),
82  tokenSoAHitsReference_(consumes(iConfig.getParameter<edm::InputTag>("pixelHitsReferenceSoA"))),
83  tokenSoAHitsTarget_(consumes(iConfig.getParameter<edm::InputTag>("pixelHitsTargetSoA"))),
84  topFolderName_(iConfig.getParameter<std::string>("topFolderName")),
85  mind2cut_(iConfig.getParameter<double>("minD2cut")) {}
86 
87 //
88 // Begin Run
89 //
90 template <typename T>
92  tkGeom_ = &iSetup.getData(geomToken_);
93  tTopo_ = &iSetup.getData(topoToken_);
94 }
95 
96 template <typename T>
97 template <typename U, typename V>
99  const auto& rhsoaHandleRef = iEvent.getHandle(tokenRef);
100  const auto& rhsoaHandleTar = iEvent.getHandle(tokenTar);
101 
102  if (not rhsoaHandleRef or not rhsoaHandleTar) {
103  edm::LogWarning out("SiPixelCompareRecHits");
104  if (not rhsoaHandleRef) {
105  out << "reference rechits not found; ";
106  }
107  if (not rhsoaHandleTar) {
108  out << "target rechits not found; ";
109  }
110  out << "the comparison will not run.";
111  return;
112  }
113 
114  auto const& rhsoaRef = *rhsoaHandleRef;
115  auto const& rhsoaTar = *rhsoaHandleTar;
116 
117  auto const& soa2dRef = rhsoaRef.const_view();
118  auto const& soa2dTar = rhsoaTar.const_view();
119 
120  uint32_t nHitsRef = soa2dRef.metadata().size();
121  uint32_t nHitsTar = soa2dTar.metadata().size();
122 
123  hnHits_->Fill(nHitsRef, nHitsTar);
124  auto detIds = tkGeom_->detUnitIds();
125  for (uint32_t i = 0; i < nHitsRef; i++) {
126  float minD = mind2cut_;
127  uint32_t matchedHit = invalidHit_;
128  uint16_t indRef = soa2dRef[i].detectorIndex();
129  float xLocalRef = soa2dRef[i].xLocal();
130  float yLocalRef = soa2dRef[i].yLocal();
131  for (uint32_t j = 0; j < nHitsTar; j++) {
132  if (soa2dTar.detectorIndex(j) == indRef) {
133  float dx = xLocalRef - soa2dTar[j].xLocal();
134  float dy = yLocalRef - soa2dTar[j].yLocal();
135  float distance = dx * dx + dy * dy;
136  if (distance < minD) {
137  minD = distance;
138  matchedHit = j;
139  }
140  }
141  }
142  DetId id = detIds[indRef];
143  uint32_t chargeRef = soa2dRef[i].chargeAndStatus().charge;
144  int16_t sizeXRef = std::ceil(float(std::abs(soa2dRef[i].clusterSizeX()) / 8.));
145  int16_t sizeYRef = std::ceil(float(std::abs(soa2dRef[i].clusterSizeY()) / 8.));
146  uint32_t chargeTar = 0;
147  int16_t sizeXTar = -99;
148  int16_t sizeYTar = -99;
149  float xLocalTar = -999.;
150  float yLocalTar = -999.;
151  if (matchedHit != invalidHit_) {
152  chargeTar = soa2dTar[matchedHit].chargeAndStatus().charge;
153  sizeXTar = std::ceil(float(std::abs(soa2dTar[matchedHit].clusterSizeX()) / 8.));
154  sizeYTar = std::ceil(float(std::abs(soa2dTar[matchedHit].clusterSizeY()) / 8.));
155  xLocalTar = soa2dTar[matchedHit].xLocal();
156  yLocalTar = soa2dTar[matchedHit].yLocal();
157  }
158  switch (id.subdetId()) {
160  hBchargeL_[tTopo_->pxbLayer(id) - 1]->Fill(chargeRef, chargeTar);
161  hBsizexL_[tTopo_->pxbLayer(id) - 1]->Fill(sizeXRef, sizeXTar);
162  hBsizeyL_[tTopo_->pxbLayer(id) - 1]->Fill(sizeYRef, sizeYTar);
163  hBposxL_[tTopo_->pxbLayer(id) - 1]->Fill(xLocalRef, xLocalTar);
164  hBposyL_[tTopo_->pxbLayer(id) - 1]->Fill(yLocalRef, yLocalTar);
165  hBchargeDiff_->Fill(chargeRef - chargeTar);
166  hBsizeXDiff_->Fill(sizeXRef - sizeXTar);
167  hBsizeYDiff_->Fill(sizeYRef - sizeYTar);
168  hBposXDiff_->Fill(micron_ * (xLocalRef - xLocalTar));
169  hBposYDiff_->Fill(micron_ * (yLocalRef - yLocalTar));
170  break;
172  hFchargeD_[tTopo_->pxfSide(id) - 1][tTopo_->pxfDisk(id) - 1]->Fill(chargeRef, chargeTar);
173  hFsizexD_[tTopo_->pxfSide(id) - 1][tTopo_->pxfDisk(id) - 1]->Fill(sizeXRef, sizeXTar);
174  hFsizeyD_[tTopo_->pxfSide(id) - 1][tTopo_->pxfDisk(id) - 1]->Fill(sizeYRef, sizeYTar);
175  hFposxD_[tTopo_->pxfSide(id) - 1][tTopo_->pxfDisk(id) - 1]->Fill(xLocalRef, xLocalTar);
176  hFposyD_[tTopo_->pxfSide(id) - 1][tTopo_->pxfDisk(id) - 1]->Fill(yLocalRef, yLocalTar);
177  hFchargeDiff_->Fill(chargeRef - chargeTar);
178  hFsizeXDiff_->Fill(sizeXRef - sizeXTar);
179  hFsizeYDiff_->Fill(sizeYRef - sizeYTar);
180  hFposXDiff_->Fill(micron_ * (xLocalRef - xLocalTar));
181  hFposYDiff_->Fill(micron_ * (yLocalRef - yLocalTar));
182  break;
183  }
184  }
185 }
186 
187 //
188 // -- Analyze
189 //
190 template <typename T>
192  // The default use case is to use vertices from Alpaka reconstructed on CPU and GPU;
193  // The function is left templated if any other cases need to be added
194  analyzeSeparate(tokenSoAHitsReference_, tokenSoAHitsTarget_, iEvent);
195 }
196 
197 //
198 // -- Book Histograms
199 //
200 template <typename T>
202  edm::Run const& iRun,
203  edm::EventSetup const& iSetup) {
204  iBook.cd();
205  iBook.setCurrentFolder(topFolderName_);
206 
207  // clang-format off
208  //Global
209  hnHits_ = iBook.book2I("nHits", "ReferencevsTarget RecHits per event;#Reference RecHits;#Target RecHits", 200, 0, 5000,200, 0, 5000);
210  //Barrel Layer
211  for(unsigned int il=0;il<tkGeom_->numberOfLayers(PixelSubdetector::PixelBarrel);il++){
212  hBchargeL_[il] = iBook.book2I(Form("recHitsBLay%dCharge",il+1), Form("ReferencevsTarget RecHits Charge Barrel Layer%d;Reference Charge;Target Charge",il+1), 250, 0, 100000, 250, 0, 100000);
213  hBsizexL_[il] = iBook.book2I(Form("recHitsBLay%dSizex",il+1), Form("ReferencevsTarget RecHits SizeX Barrel Layer%d;Reference SizeX;Target SizeX",il+1), 30, 0, 30, 30, 0, 30);
214  hBsizeyL_[il] = iBook.book2I(Form("recHitsBLay%dSizey",il+1), Form("ReferencevsTarget RecHits SizeY Barrel Layer%d;Reference SizeY;Target SizeY",il+1), 30, 0, 30, 30, 0, 30);
215  hBposxL_[il] = iBook.book2D(Form("recHitsBLay%dPosx",il+1), Form("ReferencevsTarget RecHits x-pos in Barrel Layer%d;Reference pos x;Target pos x",il+1), 200, -5, 5, 200,-5,5);
216  hBposyL_[il] = iBook.book2D(Form("recHitsBLay%dPosy",il+1), Form("ReferencevsTarget RecHits y-pos in Barrel Layer%d;Reference pos y;Target pos y",il+1), 200, -5, 5, 200,-5,5);
217  }
218  //Endcaps
219  //Endcaps Disk
220  for(int is=0;is<2;is++){
221  int sign=is==0? -1:1;
222  for(unsigned int id=0;id<tkGeom_->numberOfLayers(PixelSubdetector::PixelEndcap);id++){
223  hFchargeD_[is][id] = iBook.book2I(Form("recHitsFDisk%+dCharge",id*sign+sign), Form("ReferencevsTarget RecHits Charge Endcaps Disk%+d;Reference Charge;Target Charge",id*sign+sign), 250, 0, 100000, 250, 0, 100000);
224  hFsizexD_[is][id] = iBook.book2I(Form("recHitsFDisk%+dSizex",id*sign+sign), Form("ReferencevsTarget RecHits SizeX Endcaps Disk%+d;Reference SizeX;Target SizeX",id*sign+sign), 30, 0, 30, 30, 0, 30);
225  hFsizeyD_[is][id] = iBook.book2I(Form("recHitsFDisk%+dSizey",id*sign+sign), Form("ReferencevsTarget RecHits SizeY Endcaps Disk%+d;Reference SizeY;Target SizeY",id*sign+sign), 30, 0, 30, 30, 0, 30);
226  hFposxD_[is][id] = iBook.book2D(Form("recHitsFDisk%+dPosx",id*sign+sign), Form("ReferencevsTarget RecHits x-pos Endcaps Disk%+d;Reference pos x;Target pos x",id*sign+sign), 200, -5, 5, 200, -5, 5);
227  hFposyD_[is][id] = iBook.book2D(Form("recHitsFDisk%+dPosy",id*sign+sign), Form("ReferencevsTarget RecHits y-pos Endcaps Disk%+d;Reference pos y;Target pos y",id*sign+sign), 200, -5, 5, 200, -5, 5);
228  }
229  }
230  //1D differences
231  hBchargeDiff_ = iBook.book1D("rechitChargeDiffBpix","Charge differnce of rechits in BPix; rechit charge difference (Reference - Target)", 101, -50.5, 50.5);
232  hFchargeDiff_ = iBook.book1D("rechitChargeDiffFpix","Charge differnce of rechits in FPix; rechit charge difference (Reference - Target)", 101, -50.5, 50.5);
233  hBsizeXDiff_ = iBook.book1D("rechitsizeXDiffBpix","SizeX difference of rechits in BPix; rechit sizex difference (Reference - Target)", 21, -10.5, 10.5);
234  hFsizeXDiff_ = iBook.book1D("rechitsizeXDiffFpix","SizeX difference of rechits in FPix; rechit sizex difference (Reference - Target)", 21, -10.5, 10.5);
235  hBsizeYDiff_ = iBook.book1D("rechitsizeYDiffBpix","SizeY difference of rechits in BPix; rechit sizey difference (Reference - Target)", 21, -10.5, 10.5);
236  hFsizeYDiff_ = iBook.book1D("rechitsizeYDiffFpix","SizeY difference of rechits in FPix; rechit sizey difference (Reference - Target)", 21, -10.5, 10.5);
237  hBposXDiff_ = iBook.book1D("rechitsposXDiffBpix","x-position difference of rechits in BPix; rechit x-pos difference (Reference - Target)", 1000, -10, 10);
238  hFposXDiff_ = iBook.book1D("rechitsposXDiffFpix","x-position difference of rechits in FPix; rechit x-pos difference (Reference - Target)", 1000, -10, 10);
239  hBposYDiff_ = iBook.book1D("rechitsposYDiffBpix","y-position difference of rechits in BPix; rechit y-pos difference (Reference - Target)", 1000, -10, 10);
240  hFposYDiff_ = iBook.book1D("rechitsposYDiffFpix","y-position difference of rechits in FPix; rechit y-pos difference (Reference - Target)", 1000, -10, 10);
241 }
242 
243 template<typename T>
245  // monitorpixelRecHitsSoAAlpaka
247  desc.add<edm::InputTag>("pixelHitsReferenceSoA", edm::InputTag("siPixelRecHitsPreSplittingAlpakaSerial"));
248  desc.add<edm::InputTag>("pixelHitsTargetSoA", edm::InputTag("siPixelRecHitsPreSplittingAlpaka"));
249  desc.add<std::string>("topFolderName", "SiPixelHeterogeneous/PixelRecHitsCompareDeviceVSHost");
250  desc.add<double>("minD2cut", 0.0001);
251  descriptions.addWithDefaultLabel(desc);
252 }
253 
257 
259 // TODO: change module names to SiPixel*CompareRecHitsSoA when CUDA code is removed
263 
constexpr int32_t ceil(float num)
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > topoToken_
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
MonitorElement * hFsizexD_[2][12]
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
MonitorElement * hBchargeDiff_
MonitorElement * hBchargeL_[4]
MonitorElement * hBposxL_[4]
const std::string topFolderName_
const edm::EDGetTokenT< HitsSoA > tokenSoAHitsReference_
MonitorElement * hBsizeyL_[4]
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > geomToken_
void analyzeSeparate(U tokenRef, V tokenTar, const edm::Event &iEvent)
void bookHistograms(DQMStore::IBooker &ibooker, edm::Run const &iRun, edm::EventSetup const &iSetup) override
const edm::EDGetTokenT< HitsSoA > tokenSoAHitsTarget_
MonitorElement * hFchargeD_[2][12]
~SiPixelCompareRecHits() override=default
int iEvent
Definition: GenABIO.cc:224
uint32_t T const *__restrict__ uint32_t const *__restrict__ int32_t int Histo::index_type cudaStream_t V
MonitorElement * hFchargeDiff_
MonitorElement * hFsizeyD_[2][12]
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
const TrackerTopology * tTopo_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Transition
Definition: Transition.h:12
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
MonitorElement * hBposyL_[4]
MonitorElement * hFposxD_[2][12]
Definition: DetId.h:17
MonitorElement * hBsizexL_[4]
SiPixelCompareRecHits(const edm::ParameterSet &)
static constexpr uint32_t invalidHit_
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:221
HLT enums.
const TrackerGeometry * tkGeom_
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 * book2I(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:305
static constexpr float micron_
void dqmBeginRun(const edm::Run &, const edm::EventSetup &) override
Definition: Run.h:45
MonitorElement * hFposyD_[2][12]