CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Attributes | Private Attributes
ETLElectronicsSim Class Reference

#include <ETLElectronicsSim.h>

Public Member Functions

 ETLElectronicsSim (const edm::ParameterSet &pset, edm::ConsumesCollector iC)
 
void getEvent (const edm::Event &evt)
 
void getEventSetup (const edm::EventSetup &evt)
 
void run (const mtd::MTDSimHitDataAccumulator &input, ETLDigiCollection &output, CLHEP::HepRandomEngine *hre) const
 
void runTrivialShaper (ETLDataFrame &dataFrame, const mtd::MTDSimHitData &chargeColl, const mtd::MTDSimHitData &toa1, const mtd::MTDSimHitData &toa2, const uint8_t row, const uint8_t column) const
 
void updateOutput (ETLDigiCollection &coll, const ETLDataFrame &rawDataFrame) const
 

Static Public Attributes

static constexpr int dfSIZE = 5
 

Private Attributes

const uint32_t adcBitSaturation_
 
const float adcLSB_MIP_
 
const uint32_t adcNbits_
 
const float adcSaturation_MIP_
 
const float adcThreshold_MIP_
 
const float bxTime_
 
const bool debug_
 
const ETLPulseShape etlPulseShape_
 
const reco::FormulaEvaluator formulaLandauNoise_
 
const MTDGeometrygeom_
 
const edm::ESGetToken< MTDGeometry, MTDDigiGeometryRecordgeomToken_
 
const float integratedLum_
 
const float iThreshold_MIP_
 
const float noiseLevel_
 
const float referenceChargeColl_
 
const float sigmaDistorsion_
 
const float sigmaTDC_
 
const uint32_t tdcBitSaturation_
 
const uint32_t tdcNbits_
 
const float toaLSB_ns_
 

Detailed Description

Definition at line 27 of file ETLElectronicsSim.h.

Constructor & Destructor Documentation

◆ ETLElectronicsSim()

ETLElectronicsSim::ETLElectronicsSim ( const edm::ParameterSet pset,
edm::ConsumesCollector  iC 
)

Definition at line 9 of file ETLElectronicsSim.cc.

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  adcNbits_(pset.getParameter<uint32_t>("adcNbits")),
16  tdcNbits_(pset.getParameter<uint32_t>("tdcNbits")),
17  adcSaturation_MIP_(pset.getParameter<double>("adcSaturation_MIP")),
20  adcThreshold_MIP_(pset.getParameter<double>("adcThreshold_MIP")),
21  iThreshold_MIP_(pset.getParameter<double>("iThreshold_MIP")),
22  toaLSB_ns_(pset.getParameter<double>("toaLSB_ns")),
24  referenceChargeColl_(pset.getParameter<double>("referenceChargeColl")),
25  noiseLevel_(pset.getParameter<double>("noiseLevel")),
26  sigmaDistorsion_(pset.getParameter<double>("sigmaDistorsion")),
27  sigmaTDC_(pset.getParameter<double>("sigmaTDC")),
28  formulaLandauNoise_(pset.getParameter<std::string>("formulaLandauNoise")) {}
const edm::ESGetToken< MTDGeometry, MTDDigiGeometryRecord > geomToken_
const reco::FormulaEvaluator formulaLandauNoise_
const float integratedLum_
const float toaLSB_ns_
const uint32_t tdcNbits_
const float referenceChargeColl_
const uint32_t adcNbits_
const float sigmaDistorsion_
const uint32_t adcBitSaturation_
const float iThreshold_MIP_
const uint32_t tdcBitSaturation_
const float adcLSB_MIP_
const float noiseLevel_
const float adcSaturation_MIP_
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
const MTDGeometry * geom_
const float adcThreshold_MIP_

Member Function Documentation

◆ getEvent()

void ETLElectronicsSim::getEvent ( const edm::Event evt)
inline

Definition at line 31 of file ETLElectronicsSim.h.

31 {}

◆ getEventSetup()

void ETLElectronicsSim::getEventSetup ( const edm::EventSetup evt)

Definition at line 30 of file ETLElectronicsSim.cc.

