testSiStripHashedDetId

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206
// system includes
#include <algorithm>
#include <iostream>
#include <sstream>
#include <time.h>
#include <vector>
#include <cstdint>

// user includes
#include "CalibFormats/SiStripObjects/interface/SiStripHashedDetId.h"
#include "DataFormats/SiStripCommon/interface/SiStripConstants.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/Framework/interface/ESHandle.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/one/EDAnalyzer.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "Geometry/CommonDetUnit/interface/GeomDet.h"
#include "Geometry/CommonDetUnit/interface/GeomDetType.h"
#include "Geometry/CommonTopologies/interface/StripTopology.h"
#include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
#include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetType.h"
#include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetUnit.h"
#include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"

/**
   @class test_SiStripHashedDetId
   @author R.Bainbridge
   @brief Simple class that tests SiStripHashedDetId.
*/
class testSiStripHashedDetId : public edm::one::EDAnalyzer<> {
public:
  testSiStripHashedDetId(const edm::ParameterSet &);
  ~testSiStripHashedDetId();

  void initialize(edm::EventSetup const &);
  void analyze(const edm::Event &, const edm::EventSetup &);

private:
  //-------------- member data ----------------//
  const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> geomToken_;
};

using namespace sistrip;

// -----------------------------------------------------------------------------
//
testSiStripHashedDetId::testSiStripHashedDetId(const edm::ParameterSet &pset) : geomToken_(esConsumes()) {
  edm::LogVerbatim(mlDqmCommon_) << "[testSiStripHashedDetId::" << __func__ << "]"
                                 << " Constructing object...";
}

// -----------------------------------------------------------------------------
//
testSiStripHashedDetId::~testSiStripHashedDetId() {
  edm::LogVerbatim(mlDqmCommon_) << "[testSiStripHashedDetId::" << __func__ << "]"
                                 << " Destructing object...";
}

