CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HcalElectronicsSim.cc
Go to the documentation of this file.
9 #include "CLHEP/Random/RandFlat.h"
10 #include <math.h>
11 
12 HcalElectronicsSim::HcalElectronicsSim(HcalAmplifier * amplifier, const HcalCoderFactory * coderFactory, bool PreMixing)
13  : theAmplifier(amplifier),
14  theCoderFactory(coderFactory),
15  theStartingCapId(0),
16  theStartingCapIdIsRandom(true),
17  PreMixDigis(PreMixing)
18 {
19 }
20 
21 
23 }
24 
25 
27  // theAmplifier->setDbService(service);
28  theTDC.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);
35  theCoderFactory->coder(frame.id())->fC2adc(frame, result, theStartingCapId);
36 }
37 
38 
39 void HcalElectronicsSim::analogToDigital(CLHEP::HepRandomEngine* engine, CaloSamples & lf, HBHEDataFrame & result) {
40  convert<HBHEDataFrame>(lf, result, engine);
41  if(PreMixDigis) {
42  for(int isample = 0; isample !=lf.size(); ++isample) {
43  uint16_t theADC = round(10.0*lf[isample]);
44  unsigned capId = result[isample].capid();
45 
46  if(theADC > 126) {
47  uint16_t keepADC = result[isample].adc();
48  HcalQIESample mysamp(keepADC, capId, 0, 0, true, true); // set error bit as a flag
49  result.setSample(isample, HcalQIESample(keepADC, capId, 0, 0, true, true) );
50  }
51  else {
52  result.setSample(isample, HcalQIESample(theADC, capId, 0, 0) ); // preserve fC, no noise
53  HcalQIESample mysamp(theADC, capId, 0, 0);
54  }
55  }
56  }
57 }
58 
59 
60 void HcalElectronicsSim::analogToDigital(CLHEP::HepRandomEngine* engine, CaloSamples & lf, HODataFrame & result) {
61  convert<HODataFrame>(lf, result, engine);
62  if(PreMixDigis) {
63  for(int isample = 0; isample !=lf.size(); ++isample) {
64  uint16_t theADC = round(10.0*lf[isample]);
65  unsigned capId = result[isample].capid();
66 
67  if(theADC > 126) {
68  uint16_t keepADC = result[isample].adc();
69  HcalQIESample mysamp(keepADC, capId, 0, 0, true, true);// set error bit as a flag
70  result.setSample(isample, HcalQIESample(keepADC, capId, 0, 0, true, true) );
71  }
72  else {
73  result.setSample(isample, HcalQIESample(theADC, capId, 0, 0) ); // preserve fC, no noise
74  HcalQIESample mysamp(theADC, capId, 0, 0);
75  }
76  }
77  }
78 }
79 
80 
81 void HcalElectronicsSim::analogToDigital(CLHEP::HepRandomEngine* engine, CaloSamples & lf, HFDataFrame & result) {
82  convert<HFDataFrame>(lf, result, engine);
83  if(PreMixDigis) {
84  for(int isample = 0; isample !=lf.size(); ++isample) {
85  uint16_t theADC = round(10.0*lf[isample]);
86  unsigned capId = result[isample].capid();
87 
88  if(theADC > 126) {
89  uint16_t keepADC = result[isample].adc();
90  HcalQIESample mysamp(keepADC, capId, 0, 0, true, true);// set error bit as a flag
91  result.setSample(isample, HcalQIESample(keepADC, capId, 0, 0, true, true) );
92  }
93  else {
94  result.setSample(isample, HcalQIESample(theADC, capId, 0, 0) ); // preserve fC, no noise
95  HcalQIESample mysamp(theADC, capId, 0, 0);
96  }
97  }
98  }
99 }
100 
101 void HcalElectronicsSim::analogToDigital(CLHEP::HepRandomEngine* engine, CaloSamples & lf, ZDCDataFrame & result) {
102  convert<ZDCDataFrame>(lf, result, engine);
103  if(PreMixDigis) {
104  for(int isample = 0; isample !=lf.size(); ++isample) {
105  uint16_t theADC = round(10.0*lf[isample]);
106  unsigned capId = result[isample].capid();
107 
108  if(theADC > 126) {
109  uint16_t keepADC = result[isample].adc();
110  HcalQIESample mysamp(keepADC, capId, 0, 0, true, true);//set error bit as a flag
111  result.setSample(isample, HcalQIESample(keepADC, capId, 0, 0, true, true) );
112  }
113  else {
114  result.setSample(isample, HcalQIESample(theADC, capId, 0, 0) ); // preserve fC, no noise
115  HcalQIESample mysamp(theADC, capId, 0, 0);
116  }
117  }
118  }
119 }
120 
121 
122 void HcalElectronicsSim::analogToDigital(CLHEP::HepRandomEngine* engine, CaloSamples & lf,
124  convert<HcalUpgradeDataFrame>(lf, result, engine);
125 // std::cout << HcalDetId(lf.id()) << ' ' << lf;
126  theTDC.timing(lf, result, engine);
127 }
128 
129 void HcalElectronicsSim::newEvent(CLHEP::HepRandomEngine* engine) {
130  // pick a new starting Capacitor ID
132  {
133  theStartingCapId = CLHEP::RandFlat::shootInt(engine, 4);
135  }
136 }
137 
138 
140 {
141  theStartingCapId = startingCapId;
143  // turns off random capIDs forever for this instance
144  theStartingCapIdIsRandom = false;
145 }
146 
HcalAmplifier * theAmplifier
void setStartingCapId(int startingCapId)
void analogToDigital(CLHEP::HepRandomEngine *, CaloSamples &linearFrame, HBHEDataFrame &result)
void setStartingCapId(int capId)
Definition: HcalAmplifier.h:44
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 *)
tuple result
Definition: query.py:137
const HcalCoderFactory * theCoderFactory
virtual void amplify(CaloSamples &linearFrame, CLHEP::HepRandomEngine *) const
int size() const
get the size
Definition: CaloSamples.h:24
void timing(const CaloSamples &lf, HcalUpgradeDataFrame &digi, CLHEP::HepRandomEngine *) const
adds timing information to the digi
Definition: HcalTDC.cc:16
DetId id() const
get the (generic) id
Definition: CaloSamples.h:21
void setDbService(const HcalDbService *service)
the Producer will probably update this every event
Definition: HcalTDC.cc:132
std::auto_ptr< HcalCoder > coder(const DetId &detId) const
user gets control of the pointer