References geom_, geomToken_, and edm::EventSetup::getData().

30 { geom_ = &evs.getData(geomToken_); }
const edm::ESGetToken< MTDGeometry, MTDDigiGeometryRecord > geomToken_
const MTDGeometry * geom_

◆ run()

void ETLElectronicsSim::run ( const mtd::MTDSimHitDataAccumulator input,
ETLDigiCollection output,
CLHEP::HepRandomEngine *  hre 
) const

Definition at line 32 of file ETLElectronicsSim.cc.

References adcThreshold_MIP_, bxTime_, TauDecayModes::dec, hcalRecHitTable_cff::detId, etlPulseShape_, reco::FormulaEvaluator::evaluate(), Exception, f, MTDShapeBase::fallTime(), formulaLandauNoise_, geom_, mps_fire::i, MTDGeometry::idToDet(), input, ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, iThreshold_MIP_, PixelTopology::localX(), PixelTopology::localY(), MTDShapeBase::maximum(), noiseLevel_, CosmicsPD_Skims::radius, DetId::rawId(), referenceChargeColl_, runTrivialShaper(), mtdDigitizer_cfi::sigmaDistorsion, sigmaDistorsion_, mtdDigitizer_cfi::sigmaTDC, sigmaTDC_, mathSSE::sqrt(), MTDShapeBase::timeAtThr(), MTDShapeBase::timeOfMax(), GeomDet::topology(), compareTotals::tot, and updateOutput().

