Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-10-25 09:33:58

0001 // system includes
0002 #include <algorithm>
0003 #include <iostream>
0004 #include <sstream>
0005 #include <time.h>
0006 #include <vector>
0007 #include <cstdint>
0008 
0009 // user includes
0010 #include "CalibFormats/SiStripObjects/interface/SiStripHashedDetId.h"
0011 #include "DataFormats/SiStripCommon/interface/SiStripConstants.h"
0012 #include "FWCore/Framework/interface/MakerMacros.h"
0013 #include "FWCore/Framework/interface/ESHandle.h"
0014 #include "FWCore/Framework/interface/Event.h"
0015 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0017 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
0018 #include "Geometry/CommonDetUnit/interface/GeomDetType.h"
0019 #include "Geometry/CommonTopologies/interface/StripTopology.h"
0020 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0021 #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetType.h"
0022 #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetUnit.h"
0023 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0024 
0025 /**
0026    @class test_SiStripHashedDetId
0027    @author R.Bainbridge
0028    @brief Simple class that tests SiStripHashedDetId.
0029 */
0030 class testSiStripHashedDetId : public edm::one::EDAnalyzer<> {
0031 public:
0032   testSiStripHashedDetId(const edm::ParameterSet &);
0033   ~testSiStripHashedDetId();
0034 
0035   void initialize(edm::EventSetup const &);
0036   void analyze(const edm::Event &, const edm::EventSetup &);
0037 
0038 private:
0039   //-------------- member data ----------------//
0040   const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> geomToken_;
0041 };
0042 
0043 using namespace sistrip;
0044 
0045 // -----------------------------------------------------------------------------
0046 //
0047 testSiStripHashedDetId::testSiStripHashedDetId(const edm::ParameterSet &pset) : geomToken_(esConsumes()) {
0048   edm::LogVerbatim(mlDqmCommon_) << "[testSiStripHashedDetId::" << __func__ << "]"
0049                                  << " Constructing object...";
0050 }
0051 
0052 // -----------------------------------------------------------------------------
0053 //
0054 testSiStripHashedDetId::~testSiStripHashedDetId() {
0055   edm::LogVerbatim(mlDqmCommon_) << "[testSiStripHashedDetId::" << __func__ << "]"
0056                                  << " Destructing object...";
0057 }
0058 
0059 // -----------------------------------------------------------------------------
0060 //
0061 void testSiStripHashedDetId::initialize(const edm::EventSetup &setup) {
0062   edm::LogVerbatim(mlDqmCommon_) << "[SiStripHashedDetId::" << __func__ << "]"
0063                                  << " Tests the generation of DetId hash map...";
0064 
0065   // Retrieve geometry
0066   edm::ESHandle geom = setup.getHandle(geomToken_);
0067 
0068   // Build list of DetIds
0069   std::vector<uint32_t> dets;
0070   dets.reserve(16000);
0071   TrackerGeometry::DetContainer::const_iterator iter = geom->detUnits().begin();
0072   for (; iter != geom->detUnits().end(); ++iter) {
0073     const auto strip = dynamic_cast<const StripGeomDetUnit *>(*iter);
0074     if (strip) {
0075       dets.push_back((strip->geographicalId()).rawId());
0076     }
0077   }
0078   edm::LogVerbatim(mlDqmCommon_) << "[testSiStripHashedDetId::" << __func__ << "]"
0079                                  << " Retrieved " << dets.size() << " strip DetIds from geometry!";
0080 
0081   // Sorted DetId list gives max performance, anything else is worse
0082   if (true) {
0083     std::sort(dets.begin(), dets.end());
0084   } else {
0085     std::reverse(dets.begin(), dets.end());
0086   }
0087 
0088   // Manipulate DetId list
0089   if (false) {
0090     if (dets.size() > 4) {
0091       uint32_t temp = dets.front();
0092       dets.front() = dets.back();             // swapped
0093       dets.back() = temp;                     // swapped
0094       dets.at(1) = 0x00000001;                // wrong
0095       dets.at(dets.size() - 2) = 0xFFFFFFFF;  // wrong
0096     }
0097   }
0098 
0099   // Create hash map
0100   SiStripHashedDetId hash(dets);
0101   LogTrace(mlDqmCommon_) << "[testSiStripHashedDetId::" << __func__ << "]"
0102                          << " DetId hash map: " << std::endl
0103                          << hash;
0104 
0105   // Manipulate DetId list
0106   if (false) {
0107     if (dets.size() > 4) {
0108       uint32_t temp = dets.front();
0109       dets.front() = dets.back();             // swapped
0110       dets.back() = temp;                     // swapped
0111       dets.at(1) = 0x00000001;                // wrong
0112       dets.at(dets.size() - 2) = 0xFFFFFFFF;  // wrong
0113     }
0114   }
0115 
0116   // Retrieve hashed indices
0117   std::vector<uint32_t> hashes;
0118   uint32_t istart = time(NULL);
0119   for (uint16_t tt = 0; tt < 10000; ++tt) {  // 10000 loops just to see some non-negligible time meaasurement!
0120     hashes.clear();
0121     hashes.reserve(dets.size());
0122     std::vector<uint32_t>::const_iterator idet = dets.begin();
0123     for (; idet != dets.end(); ++idet) {
0124       hashes.push_back(hash.hashedIndex(*idet));
0125     }
0126   }
0127 
0128   // Some debug
0129   std::stringstream ss;
0130   ss << "[testSiStripHashedDetId::" << __func__ << "]";
0131   std::vector<uint32_t>::const_iterator ii = hashes.begin();
0132   uint16_t cntr1 = 0;
0133   for (; ii != hashes.end(); ++ii) {
0134     if (*ii == sistrip::invalid32_) {
0135       cntr1++;
0136       ss << std::endl << " Invalid index " << *ii;
0137       continue;
0138     }
0139     uint32_t detid = hash.unhashIndex(*ii);
0140     std::vector<uint32_t>::const_iterator iter = find(dets.begin(), dets.end(), detid);
0141     if (iter == dets.end()) {
0142       cntr1++;
0143       ss << std::endl << " Did not find value " << detid << " at index " << ii - hashes.begin() << " in vector!";
0144     } else if (*ii != static_cast<uint32_t>(iter - dets.begin())) {
0145       cntr1++;
0146       ss << std::endl
0147          << " Found same value " << detid << " at different indices " << *ii << " and " << iter - dets.begin();
0148     }
0149   }
0150   if (cntr1) {
0151     ss << std::endl << " Found " << cntr1 << " incompatible values!";
0152   } else {
0153     ss << " Found no incompatible values!";
0154   }
0155   LogTrace(mlDqmCommon_) << ss.str();
0156 
0157   edm::LogVerbatim(mlDqmCommon_) << "[testSiStripHashedDetId::" << __func__ << "]"
0158                                  << " Processed " << hashes.size() << " DetIds in " << (time(NULL) - istart)
0159                                  << " seconds";
0160 
0161   // Retrieve DetIds
0162   std::vector<uint32_t> detids;
0163   uint32_t jstart = time(NULL);
0164   for (uint16_t ttt = 0; ttt < 10000; ++ttt) {  // 10000 loops just to see some non-negligible time
0165                                                 // meaasurement!
0166     detids.clear();
0167     detids.reserve(dets.size());
0168     for (uint16_t idet = 0; idet < dets.size(); ++idet) {
0169       detids.push_back(hash.unhashIndex(idet));
0170     }
0171   }
0172 
0173   // Some debug
0174   std::stringstream sss;
0175   sss << "[testSiStripHashedDetId::" << __func__ << "]";
0176   uint16_t cntr2 = 0;
0177   std::vector<uint32_t>::const_iterator iii = detids.begin();
0178   for (; iii != detids.end(); ++iii) {
0179     if (*iii != dets.at(iii - detids.begin())) {
0180       cntr2++;
0181       sss << std::endl
0182           << " Diff values " << *iii << " and " << dets.at(iii - detids.begin()) << " found at index "
0183           << iii - detids.begin() << " ";
0184     }
0185   }
0186   if (cntr2) {
0187     sss << std::endl << " Found " << cntr2 << " incompatible values!";
0188   } else {
0189     sss << " Found no incompatible values!";
0190   }
0191   LogTrace(mlDqmCommon_) << sss.str();
0192 
0193   edm::LogVerbatim(mlDqmCommon_) << "[testSiStripHashedDetId::" << __func__ << "]"
0194                                  << " Processed " << detids.size() << " hashed indices in " << (time(NULL) - jstart)
0195                                  << " seconds";
0196 }
0197 
0198 // -----------------------------------------------------------------------------
0199 //
0200 void testSiStripHashedDetId::analyze(const edm::Event &event, const edm::EventSetup &setup) {
0201   initialize(setup);
0202   LogTrace(mlDqmCommon_) << "[testSiStripHashedDetId::" << __func__ << "]"
0203                          << " Analyzing run/event " << event.id().run() << "/" << event.id().event();
0204 }
0205 
0206 DEFINE_FWK_MODULE(testSiStripHashedDetId);