CMS 3D CMS Logo

MTDTimeCalib.cc
Go to the documentation of this file.
2 
5 
10 
12  : geom_(geom),
13  topo_(topo),
14  btlTimeOffset_(conf.getParameter<double>("BTLTimeOffset")),
15  etlTimeOffset_(conf.getParameter<double>("ETLTimeOffset")),
16  btlLightCollTime_(conf.getParameter<double>("BTLLightCollTime")),
17  btlLightCollSlope_(conf.getParameter<double>("BTLLightCollSlope")) {}
18 
19 float MTDTimeCalib::getTimeCalib(const MTDDetId& id) const {
20  if (id.subDetector() != MTDDetId::FastTime) {
21  throw cms::Exception("MTDTimeCalib") << "MTDDetId: " << std::hex << id.rawId() << " is invalid!" << std::dec
22  << std::endl;
23  }
24 
25  float time_calib = 0.;
26 
27  if (id.mtdSubDetector() == MTDDetId::BTL) {
28  time_calib += btlTimeOffset_;
29  BTLDetId hitId(id);
30  DetId geoId = hitId.geographicalId(
31  (BTLDetId::CrysLayout)topo_->getMTDTopologyMode()); //for BTL topology gives different layout id
32  const MTDGeomDet* thedet = geom_->idToDet(geoId);
33 
34  if (thedet == nullptr) {
35  throw cms::Exception("MTDTimeCalib") << "GeographicalID: " << std::hex << geoId.rawId() << " (" << id.rawId()
36  << ") is invalid!" << std::dec << std::endl;
37  }
38  const ProxyMTDTopology& topoproxy = static_cast<const ProxyMTDTopology&>(thedet->topology());
39  const RectangularMTDTopology& topo = static_cast<const RectangularMTDTopology&>(topoproxy.specificTopology());
40 
42  time_calib -= btlLightCollTime_; //simply remove the offset introduced at sim level
45  //for bars in phi
46  time_calib -= 0.5 * topo.pitch().first * btlLightCollSlope_; //time offset for bar time is L/2v
48  //for bars in z
49  time_calib -= 0.5 * topo.pitch().second * btlLightCollSlope_; //time offset for bar time is L/2v
50  }
51  } else if (id.mtdSubDetector() == MTDDetId::ETL) {
52  time_calib += etlTimeOffset_;
53  } else {
54  throw cms::Exception("MTDTimeCalib") << "MTDDetId: " << std::hex << id.rawId() << " is invalid!" << std::dec
55  << std::endl;
56  }
57 
58  return time_calib;
59 }
60 
62 
63 //--- Now use the Framework macros to set it all up:
float etlTimeOffset_
Definition: MTDTimeCalib.h:24
CrysLayout
Definition: BTLDetId.h:70
virtual const Topology & topology() const
Definition: GeomDet.cc:67
float btlLightCollSlope_
Definition: MTDTimeCalib.h:28
int getMTDTopologyMode() const
Definition: MTDTopology.h:73
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
float btlLightCollTime_
Definition: MTDTimeCalib.h:27
std::pair< float, float > pitch() const override
Detector identifier base class for the MIP Timing Layer.
Definition: MTDDetId.h:21
const MTDGeomDet * idToDet(DetId) const override
Definition: MTDGeometry.cc:160
BTLDetId geographicalId(CrysLayout lay) const
Definition: BTLDetId.cc:163
MTDTimeCalib(edm::ParameterSet const &conf, const MTDGeometry *geom, const MTDTopology *topo)
Definition: MTDTimeCalib.cc:11
virtual const PixelTopology & specificTopology() const
const MTDTopology * topo_
Definition: MTDTimeCalib.h:22
float btlTimeOffset_
Definition: MTDTimeCalib.h:23
#define TYPELOOKUP_DATA_REG(_dataclass_)
Definition: typelookup.h:102
Definition: DetId.h:17
float getTimeCalib(const MTDDetId &id) const
Definition: MTDTimeCalib.cc:19
const MTDGeometry * geom_
Definition: MTDTimeCalib.h:21
Detector identifier class for the Barrel Timing Layer. The crystal count must start from 0...
Definition: BTLDetId.h:18