34  {
35  MTDSimHitData chargeColl, toa1, toa2, tot;
36 
37  std::vector<double> emptyV;
38  std::vector<double> radius(1);
39  std::vector<double> fluence(1);
40  std::vector<double> chOverMPV(1);
41 
42  for (MTDSimHitDataAccumulator::const_iterator it = input.begin(); it != input.end(); it++) {
43  chargeColl.fill(0.f);
44  toa1.fill(0.f);
45  toa2.fill(0.f);
46  tot.fill(0.f);
47 
48  ETLDetId detId = it->first.detid_;
49  DetId geoId = detId.geographicalId();
50  const MTDGeomDet* thedet = geom_->idToDet(geoId);
51  if (thedet == nullptr)
52  throw cms::Exception("EtlElectronicsSim") << "GeographicalID: " << std::hex << geoId.rawId() << " ("
53  << detId.rawId() << ") is invalid!" << std::dec << std::endl;
54 
55  const PixelTopology& topo = static_cast<const PixelTopology&>(thedet->topology());
56  Local3DPoint local_point(topo.localX(it->first.row_), topo.localY(it->first.column_), 0.);
57 
58  //Here we are sampling over all the buckets. However for the ETL there should be only one bucket at i = 0
59  for (size_t i = 0; i < it->second.hit_info[0].size(); i++) {
60  if ((it->second).hit_info[0][i] == 0) {
61  continue;
62  }
63  if ((it->second).hit_info[0][i] < adcThreshold_MIP_)
64  continue;
65  // time of arrival
66  double finalToA = (it->second).hit_info[1][i];
67  double finalToC = (it->second).hit_info[1][i];
68 
69  // fill the time and charge arrays
70  const unsigned int ibucket = std::floor(finalToA / bxTime_);
71  if ((i + ibucket) >= chargeColl.size())
72  continue;
73 
74  //Calculate the jitter
75  double SignalToNoise =
76  etlPulseShape_.maximum() * ((it->second).hit_info[0][i] / referenceChargeColl_) / noiseLevel_;
77  double sigmaJitter1 = etlPulseShape_.timeOfMax() / SignalToNoise;
78  double sigmaJitter2 = (etlPulseShape_.fallTime() - etlPulseShape_.timeOfMax()) / SignalToNoise;
79  //Calculate the distorsion
81  //Calculate the TDC
82  double sigmaTDC = sigmaTDC_;
83  //Calculate the Landau Noise
84  chOverMPV[0] = (it->second).hit_info[0][i] / (it->second).hit_info[2][i];
85  double sigmaLN = formulaLandauNoise_.evaluate(chOverMPV, emptyV);
86  double sigmaToA = sqrt(sigmaJitter1 * sigmaJitter1 + sigmaDistorsion * sigmaDistorsion + sigmaTDC * sigmaTDC +
87  sigmaLN * sigmaLN);
88  double sigmaToC = sqrt(sigmaJitter2 * sigmaJitter2 + sigmaDistorsion * sigmaDistorsion + sigmaTDC * sigmaTDC +
89  sigmaLN * sigmaLN);
90  double smearing1 = 0.0;
91  double smearing2 = 0.0;
92 
93  if (sigmaToA > 0. && sigmaToC > 0.) {
94  smearing1 = CLHEP::RandGaussQ::shoot(hre, 0., sigmaToA);
95  smearing2 = CLHEP::RandGaussQ::shoot(hre, 0., sigmaToC);
96  }
97 
98  finalToA += smearing1;
99  finalToC += smearing1 + smearing2;
100 
101  std::array<float, 3> times = etlPulseShape_.timeAtThr(
102  (it->second).hit_info[0][i] / referenceChargeColl_, iThreshold_MIP_, iThreshold_MIP_);
103 
104  //The signal is below the threshold
105  if (times[0] == 0 && times[1] == 0 && times[2] == 0) {
106  continue;
107  }
108  //The signal is considered to be below the threshold
109  finalToA += times[0];
110  finalToC += times[2];
111  if (finalToA >= finalToC)
112  continue;
113  chargeColl[i + ibucket] += (it->second).hit_info[0][i];
114 
115  if (toa1[i + ibucket] == 0. || (finalToA - ibucket * bxTime_) < toa1[i + ibucket]) {
116  toa1[i + ibucket] = finalToA - ibucket * bxTime_;
117  toa2[i + ibucket] = finalToC - ibucket * bxTime_;
118  }
119 
120  tot[i + ibucket] = finalToC - finalToA;
121  }
122  // Run the shaper to create a new data frame
123  ETLDataFrame rawDataFrame(it->first.detid_);
124  runTrivialShaper(rawDataFrame, chargeColl, toa1, tot, it->first.row_, it->first.column_);
125  updateOutput(output, rawDataFrame);
126  }
127 }
std::array< float, 3 > timeAtThr(const float scale, const float threshold1, const float threshold2) const
Definition: MTDShapeBase.cc:13
const reco::FormulaEvaluator formulaLandauNoise_
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
float fallTime() const
Definition: MTDShapeBase.cc:79
const ETLPulseShape etlPulseShape_
const float referenceChargeColl_
static std::string const input
Definition: EdmProvDump.cc:50
virtual float localX(float mpX) const =0
T sqrt(T t)
Definition: SSEVec.h:19
const MTDGeomDet * idToDet(DetId) const override
Definition: MTDGeometry.cc:171
double f[11][100]
float maximum() const
Definition: MTDShapeBase.cc:77
const float sigmaDistorsion_
Definition: DetId.h:17
double timeOfMax() const
Definition: MTDShapeBase.cc:75
virtual float localY(float mpY) const =0
const float iThreshold_MIP_
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
void runTrivialShaper(ETLDataFrame &dataFrame, const mtd::MTDSimHitData &chargeColl, const mtd::MTDSimHitData &toa1, const mtd::MTDSimHitData &toa2, const uint8_t row, const uint8_t column) const
Detector identifier class for the Endcap Timing Layer.
Definition: ETLDetId.h:16
const float noiseLevel_
Readout digi for HGC.
Definition: FTLDataFrameT.h:14
Definition: output.py:1
const MTDGeometry * geom_
const float adcThreshold_MIP_

◆ runTrivialShaper()

void ETLElectronicsSim::runTrivialShaper ( ETLDataFrame dataFrame,
const mtd::MTDSimHitData chargeColl,
const mtd::MTDSimHitData toa1,
const mtd::MTDSimHitData toa2,
const uint8_t  row,
const uint8_t  column 
) const

Definition at line 129 of file ETLElectronicsSim.cc.

