File indexing completed on 2024-04-06 12:29:31
0001 #include "SimCalorimetry/HcalSimAlgos/interface/HcalElectronicsSim.h"
0002 #include "SimCalorimetry/HcalSimAlgos/interface/HcalSimParameters.h"
0003 #include "SimCalorimetry/HcalSimAlgos/interface/HcalTDC.h"
0004 #include "DataFormats/HcalDigi/interface/HBHEDataFrame.h"
0005 #include "DataFormats/HcalDigi/interface/HODataFrame.h"
0006 #include "DataFormats/HcalDigi/interface/HFDataFrame.h"
0007 #include "DataFormats/HcalDigi/interface/ZDCDataFrame.h"
0008 #include "DataFormats/HcalDigi/interface/QIE10DataFrame.h"
0009 #include "DataFormats/HcalDigi/interface/QIE11DataFrame.h"
0010 #include "CLHEP/Random/RandFlat.h"
0011 #include <cmath>
0012
0013 HcalElectronicsSim::HcalElectronicsSim(const HcalSimParameterMap* parameterMap,
0014 HcalAmplifier* amplifier,
0015 const HcalCoderFactory* coderFactory,
0016 bool PreMixing)
0017 : theParameterMap(parameterMap),
0018 theAmplifier(amplifier),
0019 theCoderFactory(coderFactory),
0020 theStartingCapId(0),
0021 theStartingCapIdIsRandom(true),
0022 PreMixDigis(PreMixing) {}
0023
0024 HcalElectronicsSim::~HcalElectronicsSim() {}
0025
0026 void HcalElectronicsSim::setDbService(const HcalDbService* service) {
0027
0028 theTDC.setDbService(service);
0029 }
0030
0031 template <class Digi>
0032 void HcalElectronicsSim::convert(CaloSamples& frame, Digi& result, CLHEP::HepRandomEngine* engine) {
0033 result.setSize(frame.size());
0034 theAmplifier->amplify(frame, engine);
0035 theCoderFactory->coder(frame.id())->fC2adc(frame, result, theStartingCapId);
0036 }
0037
0038 template <>
0039 void HcalElectronicsSim::convert<QIE10DataFrame>(CaloSamples& frame,
0040 QIE10DataFrame& result,
0041 CLHEP::HepRandomEngine* engine) {
0042 theAmplifier->amplify(frame, engine);
0043 theCoderFactory->coder(frame.id())->fC2adc(frame, result, theStartingCapId);
0044 }
0045
0046 template <>
0047 void HcalElectronicsSim::convert<QIE11DataFrame>(CaloSamples& frame,
0048 QIE11DataFrame& result,
0049 CLHEP::HepRandomEngine* engine) {
0050 theAmplifier->amplify(frame, engine);
0051 theCoderFactory->coder(frame.id())->fC2adc(frame, result, theStartingCapId);
0052 }
0053
0054 template <class Digi>
0055 void HcalElectronicsSim::premix(CaloSamples& frame, Digi& result, double preMixFactor, unsigned preMixBits) {
0056 for (int isample = 0; isample != frame.size(); ++isample) {
0057 uint16_t theADC = round(preMixFactor * frame[isample]);
0058 unsigned capId = result[isample].capid();
0059
0060 if (theADC > preMixBits) {
0061 uint16_t keepADC = result[isample].adc();
0062 result.setSample(isample, HcalQIESample(keepADC, capId, 0, 0, true, true));
0063 } else {
0064 result.setSample(isample, HcalQIESample(theADC, capId, 0, 0));
0065 }
0066 }
0067 }
0068
0069 template <>
0070 void HcalElectronicsSim::premix<QIE10DataFrame>(CaloSamples& frame,
0071 QIE10DataFrame& result,
0072 double preMixFactor,
0073 unsigned preMixBits) {
0074 for (int isample = 0; isample != frame.size(); ++isample) {
0075 uint16_t theADC = round(preMixFactor * frame[isample]);
0076 unsigned capId = result[isample].capid();
0077 bool ok = true;
0078
0079 if (theADC > preMixBits) {
0080 theADC = result[isample].adc();
0081 ok = false;
0082 }
0083
0084 result.setSample(
0085 isample, theADC, result[isample].le_tdc(), result[isample].te_tdc(), capId, result[isample].soi(), ok);
0086 }
0087 }
0088
0089 template <>
0090 void HcalElectronicsSim::premix<QIE11DataFrame>(CaloSamples& frame,
0091 QIE11DataFrame& result,
0092 double preMixFactor,
0093 unsigned preMixBits) {
0094 for (int isample = 0; isample != frame.size(); ++isample) {
0095 uint16_t theADC = round(preMixFactor * frame[isample]);
0096 int tdcErrorBit = 0;
0097
0098 if (theADC > preMixBits) {
0099 theADC = result[isample].adc();
0100 tdcErrorBit = 1;
0101 }
0102
0103 result.setSample(isample, theADC, tdcErrorBit, result[isample].soi());
0104 }
0105 }
0106
0107 template <class Digi>
0108 void HcalElectronicsSim::analogToDigitalImpl(
0109 CLHEP::HepRandomEngine* engine, CaloSamples& lf, Digi& result, double preMixFactor, unsigned preMixBits) {
0110 convert<Digi>(lf, result, engine);
0111 if (PreMixDigis)
0112 premix(lf, result, preMixFactor, preMixBits);
0113 }
0114
0115
0116
0117
0118 void HcalElectronicsSim::analogToDigital(
0119 CLHEP::HepRandomEngine* engine, CaloSamples& lf, HBHEDataFrame& result, double preMixFactor, unsigned preMixBits) {
0120 analogToDigitalImpl(engine, lf, result, preMixFactor, preMixBits);
0121 }
0122
0123 void HcalElectronicsSim::analogToDigital(
0124 CLHEP::HepRandomEngine* engine, CaloSamples& lf, HODataFrame& result, double preMixFactor, unsigned preMixBits) {
0125 analogToDigitalImpl(engine, lf, result, preMixFactor, preMixBits);
0126 }
0127
0128 void HcalElectronicsSim::analogToDigital(
0129 CLHEP::HepRandomEngine* engine, CaloSamples& lf, HFDataFrame& result, double preMixFactor, unsigned preMixBits) {
0130 analogToDigitalImpl(engine, lf, result, preMixFactor, preMixBits);
0131 }
0132
0133 void HcalElectronicsSim::analogToDigital(
0134 CLHEP::HepRandomEngine* engine, CaloSamples& lf, ZDCDataFrame& result, double preMixFactor, unsigned preMixBits) {
0135 analogToDigitalImpl(engine, lf, result, preMixFactor, preMixBits);
0136 }
0137
0138 void HcalElectronicsSim::analogToDigital(
0139 CLHEP::HepRandomEngine* engine, CaloSamples& lf, QIE10DataFrame& result, double preMixFactor, unsigned preMixBits) {
0140 analogToDigitalImpl(engine, lf, result, preMixFactor, preMixBits);
0141 if (!PreMixDigis) {
0142 const HFSimParameters& pars = static_cast<const HFSimParameters&>(theParameterMap->simParameters(lf.id()));
0143 if (pars.threshold_currentTDC() > 0.) {
0144 HcalTDC theTDC((pars.threshold_currentTDC()));
0145 theTDC.timing(lf, result);
0146 }
0147 }
0148 }
0149
0150 void HcalElectronicsSim::analogToDigital(
0151 CLHEP::HepRandomEngine* engine, CaloSamples& lf, QIE11DataFrame& result, double preMixFactor, unsigned preMixBits) {
0152 analogToDigitalImpl(engine, lf, result, preMixFactor, preMixBits);
0153 if (!PreMixDigis) {
0154 const HcalSimParameters& pars = static_cast<const HcalSimParameters&>(theParameterMap->simParameters(lf.id()));
0155 if (pars.threshold_currentTDC() > 0.) {
0156 HcalTDC theTDC((pars.threshold_currentTDC()));
0157 theTDC.timing(lf, result);
0158 }
0159 }
0160 }
0161
0162 void HcalElectronicsSim::newEvent(CLHEP::HepRandomEngine* engine) {
0163
0164 if (theStartingCapIdIsRandom) {
0165 theStartingCapId = CLHEP::RandFlat::shootInt(engine, 4);
0166 theAmplifier->setStartingCapId(theStartingCapId);
0167 }
0168 }
0169
0170 void HcalElectronicsSim::setStartingCapId(int startingCapId) {
0171 theStartingCapId = startingCapId;
0172 theAmplifier->setStartingCapId(theStartingCapId);
0173
0174 theStartingCapIdIsRandom = false;
0175 }