CMS 3D CMS Logo

HcalElectronicsSim.cc
Go to the documentation of this file.
8 #include "CLHEP/Random/RandFlat.h"
9 #include <cmath>
10 
11 HcalElectronicsSim::HcalElectronicsSim(HcalAmplifier * amplifier, const HcalCoderFactory * coderFactory, bool PreMixing)
12  : theAmplifier(amplifier),
13  theCoderFactory(coderFactory),
14  theStartingCapId(0),
15  theStartingCapIdIsRandom(true),
16  PreMixDigis(PreMixing)
17 {
18 }
19 
20 
22 }
23 
25  // theAmplifier->setDbService(service);
26  theTDC.setDbService(service);
27 }
28 
29 template<class Digi>
30 void HcalElectronicsSim::convert(CaloSamples & frame, Digi & result, CLHEP::HepRandomEngine* engine) {
31  result.setSize(frame.size());
32  theAmplifier->amplify(frame, engine);
33  theCoderFactory->coder(frame.id())->fC2adc(frame, result, theStartingCapId);
34 }
35 
36 template<>
37 void HcalElectronicsSim::convert<QIE10DataFrame>(CaloSamples & frame, QIE10DataFrame & result, CLHEP::HepRandomEngine* engine) {
38  theAmplifier->amplify(frame, engine);
40 }
41 
42 template<>
43 void HcalElectronicsSim::convert<QIE11DataFrame>(CaloSamples & frame, QIE11DataFrame & result, CLHEP::HepRandomEngine* engine) {
44  theAmplifier->amplify(frame, engine);
46 }
47 
48 template<class Digi>
49 void HcalElectronicsSim::premix(CaloSamples & frame, Digi & result, double preMixFactor, unsigned preMixBits){
50  for(int isample = 0; isample !=frame.size(); ++isample) {
51  uint16_t theADC = round(preMixFactor*frame[isample]);
52  unsigned capId = result[isample].capid();
53 
54  if(theADC > preMixBits) {
55  uint16_t keepADC = result[isample].adc();
56  result.setSample(isample, HcalQIESample(keepADC, capId, 0, 0, true, true) ); // set error bit as a flag
57  }
58  else {
59  result.setSample(isample, HcalQIESample(theADC, capId, 0, 0) ); // preserve fC, no noise
60  }
61  }
62 }
63 
64 template<>
65 void HcalElectronicsSim::premix<QIE10DataFrame>(CaloSamples & frame, QIE10DataFrame & result, double preMixFactor, unsigned preMixBits){
66  for(int isample = 0; isample !=frame.size(); ++isample) {
67  uint16_t theADC = round(preMixFactor*frame[isample]);
68  unsigned capId = result[isample].capid();
69  bool ok = true;
70 
71  if(theADC > preMixBits) {
72  theADC = result[isample].adc();
73  ok = false; // set error bit as a flag
74  }
75 
76  result.setSample(isample, theADC, result[isample].le_tdc(), result[isample].te_tdc(), capId, result[isample].soi(), ok);
77  }
78 }
79 
80 template<>
81 void HcalElectronicsSim::premix<QIE11DataFrame>(CaloSamples & frame, QIE11DataFrame & result, double preMixFactor, unsigned preMixBits){
82  for(int isample = 0; isample !=frame.size(); ++isample) {
83  uint16_t theADC = round(preMixFactor*frame[isample]);
84  int tdcErrorBit = 0;
85 
86  if(theADC > preMixBits) {
87  theADC = result[isample].adc();
88  tdcErrorBit = 1; //use TDC bits for error bit
89  }
90 
91  result.setSample(isample, theADC, tdcErrorBit, result[isample].soi());
92  }
93 }
94 
95 template<class Digi>
96 void HcalElectronicsSim::analogToDigitalImpl(CLHEP::HepRandomEngine* engine, CaloSamples & lf, Digi & result, double preMixFactor, unsigned preMixBits) {
97  convert<Digi>(lf, result, engine);
98  if(PreMixDigis) premix(lf,result,preMixFactor,preMixBits);
99 }
100 
101 //TODO:
102 //HcalTDC extension for QIE10? and QIE11?
103 
104 void HcalElectronicsSim::analogToDigital(CLHEP::HepRandomEngine* engine, CaloSamples & lf, HBHEDataFrame & result, double preMixFactor, unsigned preMixBits) {
105  analogToDigitalImpl(engine,lf,result,preMixFactor,preMixBits);
106 }
107 
108 void HcalElectronicsSim::analogToDigital(CLHEP::HepRandomEngine* engine, CaloSamples & lf, HODataFrame & result, double preMixFactor, unsigned preMixBits) {
109  analogToDigitalImpl(engine,lf,result,preMixFactor,preMixBits);
110 }
111 
112 void HcalElectronicsSim::analogToDigital(CLHEP::HepRandomEngine* engine, CaloSamples & lf, HFDataFrame & result, double preMixFactor, unsigned preMixBits) {
113  analogToDigitalImpl(engine,lf,result,preMixFactor,preMixBits);
114 }
115 
116 void HcalElectronicsSim::analogToDigital(CLHEP::HepRandomEngine* engine, CaloSamples & lf, ZDCDataFrame & result, double preMixFactor, unsigned preMixBits) {
117  analogToDigitalImpl(engine,lf,result,preMixFactor,preMixBits);
118 }
119 
120 void HcalElectronicsSim::analogToDigital(CLHEP::HepRandomEngine* engine, CaloSamples & lf, QIE10DataFrame & result, double preMixFactor, unsigned preMixBits) {
121  analogToDigitalImpl(engine,lf,result,preMixFactor,preMixBits);
122 }
123 
124 void HcalElectronicsSim::analogToDigital(CLHEP::HepRandomEngine* engine, CaloSamples & lf, QIE11DataFrame & result, double preMixFactor, unsigned preMixBits) {
125  analogToDigitalImpl(engine,lf,result,preMixFactor,preMixBits);
126 }
127 
128 void HcalElectronicsSim::newEvent(CLHEP::HepRandomEngine* engine) {
129  // pick a new starting Capacitor ID
131  {
132  theStartingCapId = CLHEP::RandFlat::shootInt(engine, 4);
134  }
135 }
136 
138 {
139  theStartingCapId = startingCapId;
141  // turns off random capIDs forever for this instance
142  theStartingCapIdIsRandom = false;
143 }
HcalAmplifier * theAmplifier
void setStartingCapId(int startingCapId)
void setStartingCapId(int capId)
Definition: HcalAmplifier.h:46
std::tuple< unsigned int, int, int, DigiType, int, int, int, float > Digi
Definition: GenericDigi.h:40
void setDbService(const HcalDbService *service)
void newEvent(CLHEP::HepRandomEngine *)
HcalElectronicsSim(HcalAmplifier *amplifier, const HcalCoderFactory *coderFactory, bool PreMix)
void convert(CaloSamples &frame, Digi &result, CLHEP::HepRandomEngine *)
const HcalCoderFactory * theCoderFactory
virtual void amplify(CaloSamples &linearFrame, CLHEP::HepRandomEngine *) const
int size() const
get the size
Definition: CaloSamples.h:24
void analogToDigitalImpl(CLHEP::HepRandomEngine *, CaloSamples &linearFrame, Digi &result, double preMixFactor, unsigned preMixBits)
void analogToDigital(CLHEP::HepRandomEngine *, CaloSamples &linearFrame, HBHEDataFrame &result, double preMixFactor=10.0, unsigned preMixBits=126)
DetId id() const
get the (generic) id
Definition: CaloSamples.h:21
std::unique_ptr< HcalCoder > coder(const DetId &detId) const
user gets control of the pointer
void premix(CaloSamples &frame, Digi &result, double preMixFactor, unsigned preMixBits)
void setDbService(const HcalDbService *service)
the Producer will probably update this every event
Definition: HcalTDC.cc:133