References gpuClustering::adc, adcBitSaturation_, adcLSB_MIP_, adcThreshold_MIP_, cuy::col, debug, debug_, createfilelist::int, ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, SiStripPI::min, mps_check::msg, FTLDataFrameT< D, S, DECODE >::print(), ETLSample::set(), FTLDataFrameT< D, S, DECODE >::setSample(), tdcBitSaturation_, toaLSB_ns_, and compareTotals::tot.

Referenced by run().

134  {
135  bool debug = debug_;
136 #ifdef EDM_ML_DEBUG
137  for (int it = 0; it < (int)(chargeColl.size()); it++)
138  debug |= (chargeColl[it] > adcThreshold_MIP_);
139 #endif
140 
141  if (debug)
142  edm::LogVerbatim("ETLElectronicsSim") << "[runTrivialShaper]" << std::endl;
143 
144  //set new ADCs. Notice that we are only interested in the first element of the array for the ETL
145  for (int it = 0; it < (int)(chargeColl.size()); it++) {
146  //brute force saturation, maybe could to better with an exponential like saturation
147  const uint32_t adc = std::min((uint32_t)std::floor(chargeColl[it] / adcLSB_MIP_), adcBitSaturation_);
148  const uint32_t tdc_time1 = std::min((uint32_t)std::floor(toa[it] / toaLSB_ns_), tdcBitSaturation_);
149  const uint32_t tdc_time2 = std::min((uint32_t)std::floor(tot[it] / toaLSB_ns_), tdcBitSaturation_);
150  //If time over threshold is 0 the event is assumed to not pass the threshold
151  bool thres = true;
152  if (tdc_time2 == 0 || chargeColl[it] < adcThreshold_MIP_)
153  thres = false;
154 
155  ETLSample newSample;
156  newSample.set(thres, false, tdc_time1, tdc_time2, adc, row, col);
157 
158  //ETLSample newSample;
159  dataFrame.setSample(it, newSample);
160 
161  if (debug)
162  edm::LogVerbatim("ETLElectronicsSim") << adc << " (" << chargeColl[it] << "/" << adcLSB_MIP_ << ") ";
163  }
164 
165  if (debug) {
166  std::ostringstream msg;
167  dataFrame.print(msg);
168  edm::LogVerbatim("ETLElectronicsSim") << msg.str() << std::endl;
169  }
170 }
Log< level::Info, true > LogVerbatim
const float toaLSB_ns_
wrapper for a data word
Definition: ETLSample.h:13
void set(bool thr, bool mode, uint16_t toa, uint16_t tot, uint16_t data, uint8_t row, uint8_t col)
Definition: ETLSample.h:53
void setSample(int i, const S &sample)
Definition: FTLDataFrameT.h:58
const uint32_t adcBitSaturation_
#define debug
Definition: HDRShower.cc:19
const uint32_t tdcBitSaturation_
tuple msg
Definition: mps_check.py:286
const float adcLSB_MIP_
col
Definition: cuy.py:1009
void print(std::ostream &out=std::cout)
Definition: FTLDataFrameT.h:62
uint16_t *__restrict__ uint16_t const *__restrict__ adc
const float adcThreshold_MIP_

◆ updateOutput()

void ETLElectronicsSim::updateOutput ( ETLDigiCollection coll,
const ETLDataFrame rawDataFrame 
) const

Definition at line 172 of file ETLElectronicsSim.cc.

References dfSIZE, FTLDataFrameT< D, S, DECODE >::id(), ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, edm::SortedCollection< T, SORT >::push_back(), FTLDataFrameT< D, S, DECODE >::resize(), and FTLDataFrameT< D, S, DECODE >::size().

Referenced by run().

172  {
173  int itIdx(9);
174  if (rawDataFrame.size() <= itIdx + 2)
175  return;
176 
177  ETLDataFrame dataFrame(rawDataFrame.id());
178  dataFrame.resize(dfSIZE);
179  bool putInEvent(false);
180  for (int it = 0; it < dfSIZE; ++it) {
181  dataFrame.setSample(it, rawDataFrame[itIdx - 2 + it]);
182  if (it == 2)
183  putInEvent = rawDataFrame[itIdx - 2 + it].threshold();
184  }
185 
186  if (putInEvent) {
187  coll.push_back(dataFrame);
188  }
189 }
static constexpr int dfSIZE
void push_back(T const &t)
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
const D & id() const
det id
Definition: FTLDataFrameT.h:31
Readout digi for HGC.
Definition: FTLDataFrameT.h:14