// -----------------------------------------------------------------------------
//
void testSiStripHashedDetId::initialize(const edm::EventSetup &setup) {
  edm::LogVerbatim(mlDqmCommon_) << "[SiStripHashedDetId::" << __func__ << "]"
                                 << " Tests the generation of DetId hash map...";

  // Retrieve geometry
  edm::ESHandle geom = setup.getHandle(geomToken_);

  // Build list of DetIds
  std::vector<uint32_t> dets;
  dets.reserve(16000);
  TrackerGeometry::DetContainer::const_iterator iter = geom->detUnits().begin();
  for (; iter != geom->detUnits().end(); ++iter) {
    const auto strip = dynamic_cast<const StripGeomDetUnit *>(*iter);
    if (strip) {
      dets.push_back((strip->geographicalId()).rawId());
    }
  }
  edm::LogVerbatim(mlDqmCommon_) << "[testSiStripHashedDetId::" << __func__ << "]"
                                 << " Retrieved " << dets.size() << " strip DetIds from geometry!";

  // Sorted DetId list gives max performance, anything else is worse
  if (true) {
    std::sort(dets.begin(), dets.end());
  } else {
    std::reverse(dets.begin(), dets.end());
  }

  // Manipulate DetId list
  if (false) {
    if (dets.size() > 4) {
      uint32_t temp = dets.front();
      dets.front() = dets.back();             // swapped
      dets.back() = temp;                     // swapped
      dets.at(1) = 0x00000001;                // wrong
      dets.at(dets.size() - 2) = 0xFFFFFFFF;  // wrong
    }
  }

  // Create hash map
  SiStripHashedDetId hash(dets);
  LogTrace(mlDqmCommon_) << "[testSiStripHashedDetId::" << __func__ << "]"
                         << " DetId hash map: " << std::endl
                         << hash;

  // Manipulate DetId list
  if (false) {
    if (dets.size() > 4) {
      uint32_t temp = dets.front();
      dets.front() = dets.back();             // swapped
      dets.back() = temp;                     // swapped
      dets.at(1) = 0x00000001;                // wrong
      dets.at(dets.size() - 2) = 0xFFFFFFFF;  // wrong
    }
  }

  // Retrieve hashed indices
  std::vector<uint32_t> hashes;
  uint32_t istart = time(NULL);
  for (uint16_t tt = 0; tt < 10000; ++tt) {  // 10000 loops just to see some non-negligible time meaasurement!
    hashes.clear();
    hashes.reserve(dets.size());
    std::vector<uint32_t>::const_iterator idet = dets.begin();
    for (; idet != dets.end(); ++idet) {
      hashes.push_back(hash.hashedIndex(*idet));
    }
  }

  // Some debug
  std::stringstream ss;
  ss << "[testSiStripHashedDetId::" << __func__ << "]";
  std::vector<uint32_t>::const_iterator ii = hashes.begin();
  uint16_t cntr1 = 0;
  for (; ii != hashes.end(); ++ii) {
    if (*ii == sistrip::invalid32_) {
      cntr1++;
      ss << std::endl << " Invalid index " << *ii;
      continue;
    }
    uint32_t detid = hash.unhashIndex(*ii);
    std::vector<uint32_t>::const_iterator iter = find(dets.begin(), dets.end(), detid);
    if (iter == dets.end()) {
      cntr1++;
      ss << std::endl << " Did not find value " << detid << " at index " << ii - hashes.begin() << " in vector!";
    } else if (*ii != static_cast<uint32_t>(iter - dets.begin())) {
      cntr1++;
      ss << std::endl
         << " Found same value " << detid << " at different indices " << *ii << " and " << iter - dets.begin();
    }
  }
  if (cntr1) {
    ss << std::endl << " Found " << cntr1 << " incompatible values!";
  } else {
    ss << " Found no incompatible values!";
  }
  LogTrace(mlDqmCommon_) << ss.str();

  edm::LogVerbatim(mlDqmCommon_) << "[testSiStripHashedDetId::" << __func__ << "]"
                                 << " Processed " << hashes.size() << " DetIds in " << (time(NULL) - istart)
                                 << " seconds";

  // Retrieve DetIds
  std::vector<uint32_t> detids;
  uint32_t jstart = time(NULL);
  for (uint16_t ttt = 0; ttt < 10000; ++ttt) {  // 10000 loops just to see some non-negligible time
                                                // meaasurement!
    detids.clear();
    detids.reserve(dets.size());
    for (uint16_t idet = 0; idet < dets.size(); ++idet) {
      detids.push_back(hash.unhashIndex(idet));
    }
  }

  // Some debug
  std::stringstream sss;
  sss << "[testSiStripHashedDetId::" << __func__ << "]";
  uint16_t cntr2 = 0;
  std::vector<uint32_t>::const_iterator iii = detids.begin();
  for (; iii != detids.end(); ++iii) {
    if (*iii != dets.at(iii - detids.begin())) {
      cntr2++;
      sss << std::endl
          << " Diff values " << *iii << " and " << dets.at(iii - detids.begin()) << " found at index "
          << iii - detids.begin() << " ";
    }
  }
  if (cntr2) {
    sss << std::endl << " Found " << cntr2 << " incompatible values!";
  } else {
    sss << " Found no incompatible values!";
  }
  LogTrace(mlDqmCommon_) << sss.str();

  edm::LogVerbatim(mlDqmCommon_) << "[testSiStripHashedDetId::" << __func__ << "]"
                                 << " Processed " << detids.size() << " hashed indices in " << (time(NULL) - jstart)
                                 << " seconds";
}

// -----------------------------------------------------------------------------
//
void testSiStripHashedDetId::analyze(const edm::Event &event, const edm::EventSetup &setup) {
  initialize(setup);
  LogTrace(mlDqmCommon_) << "[testSiStripHashedDetId::" << __func__ << "]"
                         << " Analyzing run/event " << event.id().run() << "/" << event.id().event();
}

DEFINE_FWK_MODULE(testSiStripHashedDetId);