test
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 
10 #include "CLHEP/Random/RandFlat.h"
11 #include <math.h>
12 
13 HcalElectronicsSim::HcalElectronicsSim(HcalAmplifier * amplifier, const HcalCoderFactory * coderFactory, bool PreMixing)
14  : theAmplifier(amplifier),
15  theCoderFactory(coderFactory),
16  theStartingCapId(0),
17  theStartingCapIdIsRandom(true),
18  PreMixDigis(PreMixing)
19 {
20 }
21 
22 
24 }
25 
26 
28  // theAmplifier->setDbService(service);
29  theTDC.setDbService(service);
30 }
31 
32 template<class Digi>
33 void HcalElectronicsSim::convert(CaloSamples & frame, Digi & result, CLHEP::HepRandomEngine* engine) {
34  result.setSize(frame.size());
35  theAmplifier->amplify(frame, engine);
36  theCoderFactory->coder(frame.id())->fC2adc(frame, result, theStartingCapId);
37 }
38 
39 
40 void HcalElectronicsSim::analogToDigital(CLHEP::HepRandomEngine* engine, CaloSamples & lf, HBHEDataFrame & result) {
41  convert<HBHEDataFrame>(lf, result, engine);
42  if(PreMixDigis) {
43  for(int isample = 0; isample !=lf.size(); ++isample) {
44  uint16_t theADC = round(10.0*lf[isample]);
45  unsigned capId = result[isample].capid();
46 
47  if(theADC > 126) {
48  uint16_t keepADC = result[isample].adc();
49  HcalQIESample mysamp(keepADC, capId, 0, 0, true, true); // set error bit as a flag
50  result.setSample(isample, HcalQIESample(keepADC, capId, 0, 0, true, true) );
51  }
52  else {
53  result.setSample(isample, HcalQIESample(theADC, capId, 0, 0) ); // preserve fC, no noise
54  HcalQIESample mysamp(theADC, capId, 0, 0);
55  }
56  }
57  }
58 }
59 
60 
61 void HcalElectronicsSim::analogToDigital(CLHEP::HepRandomEngine* engine, CaloSamples & lf, HODataFrame & result) {
62  convert<HODataFrame>(lf, result, engine);
63  if(PreMixDigis) {
64  for(int isample = 0; isample !=lf.size(); ++isample) {
65  uint16_t theADC = round(10.0*lf[isample]);
66  unsigned capId = result[isample].capid();
67 
68  if(theADC > 126) {
69  uint16_t keepADC = result[isample].adc();
70  HcalQIESample mysamp(keepADC, capId, 0, 0, true, true);// set error bit as a flag
71  result.setSample(isample, HcalQIESample(keepADC, capId, 0, 0, true, true) );
72  }
73  else {
74  result.setSample(isample, HcalQIESample(theADC, capId, 0, 0) ); // preserve fC, no noise
75  HcalQIESample mysamp(theADC, capId, 0, 0);
76  }
77  }
78  }
79 }
80 
81 
82 void HcalElectronicsSim::analogToDigital(CLHEP::HepRandomEngine* engine, CaloSamples & lf, HFDataFrame & result) {
83  convert<HFDataFrame>(lf, result, engine);
84  if(PreMixDigis) {
85  for(int isample = 0; isample !=lf.size(); ++isample) {
86  uint16_t theADC = round(10.0*lf[isample]);
87  unsigned capId = result[isample].capid();
88 
89  if(theADC > 126) {
90  uint16_t keepADC = result[isample].adc();
91  HcalQIESample mysamp(keepADC, capId, 0, 0, true, true);// set error bit as a flag
92  result.setSample(isample, HcalQIESample(keepADC, capId, 0, 0, true, true) );
93  }
94  else {
95  result.setSample(isample, HcalQIESample(theADC, capId, 0, 0) ); // preserve fC, no noise
96  HcalQIESample mysamp(theADC, capId, 0, 0);
97  }
98  }
99  }
100 }
101 
102 void HcalElectronicsSim::analogToDigital(CLHEP::HepRandomEngine* engine, CaloSamples & lf, ZDCDataFrame & result) {
103  convert<ZDCDataFrame>(lf, result, engine);
104  if(PreMixDigis) {
105  for(int isample = 0; isample !=lf.size(); ++isample) {
106  uint16_t theADC = round(10.0*lf[isample]);
107  unsigned capId = result[isample].capid();
108 
109  if(theADC > 126) {
110  uint16_t keepADC = result[isample].adc();
111  HcalQIESample mysamp(keepADC, capId, 0, 0, true, true);//set error bit as a flag
112  result.setSample(isample, HcalQIESample(keepADC, capId, 0, 0, true, true) );
113  }
114  else {
115  result.setSample(isample, HcalQIESample(theADC, capId, 0, 0) ); // preserve fC, no noise
116  HcalQIESample mysamp(theADC, capId, 0, 0);
117  }
118  }
119  }
120 }
121 
122 
123 void HcalElectronicsSim::analogToDigital(CLHEP::HepRandomEngine* engine, CaloSamples & lf,
125  convert<HcalUpgradeDataFrame>(lf, result, engine);
126 // std::cout << HcalDetId(lf.id()) << ' ' << lf;
127  theTDC.timing(lf, result, engine);
128 }
129 
130 void HcalElectronicsSim::newEvent(CLHEP::HepRandomEngine* engine) {
131  // pick a new starting Capacitor ID
133  {
134  theStartingCapId = CLHEP::RandFlat::shootInt(engine, 4);
136  }
137 }
138 
139 
141 {
142  theStartingCapId = startingCapId;
144  // turns off random capIDs forever for this instance
145  theStartingCapIdIsRandom = false;
146 }
147 
HcalAmplifier * theAmplifier
void setStartingCapId(int startingCapId)
void analogToDigital(CLHEP::HepRandomEngine *, CaloSamples &linearFrame, HBHEDataFrame &result)
void setStartingCapId(int capId)
Definition: HcalAmplifier.h:43
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