CMS 3D CMS Logo

ETLElectronicsSim.cc
Go to the documentation of this file.
2 
5 #include "CLHEP/Random/RandGaussQ.h"
6 
7 using namespace mtd;
8 
10  : geomToken_(iC.esConsumes()),
11  geom_(nullptr),
12  debug_(pset.getUntrackedParameter<bool>("debug", false)),
13  bxTime_(pset.getParameter<double>("bxTime")),
14  integratedLum_(pset.getParameter<double>("IntegratedLuminosity")),
15  fluence_(pset.getParameter<std::string>("FluenceVsRadius")),
16  lgadGain_(pset.getParameter<std::string>("LGADGainVsFluence")),
17  timeRes2_(pset.getParameter<std::string>("TimeResolution2")),
18  adcNbits_(pset.getParameter<uint32_t>("adcNbits")),
19  tdcNbits_(pset.getParameter<uint32_t>("tdcNbits")),
20  adcSaturation_MIP_(pset.getParameter<double>("adcSaturation_MIP")),
21  adcLSB_MIP_(adcSaturation_MIP_ / std::pow(2., adcNbits_)),
22  adcBitSaturation_(std::pow(2, adcNbits_) - 1),
23  adcThreshold_MIP_(pset.getParameter<double>("adcThreshold_MIP")),
24  toaLSB_ns_(pset.getParameter<double>("toaLSB_ns")),
25  tdcBitSaturation_(std::pow(2, tdcNbits_) - 1) {}
26 
28 
31  CLHEP::HepRandomEngine* hre) const {
32  MTDSimHitData chargeColl, toa;
33 
34  std::vector<double> emptyV;
35  std::vector<double> radius(1);
36  std::vector<double> fluence(1);
37  std::vector<double> gain(1);
38 
39  for (MTDSimHitDataAccumulator::const_iterator it = input.begin(); it != input.end(); it++) {
40  chargeColl.fill(0.f);
41  toa.fill(0.f);
42 
43  ETLDetId detId = it->first.detid_;
44  DetId geoId = detId.geographicalId();
45  const MTDGeomDet* thedet = geom_->idToDet(geoId);
46  if (thedet == nullptr)
47  throw cms::Exception("EtlElectronicsSim") << "GeographicalID: " << std::hex << geoId.rawId() << " ("
48  << detId.rawId() << ") is invalid!" << std::dec << std::endl;
49  const PixelTopology& topo = static_cast<const PixelTopology&>(thedet->topology());
50 
51  Local3DPoint local_point(topo.localX(it->first.row_), topo.localY(it->first.column_), 0.);
52  const auto& global_point = thedet->toGlobal(local_point);
53 
54  for (size_t i = 0; i < it->second.hit_info[0].size(); i++) {
55  if ((it->second).hit_info[0][i] < adcThreshold_MIP_)
56  continue;
57 
58  // time of arrival
59  float finalToA = (it->second).hit_info[1][i];
60 
61  // calculate the LGAD gain as a function of the fluence at R = radius
62  radius[0] = global_point.perp();
63  fluence[0] = integratedLum_ * fluence_.evaluate(radius, emptyV);
64  gain[0] = lgadGain_.evaluate(fluence, emptyV);
65  if (gain[0] <= 0.)
66  throw cms::Exception("EtlElectronicsSim") << "Null or negative LGAD gain!" << std::endl;
67 
68  // Gaussian smearing of the time of arrival
69  double sigmaToA = sqrt(timeRes2_.evaluate(gain, emptyV));
70 
71  if (sigmaToA > 0.)
72  finalToA += CLHEP::RandGaussQ::shoot(hre, 0., sigmaToA);
73 
74  // fill the time and charge arrays
75  const unsigned int ibucket = std::floor(finalToA / bxTime_);
76  if ((i + ibucket) >= chargeColl.size())
77  continue;
78 
79  chargeColl[i + ibucket] += (it->second).hit_info[0][i];
80 
81  if (toa[i + ibucket] == 0. || (finalToA - ibucket * bxTime_) < toa[i + ibucket])
82  toa[i + ibucket] = finalToA - ibucket * bxTime_;
83  }
84 
85  // run the shaper to create a new data frame
86  ETLDataFrame rawDataFrame(it->first.detid_);
87  runTrivialShaper(rawDataFrame, chargeColl, toa, it->first.row_, it->first.column_);
88  updateOutput(output, rawDataFrame);
89  }
90 }
91 
93  const mtd::MTDSimHitData& chargeColl,
94  const mtd::MTDSimHitData& toa,
95  const uint8_t row,
96  const uint8_t col) const {
97  bool debug = debug_;
98 #ifdef EDM_ML_DEBUG
99  for (int it = 0; it < (int)(chargeColl.size()); it++)
100  debug |= (chargeColl[it] > adcThreshold_MIP_);
101 #endif
102 
103  if (debug)
104  edm::LogVerbatim("ETLElectronicsSim") << "[runTrivialShaper]" << std::endl;
105 
106  //set new ADCs
107  for (int it = 0; it < (int)(chargeColl.size()); it++) {
108  //brute force saturation, maybe could to better with an exponential like saturation
109  const uint32_t adc = std::min((uint32_t)std::floor(chargeColl[it] / adcLSB_MIP_), adcBitSaturation_);
110  const uint32_t tdc_time = std::min((uint32_t)std::floor(toa[it] / toaLSB_ns_), tdcBitSaturation_);
111  ETLSample newSample;
112  newSample.set(chargeColl[it] > adcThreshold_MIP_, false, tdc_time, adc, row, col);
113  dataFrame.setSample(it, newSample);
114 
115  if (debug)
116  edm::LogVerbatim("ETLElectronicsSim") << adc << " (" << chargeColl[it] << "/" << adcLSB_MIP_ << ") ";
117  }
118 
119  if (debug) {
120  std::ostringstream msg;
121  dataFrame.print(msg);
122  edm::LogVerbatim("ETLElectronicsSim") << msg.str() << std::endl;
123  }
124 }
125 
126 void ETLElectronicsSim::updateOutput(ETLDigiCollection& coll, const ETLDataFrame& rawDataFrame) const {
127  int itIdx(9);
128  if (rawDataFrame.size() <= itIdx + 2)
129  return;
130 
131  ETLDataFrame dataFrame(rawDataFrame.id());
132  dataFrame.resize(dfSIZE);
133  bool putInEvent(false);
134  for (int it = 0; it < dfSIZE; ++it) {
135  dataFrame.setSample(it, rawDataFrame[itIdx - 2 + it]);
136  if (it == 2)
137  putInEvent = rawDataFrame[itIdx - 2 + it].threshold();
138  }
139 
140  if (putInEvent) {
141  coll.push_back(dataFrame);
142  }
143 }
Log< level::Info, true > LogVerbatim
const edm::ESGetToken< MTDGeometry, MTDDigiGeometryRecord > geomToken_
static constexpr int dfSIZE
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
double evaluate(V const &iVariables, P const &iParameters) const
void updateOutput(ETLDigiCollection &coll, const ETLDataFrame &rawDataFrame) const
std::array< MTDSimData_t, nSamples > MTDSimHitData
virtual const Topology & topology() const
Definition: GeomDet.cc:67
const float integratedLum_
void push_back(T const &t)
std::unordered_map< MTDCellId, MTDCellInfo > MTDSimHitDataAccumulator
const float toaLSB_ns_
wrapper for a data word
Definition: ETLSample.h:13
ETLDetId geographicalId() const
Definition: ETLDetId.h:108
static std::string const input
Definition: EdmProvDump.cc:47
const reco::FormulaEvaluator fluence_
virtual float localX(float mpX) const =0
int size() const
total number of samples in the digi
Definition: FTLDataFrameT.h:46
void resize(size_t s)
allow to set size
Definition: FTLDataFrameT.h:51
ETLElectronicsSim(const edm::ParameterSet &pset, edm::ConsumesCollector iC)
T sqrt(T t)
Definition: SSEVec.h:19
const D & id() const
det id
Definition: FTLDataFrameT.h:31
const MTDGeomDet * idToDet(DetId) const override
Definition: MTDGeometry.cc:171
double f[11][100]
void setSample(int i, const S &sample)
Definition: FTLDataFrameT.h:58
bool getData(T &iHolder) const
Definition: EventSetup.h:122
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
Definition: DetId.h:17
const uint32_t adcBitSaturation_
#define debug
Definition: HDRShower.cc:19
virtual float localY(float mpY) const =0
void runTrivialShaper(ETLDataFrame &dataFrame, const mtd::MTDSimHitData &chargeColl, const mtd::MTDSimHitData &toa, const uint8_t row, const uint8_t column) const
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
const uint32_t tdcBitSaturation_
const reco::FormulaEvaluator lgadGain_
tuple msg
Definition: mps_check.py:285
const float adcLSB_MIP_
void set(bool thr, bool mode, uint16_t toa, uint16_t data, uint8_t row, uint8_t col)
Definition: ETLSample.h:48
Detector identifier class for the Endcap Timing Layer.
Definition: ETLDetId.h:15
col
Definition: cuy.py:1009
Readout digi for HGC.
Definition: FTLDataFrameT.h:14
void getEventSetup(const edm::EventSetup &evt)
const reco::FormulaEvaluator timeRes2_
void run(const mtd::MTDSimHitDataAccumulator &input, ETLDigiCollection &output, CLHEP::HepRandomEngine *hre) const
void print(std::ostream &out=std::cout)
Definition: FTLDataFrameT.h:62
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
const MTDGeometry * geom_
uint16_t *__restrict__ uint16_t const *__restrict__ adc
const float adcThreshold_MIP_