File indexing completed on 2023-10-25 10:03:53
0001 #include "CLHEP/Units/GlobalPhysicalConstants.h"
0002 #include "SimFastTiming/FastTimingCommon/interface/ETLDeviceSim.h"
0003 #include "DataFormats/Math/interface/GeantUnits.h"
0004
0005 #include "DataFormats/DetId/interface/DetId.h"
0006 #include "DataFormats/ForwardDetId/interface/MTDDetId.h"
0007 #include "DataFormats/ForwardDetId/interface/ETLDetId.h"
0008 #include "FWCore/Framework/interface/ConsumesCollector.h"
0009
0010 #include "Geometry/CommonDetUnit/interface/GeomDetType.h"
0011 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0012 #include "Geometry/MTDGeometryBuilder/interface/ProxyMTDTopology.h"
0013 #include "Geometry/MTDGeometryBuilder/interface/RectangularMTDTopology.h"
0014
0015 ETLDeviceSim::ETLDeviceSim(const edm::ParameterSet& pset, edm::ConsumesCollector iC)
0016 : geomToken_(iC.esConsumes()),
0017 geom_(nullptr),
0018 MIPPerMeV_(1.0 / pset.getParameter<double>("meVPerMIP")),
0019 bxTime_(pset.getParameter<double>("bxTime")),
0020 tofDelay_(pset.getParameter<double>("tofDelay")) {}
0021
0022 void ETLDeviceSim::getEventSetup(const edm::EventSetup& evs) { geom_ = &evs.getData(geomToken_); }
0023
0024 void ETLDeviceSim::getHitsResponse(const std::vector<std::tuple<int, uint32_t, float> >& hitRefs,
0025 const edm::Handle<edm::PSimHitContainer>& hits,
0026 mtd_digitizer::MTDSimHitDataAccumulator* simHitAccumulator,
0027 CLHEP::HepRandomEngine* hre) {
0028 using namespace geant_units::operators;
0029
0030
0031 const int nchits = hitRefs.size();
0032 LogTrace("ETLDeviceSim") << "Processing " << nchits << " SIM hits";
0033 for (int i = 0; i < nchits; ++i) {
0034 const int hitidx = std::get<0>(hitRefs[i]);
0035 const uint32_t id = std::get<1>(hitRefs[i]);
0036 const MTDDetId detId(id);
0037
0038
0039 if (detId.det() != DetId::Forward || detId.mtdSubDetector() != 2) {
0040 throw cms::Exception("ETLDeviceSim")
0041 << "got a DetId that was not ETL: Det = " << detId.det() << " subDet = " << detId.mtdSubDetector();
0042 }
0043
0044 if (id == 0)
0045 continue;
0046
0047 ETLDetId etlid(detId);
0048 DetId geoId = ETLDetId(etlid.mtdSide(), etlid.mtdRR(), etlid.module(), etlid.modType());
0049 const MTDGeomDet* thedet = geom_->idToDet(geoId);
0050 if (thedet == nullptr) {
0051 throw cms::Exception("ETLDeviceSim") << "GeographicalID: " << std::hex << geoId.rawId() << " (" << detId.rawId()
0052 << ") is invalid!" << std::dec << std::endl;
0053 }
0054 const ProxyMTDTopology& topoproxy = static_cast<const ProxyMTDTopology&>(thedet->topology());
0055 const RectangularMTDTopology& topo = static_cast<const RectangularMTDTopology&>(topoproxy.specificTopology());
0056
0057 const float toa = std::get<2>(hitRefs[i]) + tofDelay_;
0058 const PSimHit& hit = hits->at(hitidx);
0059 const float charge = convertGeVToMeV(hit.energyLoss()) * MIPPerMeV_;
0060
0061
0062 const auto& position = hit.localPosition();
0063
0064 Local3DPoint simscaled(convertMmToCm(position.x()), convertMmToCm(position.y()), convertMmToCm(position.z()));
0065
0066
0067 if (!topo.isInPixel(simscaled)) {
0068 LogDebug("ETLDeviceSim") << "Skipping hit ouf of pixel # " << hitidx << " DetId " << etlid.rawId() << " tof "
0069 << toa;
0070 continue;
0071 }
0072 const auto& thepixel = topo.pixel(simscaled);
0073 const uint8_t row(thepixel.first), col(thepixel.second);
0074
0075 LogDebug("ETLDeviceSim") << "Processing hit in pixel # " << hitidx << " DetId " << etlid.rawId() << " row/col "
0076 << (uint32_t)row << " " << (uint32_t)col << " tof " << toa;
0077
0078 auto simHitIt =
0079 simHitAccumulator->emplace(mtd_digitizer::MTDCellId(id, row, col), mtd_digitizer::MTDCellInfo()).first;
0080
0081
0082 const int itime = std::floor(toa / bxTime_) + 9;
0083 if (itime < 0 || itime > 14)
0084 continue;
0085
0086
0087 if (itime >= (int)simHitIt->second.hit_info[0].size())
0088 continue;
0089
0090 (simHitIt->second).hit_info[0][itime] += charge;
0091
0092
0093 const float tof = toa - (itime - 9) * bxTime_;
0094
0095 if ((simHitIt->second).hit_info[1][itime] == 0. || tof < (simHitIt->second).hit_info[1][itime]) {
0096 (simHitIt->second).hit_info[1][itime] = tof;
0097 }
0098 }
0099 }