CMS 3D CMS Logo

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