CMS 3D CMS Logo

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