Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-07 04:35:04

0001 // -*- C++ -*-
0002 //
0003 // Package:    Calibration/TkAlCaRecoProducers
0004 // Class:      NearbyPixelClustersAnalyzer
0005 //
0006 /**\class NearbyPixelClustersAnalyzer NearbyPixelClustersAnalyzer.cc Calibration/TkAlCaRecoProducers/plugins/NearbyPixelClustersAnalyzer.cc
0007 
0008  Description: Analysis of the closebyPixelClusters collections
0009 */
0010 //
0011 // Original Author:  Marco Musich
0012 //         Created:  Tue, 30 Mar 2021 09:22:07 GMT
0013 //
0014 //
0015 
0016 // system include files
0017 #include <memory>
0018 
0019 // user include files
0020 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0021 #include "CalibTracker/SiPixelESProducers/interface/SiPixelDetInfoFileReader.h"
0022 #include "FWCore/ServiceRegistry/interface/Service.h"
0023 #include "DataFormats/Common/interface/DetSetVector.h"
0024 #include "DataFormats/Common/interface/DetSetVectorNew.h"
0025 #include "DataFormats/DetId/interface/DetId.h"
0026 #include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
0027 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0028 #include "DataFormats/TrackReco/interface/Track.h"
0029 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0030 #include "DataFormats/TrackerCommon/interface/PixelBarrelName.h"
0031 #include "DataFormats/TrackerCommon/interface/PixelEndcapName.h"
0032 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0033 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
0034 #include "FWCore/Framework/interface/Event.h"
0035 #include "FWCore/Framework/interface/Frameworkfwd.h"
0036 #include "FWCore/Framework/interface/MakerMacros.h"
0037 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0038 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0039 #include "FWCore/Utilities/interface/StreamID.h"
0040 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0041 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0042 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0043 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0044 #include "Geometry/TrackerNumberingBuilder/interface/GeometricDet.h"
0045 #include "RecoLocalTracker/ClusterParameterEstimator/interface/PixelClusterParameterEstimator.h"
0046 #include "RecoLocalTracker/Records/interface/TkPixelCPERecord.h"
0047 #include "TrackingTools/PatternTools/interface/TrajTrackAssociation.h"
0048 #include "TrackingTools/TrackFitters/interface/TrajectoryStateCombiner.h"
0049 #include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
0050 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
0051 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHit.h"
0052 
0053 #include "TH2F.h"
0054 //
0055 // class declaration
0056 //
0057 
0058 class NearbyPixelClustersAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources, edm::one::WatchRuns> {
0059 public:
0060   explicit NearbyPixelClustersAnalyzer(const edm::ParameterSet&);
0061   ~NearbyPixelClustersAnalyzer() override = default;
0062 
0063   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0064 
0065 private:
0066   void beginJob() override;
0067   void analyze(const edm::Event&, const edm::EventSetup&) override;
0068   void beginRun(edm::Run const&, edm::EventSetup const&) override;
0069   void endRun(edm::Run const&, edm::EventSetup const&) override {}
0070   std::map<uint32_t, TH2F*> bookModuleHistograms(const TrackerTopology* tTopo);
0071   std::tuple<std::string, int, int, int> setTopoInfo(uint32_t detId, const TrackerTopology* tTopo);
0072   void endJob() override;
0073 
0074   void countClusters(const edm::Handle<SiPixelClusterCollectionNew>& handle,
0075                      //const TrackerGeometry* tkGeom_,
0076                      unsigned int& nClusGlobal);
0077 
0078   bool detidIsOnPixel(const DetId& detid);
0079   TrajectoryStateOnSurface getTrajectoryStateOnSurface(const TrajectoryMeasurement& measurement);
0080   std::pair<float, float> findClosestCluster(const edm::Handle<SiPixelClusterCollectionNew>& handle,
0081                                              const PixelClusterParameterEstimator* pixelCPE_,
0082                                              const TrackerGeometry* trackerGeometry_,
0083                                              uint32_t rawId,
0084                                              float traj_lx,
0085                                              float traj_ly);
0086 
0087   void fillClusterFrames(const edm::Handle<SiPixelClusterCollectionNew>& handle);
0088 
0089   // ----------member data ---------------------------
0090   const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> geomEsToken_;
0091   const edm::ESGetToken<PixelClusterParameterEstimator, TkPixelCPERecord> pixelCPEEsToken_;
0092   const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> geomEsTokenBR_;
0093   const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> topoEsTokenBR_;
0094 
0095   edm::EDGetTokenT<SiPixelClusterCollectionNew> clustersToken_;
0096   edm::EDGetTokenT<SiPixelClusterCollectionNew> nearByClustersToken_;
0097   edm::EDGetTokenT<TrajTrackAssociationCollection> trajTrackCollectionToken_;
0098   edm::EDGetTokenT<edm::ValueMap<std::vector<float>>> distanceToken_;
0099   edm::EDGetTokenT<edm::View<reco::Track>> muonTracksToken_;
0100 
0101   edm::Service<TFileService> fs;
0102 
0103   TH1I* h_nALCARECOClusters;
0104   TH1I* h_nCloseByClusters;
0105   TH1F* h_distClosestValid;
0106   TH1F* h_distClosestMissing;
0107   TH1F* h_distClosestInactive;
0108   TH1F* h_distClosestTrack;
0109 
0110   SiPixelDetInfoFileReader reader_;
0111   std::map<std::string, TFileDirectory> outputFolders_;
0112   std::map<uint32_t, TH2F*> histoMap_;
0113   bool phase_;
0114 };
0115 
0116 //
0117 // constructors and destructor
0118 //
0119 NearbyPixelClustersAnalyzer::NearbyPixelClustersAnalyzer(const edm::ParameterSet& iConfig)
0120     : geomEsToken_(esConsumes()),
0121       pixelCPEEsToken_(esConsumes(edm::ESInputTag("", "PixelCPEGeneric"))),
0122       geomEsTokenBR_(esConsumes<edm::Transition::BeginRun>()),
0123       topoEsTokenBR_(esConsumes<edm::Transition::BeginRun>()),
0124       clustersToken_(consumes<SiPixelClusterCollectionNew>(iConfig.getParameter<edm::InputTag>("clusterCollection"))),
0125       nearByClustersToken_(
0126           consumes<SiPixelClusterCollectionNew>(iConfig.getParameter<edm::InputTag>("nearByClusterCollection"))),
0127       trajTrackCollectionToken_(
0128           consumes<TrajTrackAssociationCollection>(iConfig.getParameter<edm::InputTag>("trajectoryInput"))),
0129       distanceToken_(consumes<edm::ValueMap<std::vector<float>>>(iConfig.getParameter<edm::InputTag>("distToTrack"))),
0130       muonTracksToken_(consumes<edm::View<reco::Track>>(iConfig.getParameter<edm::InputTag>("muonTracks"))),
0131       reader_(edm::FileInPath(iConfig.getParameter<std::string>("skimmedGeometryPath")).fullPath()) {
0132   usesResource(TFileService::kSharedResource);
0133 }
0134 
0135 //
0136 // member functions
0137 //
0138 
0139 // ------------ method called for each event  ------------
0140 void NearbyPixelClustersAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0141   using namespace edm;
0142 
0143   // get the Tracker geometry from event setup
0144   const TrackerGeometry* trackerGeometry_ = &iSetup.getData(geomEsToken_);
0145 
0146   // get the Pixel CPE from event setup
0147   const PixelClusterParameterEstimator* pixelCPE_ = &iSetup.getData(pixelCPEEsToken_);
0148 
0149   // get the muon track collection
0150   const auto& muonTrackCollectionHandle = iEvent.getHandle(muonTracksToken_);
0151   auto const& muonTracks = *muonTrackCollectionHandle;
0152 
0153   // get the track distances
0154   const auto& distancesToTrack = iEvent.getHandle(distanceToken_);
0155 
0156   unsigned int nMuons = muonTracks.size();
0157   for (unsigned int ij = 0; ij < nMuons; ij++) {
0158     auto muon = muonTrackCollectionHandle->ptrAt(ij);
0159     edm::RefToBase<reco::Track> trackRef = muonTrackCollectionHandle->refAt(ij);
0160     const auto& distances = (*distancesToTrack)[trackRef];
0161 
0162     LogDebug("NearbyPixelClustersAnalyzer") << "distances size: " << distances.size() << std::endl;
0163 
0164     unsigned counter = 0;
0165     double closestDR = 999.;
0166     for (const auto& distance : distances) {
0167       counter++;
0168       LogDebug("NearbyPixelClustersAnalyzer")
0169           << "track: " << counter << " distance:" << std::sqrt(distance) << std::endl;
0170       if (distance < closestDR && distance > 0) {
0171         closestDR = distance;
0172       }
0173     }
0174 
0175     h_distClosestTrack->Fill(std::sqrt(closestDR));
0176   }
0177 
0178   // Get cluster collection
0179   const auto& clusterCollectionHandle = iEvent.getHandle(clustersToken_);
0180 
0181   unsigned int nClusGlobal = 0;
0182   countClusters(clusterCollectionHandle, nClusGlobal);
0183 
0184   h_nALCARECOClusters->Fill(nClusGlobal);
0185   edm::LogInfo("NearbyPixelClustersAnalyzer") << "total ALCARECO clusters: " << nClusGlobal << std::endl;
0186 
0187   // Get nearby cluster collection
0188   const auto& nearByClusterCollectionHandle = iEvent.getHandle(nearByClustersToken_);
0189 
0190   unsigned int nNearByClusGlobal = 0;
0191   countClusters(nearByClusterCollectionHandle, /*trackerGeometry_,*/ nNearByClusGlobal);
0192 
0193   h_nCloseByClusters->Fill(nNearByClusGlobal);
0194   edm::LogInfo("NearbyPixelClustersAnalyzer") << "total close-by clusters: " << nNearByClusGlobal << std::endl;
0195 
0196   // fill the detector frames
0197   fillClusterFrames(nearByClusterCollectionHandle);
0198 
0199   // Get Traj-Track Collection
0200   const auto& trajTrackCollectionHandle = iEvent.getHandle(trajTrackCollectionToken_);
0201 
0202   if (!trajTrackCollectionHandle.isValid())
0203     return;
0204 
0205   for (const auto& pair : *trajTrackCollectionHandle) {
0206     const edm::Ref<std::vector<Trajectory>> traj = pair.key;
0207     const reco::TrackRef track = pair.val;
0208 
0209     for (const TrajectoryMeasurement& measurement : traj->measurements()) {
0210       if (!measurement.updatedState().isValid())
0211         return;
0212 
0213       const TransientTrackingRecHit::ConstRecHitPointer& recHit = measurement.recHit();
0214 
0215       // Only looking for pixel hits
0216       DetId r_rawId = recHit->geographicalId();
0217 
0218       if (!this->detidIsOnPixel(r_rawId))
0219         continue;
0220 
0221       // Skipping hits with undeterminable positions
0222       TrajectoryStateOnSurface trajStateOnSurface = this->getTrajectoryStateOnSurface(measurement);
0223 
0224       if (!(trajStateOnSurface.isValid()))
0225         continue;
0226 
0227       // Position measurements
0228       // Looking for valid and missing hits
0229       LocalPoint localPosition = trajStateOnSurface.localPosition();
0230 
0231       const auto& traj_lx = localPosition.x();
0232       const auto& traj_ly = localPosition.y();
0233 
0234       const auto loc = this->findClosestCluster(
0235           nearByClusterCollectionHandle, pixelCPE_, trackerGeometry_, r_rawId.rawId(), traj_lx, traj_ly);
0236 
0237       float dist = (loc.first != -999.) ? std::sqrt(loc.first * loc.first + loc.second * loc.second) : -0.1;
0238 
0239       if (recHit->getType() == TrackingRecHit::valid) {
0240         edm::LogInfo("NearbyPixelClustersAnalyzer")
0241             << "RawID:" << r_rawId.rawId() << " (valid hit), distance: " << dist << std::endl;
0242         h_distClosestValid->Fill(dist);
0243       }
0244 
0245       if (recHit->getType() == TrackingRecHit::missing) {
0246         edm::LogInfo("NearbyPixelClustersAnalyzer")
0247             << "RawID:" << r_rawId.rawId() << " (missing hit), distance: " << dist << std::endl;
0248         h_distClosestMissing->Fill(dist);
0249       }
0250 
0251       if (recHit->getType() == TrackingRecHit::inactive) {
0252         edm::LogInfo("NearbyPixelClustersAnalyzer")
0253             << "RawID:" << r_rawId.rawId() << " (inactive hit), distance: " << dist << std::endl;
0254         h_distClosestInactive->Fill(dist);
0255       }
0256     }
0257   }
0258 }
0259 
0260 // ------------ method called once each job just before starting event loop  ------------
0261 void NearbyPixelClustersAnalyzer::beginJob() {
0262   TFileDirectory ClusterCounts = fs->mkdir("ClusterCounts");
0263   h_nALCARECOClusters = ClusterCounts.make<TH1I>(
0264       "h_nALCARECOClusters", "Number of Pixel clusters per event (ALCARECO) ;N_{clusters};events", 20, 0, 20);
0265   h_nCloseByClusters = ClusterCounts.make<TH1I>(
0266       "h_nCloseByClusters", "Number of Pixel clusters per event (close-by) ;N_{clusters};events", 20, 0, 20);
0267 
0268   TFileDirectory Distances = fs->mkdir("TrajDistance");
0269   h_distClosestValid = Distances.make<TH1F>(
0270       "h_distClosestValid",
0271       "Distance of Closest cluster to trajectory (valid);distance (cm); valid trajectory measurements",
0272       110,
0273       -0.105,
0274       0.995);
0275   h_distClosestMissing = Distances.make<TH1F>(
0276       "h_distClosestMissing",
0277       "Distance of Closest cluster to trajectory (missing);distance (cm);missing trajectory measurements",
0278       110,
0279       -0.105,
0280       0.995);
0281   h_distClosestInactive = Distances.make<TH1F>(
0282       "h_distClosestInactive",
0283       "Distance of Closest cluster to trajectory (inactive);distance (cm);inactive trajectory measurements",
0284       110,
0285       -0.105,
0286       0.995);
0287 
0288   TFileDirectory TkDistances = fs->mkdir("OtherTrackDistance");
0289   h_distClosestTrack = TkDistances.make<TH1F>(
0290       "h_distClosestTrack",
0291       "#DeltaR Distance of Closest track to the muon trajectory;#DeltaR distance; muon trajectories",
0292       100,
0293       0.,
0294       5.);
0295 }
0296 
0297 // ------------ method called once each job just after ending the event loop  ------------
0298 void NearbyPixelClustersAnalyzer::endJob() {
0299   // please remove this method if not needed
0300 }
0301 
0302 /*--------------------------------------------------------------------*/
0303 bool NearbyPixelClustersAnalyzer::detidIsOnPixel(const DetId& detid)
0304 /*--------------------------------------------------------------------*/
0305 {
0306   if (detid.det() != DetId::Tracker)
0307     return false;
0308   if (detid.subdetId() == PixelSubdetector::PixelBarrel)
0309     return true;
0310   if (detid.subdetId() == PixelSubdetector::PixelEndcap)
0311     return true;
0312   return false;
0313 }
0314 
0315 /*--------------------------------------------------------------------*/
0316 TrajectoryStateOnSurface NearbyPixelClustersAnalyzer::getTrajectoryStateOnSurface(
0317     const TrajectoryMeasurement& measurement)
0318 /*--------------------------------------------------------------------*/
0319 {
0320   const static TrajectoryStateCombiner trajStateCombiner;
0321 
0322   const auto& forwardPredictedState = measurement.forwardPredictedState();
0323   const auto& backwardPredictedState = measurement.backwardPredictedState();
0324 
0325   if (forwardPredictedState.isValid() && backwardPredictedState.isValid())
0326     return trajStateCombiner(forwardPredictedState, backwardPredictedState);
0327 
0328   else if (backwardPredictedState.isValid())
0329     return backwardPredictedState;
0330 
0331   else if (forwardPredictedState.isValid())
0332     return forwardPredictedState;
0333 
0334   edm::LogError("NearbyPixelClusterProducer") << "Error saving traj. measurement data."
0335                                               << " Trajectory state on surface cannot be determined." << std::endl;
0336 
0337   return TrajectoryStateOnSurface();
0338 }
0339 
0340 /*--------------------------------------------------------------------*/
0341 void NearbyPixelClustersAnalyzer::countClusters(const edm::Handle<SiPixelClusterCollectionNew>& handle,
0342                                                 //const TrackerGeometry* tkGeom_,
0343                                                 unsigned int& nClusGlobal)
0344 /*--------------------------------------------------------------------*/
0345 {
0346   for (const auto& DSVItr : *handle) {
0347     uint32_t rawid(DSVItr.detId());
0348     DetId detId(rawid);
0349     LogDebug("NearbyPixelClustersAnalyzer") << "DetId: " << detId.rawId() << " size: " << DSVItr.size() << std::endl;
0350     nClusGlobal += DSVItr.size();
0351   }
0352 }
0353 
0354 /*--------------------------------------------------------------------*/
0355 std::pair<float, float> NearbyPixelClustersAnalyzer::findClosestCluster(
0356     const edm::Handle<SiPixelClusterCollectionNew>& handle,
0357     const PixelClusterParameterEstimator* pixelCPE_,
0358     const TrackerGeometry* trackerGeometry_,
0359     uint32_t rawId,
0360     float traj_lx,
0361     float traj_ly)
0362 /*--------------------------------------------------------------------*/
0363 {
0364   const SiPixelClusterCollectionNew& clusterCollection = *handle;
0365   SiPixelClusterCollectionNew::const_iterator itClusterSet = clusterCollection.begin();
0366 
0367   float minD = 10000.;
0368 
0369   auto loc = std::make_pair(-999., -999.);
0370 
0371   for (; itClusterSet != clusterCollection.end(); itClusterSet++) {
0372     DetId detId(itClusterSet->id());
0373     if (detId.rawId() != rawId)
0374       continue;
0375 
0376     unsigned int subDetId = detId.subdetId();
0377     if (subDetId != PixelSubdetector::PixelBarrel && subDetId != PixelSubdetector::PixelEndcap) {
0378       edm::LogError("NearByPixelClustersAnalyzer")
0379           << "ERROR: not a pixel cluster!!!" << std::endl;  // should not happen
0380       continue;
0381     }
0382 
0383     const PixelGeomDetUnit* pixdet = (const PixelGeomDetUnit*)trackerGeometry_->idToDetUnit(detId);
0384     edmNew::DetSet<SiPixelCluster>::const_iterator itCluster = itClusterSet->begin();
0385     for (; itCluster != itClusterSet->end(); ++itCluster) {
0386       LocalPoint lp(itCluster->x(), itCluster->y(), 0.);
0387       PixelClusterParameterEstimator::ReturnType params = pixelCPE_->getParameters(*itCluster, *pixdet);
0388       lp = std::get<0>(params);
0389 
0390       float D = sqrt((lp.x() - traj_lx) * (lp.x() - traj_lx) + (lp.y() - traj_ly) * (lp.y() - traj_ly));
0391       if (D < minD) {
0392         minD = D;
0393         loc.first = (lp.x() - traj_lx);
0394         loc.second = (lp.y() - traj_ly);
0395       }
0396     }  // loop on cluster sets
0397   }
0398   return loc;
0399 }
0400 
0401 void NearbyPixelClustersAnalyzer::fillClusterFrames(const edm::Handle<SiPixelClusterCollectionNew>& handle) {
0402   const SiPixelClusterCollectionNew& clusterCollection = *handle;
0403   SiPixelClusterCollectionNew::const_iterator itClusterSet = clusterCollection.begin();
0404 
0405   for (; itClusterSet != clusterCollection.end(); itClusterSet++) {
0406     DetId detId(itClusterSet->id());
0407 
0408     edmNew::DetSet<SiPixelCluster>::const_iterator itCluster = itClusterSet->begin();
0409     for (; itCluster != itClusterSet->end(); ++itCluster) {
0410       const std::vector<SiPixelCluster::Pixel> pixelsVec = (*itCluster).pixels();
0411       for (unsigned int i = 0; i < pixelsVec.size(); ++i) {
0412         float pixx = pixelsVec[i].x;  // index as float=iteger, row index
0413         float pixy = pixelsVec[i].y;  // same, col index
0414         float pixel_charge = pixelsVec[i].adc;
0415         histoMap_[detId.rawId()]->Fill(pixy, pixx, pixel_charge);
0416       }
0417     }
0418   }
0419 }
0420 
0421 // ------------ method called for each run  ------------------------------------------
0422 void NearbyPixelClustersAnalyzer::beginRun(const edm::Run& iRun, edm::EventSetup const& iSetup)
0423 /*-----------------------------------------------------------------------------------*/
0424 {
0425   edm::LogInfo("NearbyPixelClustersAnalyzer")
0426       << "@SUB=NearbyPixelClustersAnalyzer::beginRun() before booking histoMap_.size(): " << histoMap_.size()
0427       << std::endl;
0428 
0429   const TrackerTopology* tTopo_ = &iSetup.getData(topoEsTokenBR_);
0430   const TrackerGeometry* pDD_ = &iSetup.getData(geomEsTokenBR_);
0431 
0432   if ((pDD_->isThere(GeomDetEnumerators::P1PXB)) || (pDD_->isThere(GeomDetEnumerators::P1PXEC))) {
0433     phase_ = true;
0434   } else {
0435     phase_ = false;
0436   }
0437 
0438   unsigned nPixelDets = 0;
0439   for (const auto& it : pDD_->detUnits()) {
0440     const PixelGeomDetUnit* mit = dynamic_cast<PixelGeomDetUnit const*>(it);
0441     if (mit != nullptr) {
0442       nPixelDets++;
0443     }
0444   }
0445 
0446   const auto& detIds = reader_.getAllDetIds();
0447   if (detIds.size() != nPixelDets) {
0448     throw cms::Exception("Inconsistent Data")
0449         << "The size of the detId list specified from file (" << detIds.size()
0450         << ") differs from the one in TrackerGeometry (" << nPixelDets << ")! Please cross-check" << std::endl;
0451   }
0452 
0453   for (const auto& it : detIds) {
0454     auto topolInfo = setTopoInfo(it, tTopo_);
0455 
0456     std::string thePart = std::get<0>(topolInfo);
0457 
0458     // book the TFileDirectory if it's not already done
0459     if (!outputFolders_.count(thePart)) {
0460       LogDebug("NearbyPixelClustersAnalyzer") << "booking " << thePart << std::endl;
0461       outputFolders_[thePart] = fs->mkdir(thePart);
0462     }
0463   }
0464 
0465   if (histoMap_.empty()) {
0466     histoMap_ = bookModuleHistograms(tTopo_);
0467   }
0468 
0469   edm::LogInfo("NearbyPixelClusterAnalyzer")
0470       << "@SUB=NearbyPixelClusterAnalyzer::beginRun() After booking histoMap_.size(): " << histoMap_.size()
0471       << std::endl;
0472 }
0473 
0474 // ------------ method called to determine the topology  ------------
0475 std::tuple<std::string, int, int, int> NearbyPixelClustersAnalyzer::setTopoInfo(uint32_t detId,
0476                                                                                 const TrackerTopology* tTopo)
0477 /*-----------------------------------------------------------------------------------*/
0478 {
0479   int subdetId_(-999), layer_(-999), side_(-999);
0480   std::string ret = "";
0481 
0482   subdetId_ = DetId(detId).subdetId();
0483   switch (subdetId_) {
0484     case PixelSubdetector::PixelBarrel:  // PXB
0485       layer_ = tTopo->pxbLayer(detId);
0486       side_ = 0;
0487       ret += Form("BPix_Layer%i", layer_);
0488       break;
0489     case PixelSubdetector::PixelEndcap:  //PXF
0490       side_ = tTopo->pxfSide(detId);
0491       layer_ = tTopo->pxfDisk(detId);
0492       ret += ("FPix_");
0493       ret += (side_ == 1) ? Form("P_disk%i", layer_) : Form("M_disk%i", layer_);
0494       break;
0495     default:
0496       edm::LogError("NearbyPixelClusterAnalyzer") << "we should never be here!" << std::endl;
0497       break;
0498   }
0499 
0500   return std::make_tuple(ret, subdetId_, layer_, side_);
0501 }
0502 
0503 /* ------------ method called once to book all the module level histograms  ---------*/
0504 std::map<uint32_t, TH2F*> NearbyPixelClustersAnalyzer::bookModuleHistograms(const TrackerTopology* tTopo_)
0505 /*-----------------------------------------------------------------------------------*/
0506 {
0507   std::map<uint32_t, TH2F*> hd;
0508 
0509   const auto& detIds = reader_.getAllDetIds();
0510   for (const auto& it : detIds) {
0511     // check if det id is correct and if it is actually cabled in the detector
0512     if (it == 0 || it == 0xFFFFFFFF) {
0513       edm::LogError("DetIdNotGood") << "@SUB=analyze"
0514                                     << "Wrong det id: " << it << "  ... neglecting!" << std::endl;
0515       continue;
0516     }
0517 
0518     auto topolInfo = setTopoInfo(it, tTopo_);
0519     std::string thePart = std::get<0>(topolInfo);
0520 
0521     unsigned int nCols = reader_.getDetUnitDimensions(it).first;
0522     unsigned int nRows = reader_.getDetUnitDimensions(it).second;
0523 
0524     int subdetId = DetId(it).subdetId();
0525 
0526     std::string moduleName = (subdetId == PixelSubdetector::PixelBarrel) ? PixelBarrelName(it, tTopo_, phase_).name()
0527                                                                          : PixelEndcapName(it, tTopo_, phase_).name();
0528 
0529     hd[it] = outputFolders_[thePart].make<TH2F>(
0530         Form("ClusterFrame_%s", moduleName.c_str()),
0531         Form("Cluster Map for module %s;n. cols;n. rows;pixel charge [ADC counts]", moduleName.c_str()),
0532         nCols,
0533         0,
0534         nCols,
0535         nRows,
0536         0,
0537         nRows);
0538   }
0539 
0540   return hd;
0541 }
0542 
0543 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0544 void NearbyPixelClustersAnalyzer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0545   edm::ParameterSetDescription desc;
0546   desc.setComment("Analysis of the closebyPixelClusters collections");
0547   desc.add<edm::InputTag>("clusterCollection", edm::InputTag("ALCARECOSiPixelCalSingleMuonTight"));
0548   desc.add<edm::InputTag>("nearByClusterCollection", edm::InputTag("closebyPixelClusters"));
0549   desc.add<edm::InputTag>("trajectoryInput", edm::InputTag("refittedTracks"));
0550   desc.add<edm::InputTag>("muonTracks", edm::InputTag("ALCARECOSiPixelCalSingleMuonTight"));
0551   desc.add<edm::InputTag>("distToTrack", edm::InputTag("trackDistances"));
0552   desc.add<std::string>("skimmedGeometryPath",
0553                         "SLHCUpgradeSimulations/Geometry/data/PhaseI/PixelSkimmedGeometry_phase1.txt");
0554   descriptions.addWithDefaultLabel(desc);
0555 }
0556 
0557 //define this as a plug-in
0558 DEFINE_FWK_MODULE(NearbyPixelClustersAnalyzer);