Member Data Documentation

◆ adcBitSaturation_

const uint32_t ETLElectronicsSim::adcBitSaturation_
private

Definition at line 64 of file ETLElectronicsSim.h.

Referenced by runTrivialShaper().

◆ adcLSB_MIP_

const float ETLElectronicsSim::adcLSB_MIP_
private

Definition at line 63 of file ETLElectronicsSim.h.

Referenced by runTrivialShaper().

◆ adcNbits_

const uint32_t ETLElectronicsSim::adcNbits_
private

Definition at line 59 of file ETLElectronicsSim.h.

◆ adcSaturation_MIP_

const float ETLElectronicsSim::adcSaturation_MIP_
private

Definition at line 62 of file ETLElectronicsSim.h.

◆ adcThreshold_MIP_

const float ETLElectronicsSim::adcThreshold_MIP_
private

Definition at line 65 of file ETLElectronicsSim.h.

Referenced by run(), and runTrivialShaper().

◆ bxTime_

const float ETLElectronicsSim::bxTime_
private

Definition at line 53 of file ETLElectronicsSim.h.

Referenced by run().

◆ debug_

const bool ETLElectronicsSim::debug_
private

Definition at line 52 of file ETLElectronicsSim.h.

Referenced by runTrivialShaper().

◆ dfSIZE

constexpr int ETLElectronicsSim::dfSIZE = 5
static

Definition at line 46 of file ETLElectronicsSim.h.

Referenced by updateOutput().

◆ etlPulseShape_

const ETLPulseShape ETLElectronicsSim::etlPulseShape_
private

Definition at line 56 of file ETLElectronicsSim.h.

Referenced by run().

◆ formulaLandauNoise_

const reco::FormulaEvaluator ETLElectronicsSim::formulaLandauNoise_
private

Definition at line 73 of file ETLElectronicsSim.h.

Referenced by run().

◆ geom_

const MTDGeometry* ETLElectronicsSim::geom_
private

Definition at line 50 of file ETLElectronicsSim.h.

Referenced by getEventSetup(), and run().

◆ geomToken_

const edm::ESGetToken<MTDGeometry, MTDDigiGeometryRecord> ETLElectronicsSim::geomToken_
private

Definition at line 49 of file ETLElectronicsSim.h.

Referenced by getEventSetup().

◆ integratedLum_

const float ETLElectronicsSim::integratedLum_
private

Definition at line 54 of file ETLElectronicsSim.h.

◆ iThreshold_MIP_

const float ETLElectronicsSim::iThreshold_MIP_
private

Definition at line 66 of file ETLElectronicsSim.h.

Referenced by run().

◆ noiseLevel_

const float ETLElectronicsSim::noiseLevel_
private

Definition at line 70 of file ETLElectronicsSim.h.

Referenced by run().

◆ referenceChargeColl_

const float ETLElectronicsSim::referenceChargeColl_
private

Definition at line 69 of file ETLElectronicsSim.h.

Referenced by run().

◆ sigmaDistorsion_

const float ETLElectronicsSim::sigmaDistorsion_
private

Definition at line 71 of file ETLElectronicsSim.h.

Referenced by run().

◆ sigmaTDC_

const float ETLElectronicsSim::sigmaTDC_
private

Definition at line 72 of file ETLElectronicsSim.h.

Referenced by run().

◆ tdcBitSaturation_

const uint32_t ETLElectronicsSim::tdcBitSaturation_
private

Definition at line 68 of file ETLElectronicsSim.h.

Referenced by runTrivialShaper().

◆ tdcNbits_

const uint32_t ETLElectronicsSim::tdcNbits_
private

Definition at line 59 of file ETLElectronicsSim.h.

◆ toaLSB_ns_

const float ETLElectronicsSim::toaLSB_ns_
private

Definition at line 67 of file ETLElectronicsSim.h.

Referenced by runTrivialShaper().