CMS 3D CMS Logo

HGCFEElectronics.cc
Go to the documentation of this file.
1 
5 
6 #include "vdt/vdtMath.h"
7 
8 using namespace hgc_digi;
9 
10 //
11 template <class DFr>
13  : fwVersion_{ps.getParameter<uint32_t>("fwVersion")},
14  adcPulse_{},
15  pulseAvgT_{},
16  tdcForToAOnset_fC_{},
17  adcSaturation_fC_{-1.0},
18  adcLSB_fC_{},
19  tdcLSB_fC_{},
20  tdcSaturation_fC_{-1.0},
21  adcThreshold_fC_{},
22  tdcOnset_fC_{},
23  toaLSB_ns_{},
24  tdcResolutionInNs_{1e-9}, // set time resolution very small by default
25  targetMIPvalue_ADC_{},
26  jitterNoise_ns_{},
27  jitterConstant_ns_{},
28  eventTimeOffset_ns_{{0.02, 0.02, 0.02}},
29  noise_fC_{},
30  toaMode_(WEIGHTEDBYE) {
31  edm::LogVerbatim("HGCFE") << "[HGCFEElectronics] running with version " << fwVersion_ << std::endl;
32  if (ps.exists("adcPulse")) {
33  auto temp = ps.getParameter<std::vector<double> >("adcPulse");
34  for (unsigned i = 0; i < temp.size(); ++i) {
35  adcPulse_[i] = (float)temp[i];
36  }
37  // normalize adc pulse
38  for (unsigned i = 0; i < adcPulse_.size(); ++i) {
39  adcPulse_[i] = adcPulse_[i] / adcPulse_[2];
40  }
41  temp = ps.getParameter<std::vector<double> >("pulseAvgT");
42  for (unsigned i = 0; i < temp.size(); ++i) {
43  pulseAvgT_[i] = (float)temp[i];
44  }
45  }
46  if (ps.exists("adcNbits")) {
47  uint32_t adcNbits = ps.getParameter<uint32_t>("adcNbits");
48  adcSaturation_fC_ = ps.getParameter<double>("adcSaturation_fC");
49  adcLSB_fC_ = adcSaturation_fC_ / pow(2., adcNbits);
50  edm::LogVerbatim("HGCFE") << "[HGCFEElectronics] " << adcNbits << " bit ADC defined"
51  << " with LSB=" << adcLSB_fC_ << " saturation to occur @ " << adcSaturation_fC_
52  << std::endl;
53  }
54 
55  if (ps.exists("tdcNbits")) {
56  tdcNbits_ = ps.getParameter<uint32_t>("tdcNbits");
57  setTDCfsc(ps.getParameter<double>("tdcSaturation_fC"));
58  edm::LogVerbatim("HGCFE") << "[HGCFEElectronics] " << tdcNbits_ << " bit TDC defined with LSB=" << tdcLSB_fC_
59  << " saturation to occur @ " << tdcSaturation_fC_
60  << " (NB lowered by 1 part in a million)" << std::endl;
61  }
62  if (ps.exists("targetMIPvalue_ADC"))
63  targetMIPvalue_ADC_ = ps.getParameter<uint32_t>("targetMIPvalue_ADC");
64  if (ps.exists("adcThreshold_fC"))
65  adcThreshold_fC_ = ps.getParameter<double>("adcThreshold_fC");
66  if (ps.exists("tdcOnset_fC"))
67  tdcOnset_fC_ = ps.getParameter<double>("tdcOnset_fC");
68  if (ps.exists("tdcForToAOnset_fC")) {
69  auto temp = ps.getParameter<std::vector<double> >("tdcForToAOnset_fC");
70  if (temp.size() == tdcForToAOnset_fC_.size()) {
71  std::copy_n(temp.begin(), temp.size(), tdcForToAOnset_fC_.begin());
72  } else {
73  throw cms::Exception("BadConfiguration") << " HGCFEElectronics wrong size for ToA thresholds ";
74  }
75  }
76  if (ps.exists("toaLSB_ns"))
77  toaLSB_ns_ = ps.getParameter<double>("toaLSB_ns");
78  if (ps.exists("tdcChargeDrainParameterisation")) {
79  for (auto val : ps.getParameter<std::vector<double> >("tdcChargeDrainParameterisation")) {
80  tdcChargeDrainParameterisation_.push_back((float)val);
81  }
82  }
83  if (ps.exists("tdcResolutionInPs"))
84  tdcResolutionInNs_ = ps.getParameter<double>("tdcResolutionInPs") * 1e-3; // convert to ns
85  if (ps.exists("toaMode"))
86  toaMode_ = ps.getParameter<uint32_t>("toaMode");
87 
88  if (ps.exists("jitterNoise_ns")) {
89  auto temp = ps.getParameter<std::vector<double> >("jitterNoise_ns");
90  if (temp.size() == jitterNoise_ns_.size()) {
91  std::copy_n(temp.begin(), temp.size(), jitterNoise_ns_.begin());
92  } else {
93  throw cms::Exception("BadConfiguration") << " HGCFEElectronics wrong size for ToA jitterNoise ";
94  }
95  }
96  if (ps.exists("jitterConstant_ns")) {
97  auto temp = ps.getParameter<std::vector<double> >("jitterConstant_ns");
98  if (temp.size() == jitterConstant_ns_.size()) {
99  std::copy_n(temp.begin(), temp.size(), jitterConstant_ns_.begin());
100  } else {
101  throw cms::Exception("BadConfiguration") << " HGCFEElectronics wrong size for ToA jitterConstant ";
102  }
103  }
104 }
105 
106 //
107 template <class DFr>
109  DFr& dataFrame, HGCSimHitData& chargeColl, uint32_t thrADC, float lsbADC, uint32_t gainIdx, float maxADC) {
110  bool debug(false);
111 
112 #ifdef EDM_ML_DEBUG
113  for (int it = 0; it < (int)(chargeColl.size()); it++)
114  debug |= (chargeColl[it] > adcThreshold_fC_);
115 #endif
116 
117  if (debug)
118  edm::LogVerbatim("HGCFE") << "[runTrivialShaper]" << std::endl;
119 
120  if (lsbADC < 0)
121  lsbADC = adcLSB_fC_;
122  if (maxADC < 0)
123  // lower adcSaturation_fC_ by one part in a million
124  // to ensure largest charge converted in bits is 0xfff==4095, not 0x1000
125  // no effect on charges loewer than; no impact on cpu time, only done once
126  maxADC = adcSaturation_fC_ * (1 - 1e-6);
127  for (int it = 0; it < (int)(chargeColl.size()); it++) {
128  //brute force saturation, maybe could to better with an exponential like saturation
129  const uint32_t adc = std::floor(std::min(chargeColl[it], maxADC) / lsbADC);
130  HGCSample newSample;
131  newSample.set(adc > thrADC, false, gainIdx, 0, adc);
132  dataFrame.setSample(it, newSample);
133 
134  if (debug)
135  edm::LogVerbatim("HGCFE") << adc << " (" << chargeColl[it] << "/" << adcLSB_fC_ << ") ";
136  }
137 
138  if (debug) {
139  std::ostringstream msg;
140  dataFrame.print(msg);
141  edm::LogVerbatim("HGCFE") << msg.str() << std::endl;
142  }
143 }
144 
145 //
146 template <class DFr>
148  HGCSimHitData& chargeColl,
149  uint32_t thrADC,
150  float lsbADC,
151  uint32_t gainIdx,
152  float maxADC,
154  //convolute with pulse shape to compute new ADCs
155  newCharge.fill(0.f);
156  bool debug(false);
157  for (int it = 0; it < (int)(chargeColl.size()); it++) {
158  const float charge(chargeColl[it]);
159  if (charge == 0.f)
160  continue;
161 
162 #ifdef EDM_ML_DEBUG
163  debug |= (charge > adcThreshold_fC_);
164 #endif
165 
166  if (debug)
167  edm::LogVerbatim("HGCFE") << "\t Redistributing SARS ADC" << charge << " @ " << it;
168 
169  for (int ipulse = -2; ipulse < (int)(adcPulse.size()) - 2; ipulse++) {
170  if (it + ipulse < 0)
171  continue;
172  if (it + ipulse >= (int)(dataFrame.size()))
173  continue;
174  const float chargeLeak = charge * adcPulse[(ipulse + 2)];
175  newCharge[it + ipulse] += chargeLeak;
176 
177  if (debug)
178  edm::LogVerbatim("HGCFE") << " | " << it + ipulse << " " << chargeLeak;
179  }
180 
181  if (debug)
182  edm::LogVerbatim("HGCFE") << std::endl;
183  }
184 
185  for (int it = 0; it < (int)(newCharge.size()); it++) {
186  //brute force saturation, maybe could to better with an exponential like saturation
187  const uint32_t adc = std::floor(std::min(newCharge[it], maxADC) / lsbADC);
188  HGCSample newSample;
189  newSample.set(adc > thrADC, false, gainIdx, 0, adc);
190  dataFrame.setSample(it, newSample);
191 
192  if (debug)
193  edm::LogVerbatim("HGCFE") << adc << " (" << std::min(newCharge[it], maxADC) << "/" << lsbADC << " ) ";
194  }
195 
196  if (debug) {
197  std::ostringstream msg;
198  dataFrame.print(msg);
199  edm::LogVerbatim("HGCFE") << msg.str() << std::endl;
200  }
201 }
202 
203 //
204 template <class DFr>
206  HGCSimHitData& chargeColl,
207  HGCSimHitData& toaColl,
208  CLHEP::HepRandomEngine* engine,
209  uint32_t thrADC,
210  float lsbADC,
211  uint32_t gainIdx,
212  float maxADC,
213  int thickness,
214  float tdcOnsetAuto,
215  float noiseWidth,
217  busyFlags.fill(false);
218  totFlags.fill(false);
219  toaFlags.fill(false);
220  newCharge.fill(0.f);
221  toaFromToT.fill(0.f);
222 
223 #ifdef EDM_ML_DEBUG
224  constexpr bool debug_state(true);
225 #else
226  constexpr bool debug_state(false);
227 #endif
228 
229  bool debug = debug_state;
230 
231  float timeToA = 0.f;
232 
233  //configure the ADC <-> TDC transition depending on the value passed as argument
234  double tdcOnset(maxADC);
235  if (maxADC < 0) {
236  maxADC = adcSaturation_fC_;
237  tdcOnset = tdcOnset_fC_;
238  }
239 
240  //configure the ADC LSB depending on the value passed as argument
241  if (lsbADC < 0)
242  lsbADC = adcLSB_fC_;
243 
244  //first look at time
245  //for pileup look only at intime signals
246  //ToA is in central BX if fired -- std::floor(BX/25.)+9;
247  int fireBX = 9;
248  //noise fluctuation on charge is added after ToA computation
249  //do not recheck the ToA firing threshold tdcForToAOnset_fC_[thickness-1] not to bias the efficiency
250  //to be done properly with realistic ToA shaper and jitter for the moment accounted in the smearing
251  if (toaColl[fireBX] != 0.f) {
252  timeToA = toaColl[fireBX];
253  float sensor_noise = noiseWidth <= 0 ? noise_fC_[thickness - 1] : noiseWidth;
254  float noise = jitterNoise_ns_[thickness - 1] * sensor_noise;
255  float jitter = chargeColl[fireBX] == 0 ? 0 : (noise / chargeColl[fireBX]);
256  if (jitter != 0)
257  timeToA = CLHEP::RandGaussQ::shoot(engine, timeToA, jitter);
258  else if (tdcResolutionInNs_ != 0)
259  timeToA = CLHEP::RandGaussQ::shoot(engine, timeToA, tdcResolutionInNs_);
260  timeToA += eventTimeOffset_ns_[thickness - 1];
261  if (timeToA >= 0.f && timeToA <= 25.f)
262  toaFlags[fireBX] = true;
263  }
264 
265  //now look at charge
266  //first identify bunches which will trigger ToT
267  //if(debug_state) edm::LogVerbatim("HGCFE") << "[runShaperWithToT]" << std::endl;
268  for (int it = 0; it < (int)(chargeColl.size()); ++it) {
269  debug = debug_state;
270  //if already flagged as busy it can't be re-used to trigger the ToT
271  if (busyFlags[it])
272  continue;
273 
274  if (tdcOnsetAuto < 0) {
275  tdcOnsetAuto = tdcOnset_fC_;
276  }
277  //if below TDC onset will be handled by SARS ADC later
278  float charge = chargeColl[it];
279  if (charge < tdcOnset) {
280  debug = false;
281  continue;
282  }
283 
284  //raise TDC mode for charge computation
285  //ToA anyway fired independently will be sorted out with realistic ToA dedicated shaper
286  float toa = timeToA;
287  totFlags[it] = true;
288 
289  if (debug)
290  edm::LogVerbatim("HGCFE") << "\t q=" << charge << " fC with <toa>=" << toa << " ns, triggers ToT @ " << it
291  << std::endl;
292 
293  //compute total charge to be integrated and integration time
294  //needs a loop as ToT will last as long as there is charge to dissipate
295  int busyBxs(0);
296  float totalCharge(charge), finalToA(toa), integTime(0);
297  while (true) {
298  //compute integration time in ns and # bunches
299  //float newIntegTime(0);
300  int poffset = 0;
301  float charge_offset = 0.f;
302  const float charge_kfC(totalCharge * 1e-3);
303  if (charge_kfC < tdcChargeDrainParameterisation_[3]) {
304  //newIntegTime=tdcChargeDrainParameterisation_[0]*pow(charge_kfC,2)+tdcChargeDrainParameterisation_[1]*charge_kfC+tdcChargeDrainParameterisation_[2];
305  } else if (charge_kfC < tdcChargeDrainParameterisation_[7]) {
306  poffset = 4;
307  charge_offset = tdcChargeDrainParameterisation_[3];
308  //newIntegTime=tdcChargeDrainParameterisation_[4]*pow(charge_kfC-tdcChargeDrainParameterisation_[3],2)+tdcChargeDrainParameterisation_[5]*(charge_kfC-tdcChargeDrainParameterisation_[3])+tdcChargeDrainParameterisation_[6];
309  } else {
310  poffset = 8;
311  charge_offset = tdcChargeDrainParameterisation_[7];
312  //newIntegTime=tdcChargeDrainParameterisation_[8]*pow(charge_kfC-tdcChargeDrainParameterisation_[7],2)+tdcChargeDrainParameterisation_[9]*(charge_kfC-tdcChargeDrainParameterisation_[7])+tdcChargeDrainParameterisation_[10];
313  }
314  const float charge_mod = charge_kfC - charge_offset;
315  const float newIntegTime =
316  ((tdcChargeDrainParameterisation_[poffset] * charge_mod + tdcChargeDrainParameterisation_[poffset + 1]) *
317  charge_mod +
318  tdcChargeDrainParameterisation_[poffset + 2]);
319 
320  const int newBusyBxs = std::floor(newIntegTime / 25.f) + 1;
321 
322  //if no update is needed regarding the number of bunches,
323  //then the ToT integration time has converged
324  integTime = newIntegTime;
325  if (newBusyBxs == busyBxs)
326  break;
327 
328  //update charge integrated during ToT
329  if (debug) {
330  if (busyBxs == 0)
331  edm::LogVerbatim("HGCFE") << "\t Intial busy estimate=" << integTime << " ns = " << newBusyBxs << " bxs"
332  << std::endl;
333  else
334  edm::LogVerbatim("HGCFE") << "\t ...integrated charge overflows initial busy estimate, interating again"
335  << std::endl;
336  }
337 
338  //update number of busy bunches
339  busyBxs = newBusyBxs;
340 
341  //reset charge to be integrated
342  totalCharge = charge;
343  if (toaMode_ == WEIGHTEDBYE)
344  finalToA = toa * charge;
345 
346  //add leakage from previous bunches in SARS ADC mode
347  for (int jt = 0; jt < it; ++jt) {
348  const unsigned int deltaT = (it - jt);
349  if ((deltaT + 2) >= adcPulse.size() || chargeColl[jt] == 0.f || totFlags[jt] || busyFlags[jt])
350  continue;
351 
352  const float leakCharge = chargeColl[jt] * adcPulse[deltaT + 2];
353  totalCharge += leakCharge;
354  if (toaMode_ == WEIGHTEDBYE)
355  finalToA += leakCharge * pulseAvgT_[deltaT + 2];
356 
357  if (debug)
358  edm::LogVerbatim("HGCFE") << "\t\t leaking " << chargeColl[jt] << " fC @ deltaT=-" << deltaT << " -> +"
359  << leakCharge << " with avgT=" << pulseAvgT_[deltaT + 2] << std::endl;
360  }
361 
362  //add contamination from posterior bunches
363  for (int jt = it + 1; jt < it + busyBxs && jt < dataFrame.size(); ++jt) {
364  //this charge will be integrated in TDC mode
365  //disable for SARS ADC
366  busyFlags[jt] = true;
367 
368  const float extraCharge = chargeColl[jt];
369  if (extraCharge == 0.f)
370  continue;
371  if (debug)
372  edm::LogVerbatim("HGCFE") << "\t\t adding " << extraCharge << " fC @ deltaT=+" << (jt - it) << std::endl;
373 
374  totalCharge += extraCharge;
375  if (toaMode_ == WEIGHTEDBYE)
376  finalToA += extraCharge * toaColl[jt];
377  }
378 
379  //finalize ToA contamination
380  if (toaMode_ == WEIGHTEDBYE)
381  finalToA /= totalCharge;
382  }
383 
384  newCharge[it] = (totalCharge - tdcOnset);
385 
386  if (debug)
387  edm::LogVerbatim("HGCFE") << "\t Final busy estimate=" << integTime << " ns = " << busyBxs << " bxs" << std::endl
388  << "\t Total integrated=" << totalCharge << " fC <toa>=" << toaFromToT[it]
389  << " (raw=" << finalToA << ") ns " << std::endl;
390 
391  //last fC (tdcOnset) are dissipated trough pulse
392  if (it + busyBxs < (int)(newCharge.size())) {
393  const float deltaT2nextBx((busyBxs * 25 - integTime));
394  const float tdcOnsetLeakage(tdcOnset * vdt::fast_expf(-deltaT2nextBx / tdcChargeDrainParameterisation_[11]));
395  if (debug)
396  edm::LogVerbatim("HGCFE") << "\t Leaking remainder of TDC onset " << tdcOnset << " fC, to be dissipated in "
397  << deltaT2nextBx << " DeltaT/tau=" << deltaT2nextBx << " / "
398  << tdcChargeDrainParameterisation_[11] << " ns, adds " << tdcOnsetLeakage << " fC @ "
399  << it + busyBxs << " bx (first free bx)" << std::endl;
400  newCharge[it + busyBxs] += tdcOnsetLeakage;
401  }
402  }
403 
404  //including the leakage from bunches in SARS ADC when not declared busy or in ToT
405  auto runChargeSharing = [&]() {
406  int ipulse = 0;
407  for (int it = 0; it < (int)(chargeColl.size()); ++it) {
408  //if busy, charge has been already integrated
409  //if(debug) edm::LogVerbatim("HGCFE") << "\t SARS ADC pulse activated @ " << it << " : ";
410  if (!totFlags[it] && !busyFlags[it]) {
411  const int start = std::max(0, 2 - it);
412  const int stop = std::min((int)adcPulse.size(), (int)newCharge.size() - it + 2);
413  for (ipulse = start; ipulse < stop; ++ipulse) {
414  const int itoffset = it + ipulse - 2;
415  //notice that if the channel is already busy,
416  //it has already been affected by the leakage of the SARS ADC
417  //if(totFlags[itoffset] || busyFlags[itoffset]) continue;
418  if (!totFlags[itoffset] && !busyFlags[itoffset]) {
419  newCharge[itoffset] += chargeColl[it] * adcPulse[ipulse];
420  }
421  //if(debug) edm::LogVerbatim("HGCFE") << " | " << itoffset << " " << chargeColl[it]*adcPulse[ipulse] << "( " << chargeColl[it] << "->";
422  //if(debug) edm::LogVerbatim("HGCFE") << newCharge[itoffset] << ") ";
423  }
424  }
425 
426  if (debug)
427  edm::LogVerbatim("HGCFE") << std::endl;
428  }
429  };
430  runChargeSharing();
431 
432  //For the future need to understand how to deal with toa for out of time signals
433  //and for that should keep track of the BX firing the ToA somewhere (also to restore the use of finalToA)
434  /*
435  float finalToA(0.);
436  for(int it=0; it<(int)(newCharge.size()); it++){
437  if(toaFlags[it]){
438  finalToA = toaFromToT[it];
439  //to avoid +=25 for small negative time taken as 0
440  while(finalToA < -1.e-5) finalToA+=25.f;
441  while(finalToA > 25.f) finalToA-=25.f;
442  toaFromToT[it] = finalToA;
443  }
444  }
445  */
446  //timeToA is already in 0-25ns range by construction
447 
448  //set new ADCs and ToA
449  if (debug)
450  edm::LogVerbatim("HGCFE") << "\t final result : ";
451  for (int it = 0; it < (int)(newCharge.size()); it++) {
452  if (debug)
453  edm::LogVerbatim("HGCFE") << chargeColl[it] << " -> " << newCharge[it] << " ";
454 
455  HGCSample newSample;
456  if (totFlags[it] || busyFlags[it]) {
457  if (totFlags[it]) {
458  //brute force saturation, maybe could to better with an exponential like saturation
459  const float saturatedCharge(std::min(newCharge[it], tdcSaturation_fC_));
460  //working version for in-time PU and signal
461  newSample.set(
462  true, true, gainIdx, (uint16_t)(timeToA / toaLSB_ns_), (uint16_t)(std::floor(saturatedCharge / tdcLSB_fC_)));
463  if (toaFlags[it])
464  newSample.setToAValid(true);
465  } else {
466  newSample.set(false, true, gainIdx, 0, 0);
467  }
468  } else {
469  //brute force saturation, maybe could to better with an exponential like saturation
470  const uint16_t adc = std::floor(std::min(newCharge[it], maxADC) / lsbADC);
471  //working version for in-time PU and signal
472  newSample.set(adc > thrADC, false, gainIdx, (uint16_t)(timeToA / toaLSB_ns_), adc);
473  if (toaFlags[it])
474  newSample.setToAValid(true);
475  }
476  dataFrame.setSample(it, newSample);
477  }
478 
479  if (debug) {
480  std::ostringstream msg;
481  dataFrame.print(msg);
482  edm::LogVerbatim("HGCFE") << msg.str() << std::endl;
483  }
484 }
485 
486 // cause the compiler to generate the appropriate code
488 template class HGCFEElectronics<HGCalDataFrame>;
Definition: start.py:1
Log< level::Info, true > LogVerbatim
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void runShaperWithToT(DFr &dataFrame, hgc::HGCSimHitData &chargeColl, hgc::HGCSimHitData &toa, CLHEP::HepRandomEngine *engine, uint32_t thrADC, float lsbADC, uint32_t gainIdx, float maxADC, int thickness, float tdcOnsetAuto, float noiseWidth, const hgc_digi::FEADCPulseShape &adcPulse)
implements pulse shape and switch to time over threshold including deadtime
wrapper for a data word
Definition: HGCSample.h:13
std::array< HGCSimData_t, nSamples > HGCSimHitData
HGCFEElectronics(const edm::ParameterSet &ps)
CTOR.
double f[11][100]
void setToAValid(bool toaFired)
Definition: HGCSample.h:47
#define debug
Definition: HDRShower.cc:19
tuple msg
Definition: mps_check.py:286
void runSimpleShaper(DFr &dataFrame, hgc::HGCSimHitData &chargeColl, uint32_t thrADC, float lsbADC, uint32_t gainIdx, float maxADC, const hgc_digi::FEADCPulseShape &adcPulse)
applies a shape to each time sample and propagates the tails to the subsequent time samples ...
float fast_expf(float x)
void runTrivialShaper(DFr &dataFrame, hgc::HGCSimHitData &chargeColl, uint32_t thrADC, float lsbADC, uint32_t gainIdx, float maxADC)
converts charge to digis without pulse shape
models the behavior of the front-end electronics
void set(bool thr, bool mode, uint16_t gain, uint16_t toa, uint16_t data)
Definition: HGCSample.h:49
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
uint16_t *__restrict__ uint16_t const *__restrict__ adc
std::array< float, 6 > FEADCPulseShape