CMS 3D CMS Logo

BTLDeviceSim.cc
Go to the documentation of this file.
7 
12 
13 #include "CLHEP/Random/RandGaussQ.h"
14 
16  : geomToken_(iC.esConsumes()),
17  topoToken_(iC.esConsumes()),
18  geom_(nullptr),
19  topo_(nullptr),
20  bxTime_(pset.getParameter<double>("bxTime")),
21  LightYield_(pset.getParameter<double>("LightYield")),
22  LightCollEff_(pset.getParameter<double>("LightCollectionEff")),
23  LightCollSlopeR_(pset.getParameter<double>("LightCollectionSlopeR")),
24  LightCollSlopeL_(pset.getParameter<double>("LightCollectionSlopeL")),
25  PDE_(pset.getParameter<double>("PhotonDetectionEff")) {}
26 
28  geom_ = &evs.getData(geomToken_);
29  topo_ = &evs.getData(topoToken_);
30 }
31 
32 void BTLDeviceSim::getHitsResponse(const std::vector<std::tuple<int, uint32_t, float> >& hitRefs,
34  mtd_digitizer::MTDSimHitDataAccumulator* simHitAccumulator,
35  CLHEP::HepRandomEngine* hre) {
36  using namespace geant_units::operators;
37 
38  //loop over sorted simHits
39  for (auto const& hitRef : hitRefs) {
40  const int hitidx = std::get<0>(hitRef);
41  const uint32_t id = std::get<1>(hitRef);
42  const MTDDetId detId(id);
43  const PSimHit& hit = hits->at(hitidx);
44 
45  // --- Safety check on the detector ID
46  if (detId.det() != DetId::Forward || detId.mtdSubDetector() != 1)
47  continue;
48 
49  if (id == 0)
50  continue; // to be ignored at RECO level
51 
52  BTLDetId btlid(detId);
54  const MTDGeomDet* thedet = geom_->idToDet(geoId);
55 
56  if (thedet == nullptr) {
57  throw cms::Exception("BTLDeviceSim") << "GeographicalID: " << std::hex << geoId.rawId() << " (" << detId.rawId()
58  << ") is invalid!" << std::dec << std::endl;
59  }
60  const ProxyMTDTopology& topoproxy = static_cast<const ProxyMTDTopology&>(thedet->topology());
61  const RectangularMTDTopology& topo = static_cast<const RectangularMTDTopology&>(topoproxy.specificTopology());
62  // calculate the simhit row and column
63  const auto& pentry = hit.entryPoint();
64  Local3DPoint simscaled(convertMmToCm(pentry.x()), convertMmToCm(pentry.y()), convertMmToCm(pentry.z()));
65  // translate from crystal-local coordinates to module-local coordinates to get the row and column
66  simscaled = topo.pixelToModuleLocalPoint(simscaled, btlid.row(topo.nrows()), btlid.column(topo.nrows()));
67  const auto& thepixel = topo.pixel(simscaled);
68  uint8_t row(thepixel.first), col(thepixel.second);
69 
70  if (btlid.row(topo.nrows()) != row || btlid.column(topo.nrows()) != col) {
71  edm::LogWarning("BTLDeviceSim") << "BTLDetId (row,column): (" << btlid.row(topo.nrows()) << ','
72  << btlid.column(topo.nrows()) << ") is not equal to "
73  << "topology (row,column): (" << uint32_t(row) << ',' << uint32_t(col)
74  << "), overriding to detid";
75  row = btlid.row(topo.nrows());
76  col = btlid.column(topo.nrows());
77  }
78 
79  // --- Store the detector element ID as a key of the MTDSimHitDataAccumulator map
80  auto simHitIt =
81  simHitAccumulator->emplace(mtd_digitizer::MTDCellId(id, row, col), mtd_digitizer::MTDCellInfo()).first;
82 
83  // --- Get the simHit energy and convert it from MeV to photo-electrons
84  float Npe = convertGeVToMeV(hit.energyLoss()) * LightYield_ * LightCollEff_ * PDE_;
85 
86  // --- Calculate the light propagation time to the crystal bases (labeled L and R)
87  double distR = 0.5 * topo.pitch().first - convertMmToCm(hit.localPosition().x());
88  double distL = 0.5 * topo.pitch().first + convertMmToCm(hit.localPosition().x());
89 
90  double tR = std::get<2>(hitRef) + LightCollSlopeR_ * distR;
91  double tL = std::get<2>(hitRef) + LightCollSlopeL_ * distL;
92 
93  // --- Accumulate in 15 buckets of 25ns (9 pre-samples, 1 in-time, 5 post-samples)
94  const int iBXR = std::floor(tR / bxTime_) + mtd_digitizer::kInTimeBX;
95  const int iBXL = std::floor(tL / bxTime_) + mtd_digitizer::kInTimeBX;
96 
97  // --- Right side
98  if (iBXR > 0 && iBXR < mtd_digitizer::kNumberOfBX) {
99  // Accumulate the energy of simHits in the same crystal
100  (simHitIt->second).hit_info[0][iBXR] += Npe;
101 
102  // Store the time of the first SimHit in the i-th BX
103  if ((simHitIt->second).hit_info[1][iBXR] == 0 || tR < (simHitIt->second).hit_info[1][iBXR])
104  (simHitIt->second).hit_info[1][iBXR] = tR - (iBXR - mtd_digitizer::kInTimeBX) * bxTime_;
105  }
106 
107  // --- Left side
108  if (iBXL > 0 && iBXL < mtd_digitizer::kNumberOfBX) {
109  // Accumulate the energy of simHits in the same crystal
110  (simHitIt->second).hit_info[2][iBXL] += Npe;
111 
112  // Store the time of the first SimHit in the i-th BX
113  if ((simHitIt->second).hit_info[3][iBXL] == 0 || tL < (simHitIt->second).hit_info[3][iBXL])
114  (simHitIt->second).hit_info[3][iBXL] = tL - (iBXL - mtd_digitizer::kInTimeBX) * bxTime_;
115  }
116 
117  } // hitRef loop
118 }
BTLDeviceSim::LightYield_
const float LightYield_
Definition: BTLDeviceSim.h:43
BTLDeviceSim::geom_
const MTDGeometry * geom_
Definition: BTLDeviceSim.h:39
BTLDetId.h
hfClusterShapes_cfi.hits
hits
Definition: hfClusterShapes_cfi.py:5
BTLDeviceSim.h
mtd_digitizer::MTDCellInfo
Definition: MTDDigitizerTypes.h:17
DetId::det
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
cuy.col
col
Definition: cuy.py:1009
TrackerGeomDet
Definition: TrackerGeomDet.h:6
BTLDeviceSim::bxTime_
const float bxTime_
Definition: BTLDeviceSim.h:42
GeomDetType.h
ProxyMTDTopology::specificTopology
virtual const PixelTopology & specificTopology() const
Definition: ProxyMTDTopology.h:90
GeomDet::topology
virtual const Topology & topology() const
Definition: GeomDet.cc:67
geant_units::operators
Definition: GeantUnits.h:18
MTDTopologyMode.h
edm::Handle< edm::PSimHitContainer >
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition: MessageLogger.h:122
MTDTopology::getMTDTopologyMode
int getMTDTopologyMode() const
Definition: MTDTopology.h:27
BTLDeviceSim::BTLDeviceSim
BTLDeviceSim(const edm::ParameterSet &pset, edm::ConsumesCollector iC)
Definition: BTLDeviceSim.cc:15
hit::x
double x
Definition: SiStripHitEffFromCalibTree.cc:89
geant_units::operators::convertGeVToMeV
constexpr NumType convertGeVToMeV(NumType gev)
Definition: GeantUnits.h:93
mtd_digitizer::kNumberOfBX
constexpr int kNumberOfBX
Definition: MTDDigitizerTypes.h:41
DetId
Definition: DetId.h:17
MTDGeometry::idToDet
const MTDGeomDet * idToDet(DetId) const override
Definition: MTDGeometry.cc:171
BTLDetId
Detector identifier class for the Barrel Timing Layer. The crystal count must start from 0,...
Definition: BTLDetId.h:18
BTLDeviceSim::LightCollSlopeR_
const float LightCollSlopeR_
Definition: BTLDeviceSim.h:46
BTLDeviceSim::topo_
const MTDTopology * topo_
Definition: BTLDeviceSim.h:40
BTLDetId::row
int row(unsigned nrows=16) const
Definition: BTLDetId.h:105
BTLDeviceSim::getHitsResponse
void getHitsResponse(const std::vector< std::tuple< int, uint32_t, float > > &hitRefs, const edm::Handle< edm::PSimHitContainer > &hits, mtd_digitizer::MTDSimHitDataAccumulator *simHitAccumulator, CLHEP::HepRandomEngine *hre)
Definition: BTLDeviceSim.cc:32
RectangularMTDTopology.h
Point3DBase< float, LocalTag >
MTDDetId::mtdSubDetector
int mtdSubDetector() const
Definition: MTDDetId.h:56
edm::ParameterSet
Definition: ParameterSet.h:47
MTDTopologyMode::crysLayoutFromTopoMode
BTLDetId::CrysLayout crysLayoutFromTopoMode(const int &topoMode)
Definition: MTDTopologyMode.h:19
MTDDetId.h
GeantUnits.h
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
mtd_digitizer::MTDCellId
Definition: MTDDigitizerTypes.h:26
edm::EventSetup
Definition: EventSetup.h:58
BTLDeviceSim::LightCollSlopeL_
const float LightCollSlopeL_
Definition: BTLDeviceSim.h:47
BTLDeviceSim::getEventSetup
void getEventSetup(const edm::EventSetup &evt)
Definition: BTLDeviceSim.cc:27
RectangularMTDTopology
Definition: RectangularMTDTopology.h:39
ProxyMTDTopology
Definition: ProxyMTDTopology.h:28
edm::EventSetup::getData
bool getData(T &iHolder) const
Definition: EventSetup.h:127
ProxyMTDTopology.h
DetId::rawId
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
DetId.h
MTDDetId
Detector identifier base class for the MIP Timing Layer.
Definition: MTDDetId.h:21
Exception
Definition: hltDiff.cc:245
mtd_digitizer::MTDSimHitDataAccumulator
std::unordered_map< MTDCellId, MTDCellInfo > MTDSimHitDataAccumulator
Definition: MTDDigitizerTypes.h:39
BTLDeviceSim::LightCollEff_
const float LightCollEff_
Definition: BTLDeviceSim.h:44
geant_units::operators::convertMmToCm
constexpr NumType convertMmToCm(NumType millimeters)
Definition: GeantUnits.h:63
ConsumesCollector.h
BTLDeviceSim::topoToken_
const edm::ESGetToken< MTDTopology, MTDTopologyRcd > topoToken_
Definition: BTLDeviceSim.h:38
PSimHit
Definition: PSimHit.h:15
DetId::Forward
Definition: DetId.h:30
mtd_digitizer::kInTimeBX
constexpr int kInTimeBX
Definition: MTDDigitizerTypes.h:42
BTLDetId::geographicalId
BTLDetId geographicalId(CrysLayout lay) const
Definition: BTLDetId.cc:171
DeDxTools::esConsumes
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
BTLDeviceSim::geomToken_
const edm::ESGetToken< MTDGeometry, MTDDigiGeometryRecord > geomToken_
Definition: BTLDeviceSim.h:37
TauDecayModes.dec
dec
Definition: TauDecayModes.py:142
edm::ConsumesCollector
Definition: ConsumesCollector.h:45
BTLDeviceSim::PDE_
const float PDE_
Definition: BTLDeviceSim.h:48
hit
Definition: SiStripHitEffFromCalibTree.cc:88
muonDTDigis_cfi.pset
pset
Definition: muonDTDigis_cfi.py:27
BTLDetId::column
int column(unsigned nrows=16) const
Definition: BTLDetId.h:110