Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

70 x hcal pre mix bugfix #4325

Merged
merged 5 commits into from Jun 20, 2014
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
2 changes: 1 addition & 1 deletion EventFilter/HcalRawToDigi/src/HcalHTRData.cc
Expand Up @@ -248,7 +248,7 @@ void HcalHTRData::pack(unsigned char* daq_lengths, unsigned short* daq_samples,
for (ichan=0; ichan<24; ichan++) {
unsigned short chanid=((ichan%3)+((ichan/3)<<2))<<11;
for (isample=0; isample<daq_lengths[ichan] && isample<MAXIMUM_SAMPLES_PER_CHANNEL; isample++) {
unsigned short basedata=daq_samples[ichan*MAXIMUM_SAMPLES_PER_CHANNEL+isample]&0x3FF;
unsigned short basedata=daq_samples[ichan*MAXIMUM_SAMPLES_PER_CHANNEL+isample]&0x7FF;
if (do_capid) basedata=(basedata&0x7F)|(0x200)|((isample%4)<<7);
ptr[daq_words_total]=chanid|basedata;
daq_words_total++;
Expand Down
3 changes: 2 additions & 1 deletion SimCalorimetry/HcalSimAlgos/interface/HcalElectronicsSim.h
Expand Up @@ -21,7 +21,7 @@ class HcalCoderFactory;
class HcalElectronicsSim {
public:
HcalElectronicsSim(HcalAmplifier * amplifier,
const HcalCoderFactory * coderFactory);
const HcalCoderFactory * coderFactory, bool PreMix);
~HcalElectronicsSim();

void setRandomEngine(CLHEP::HepRandomEngine & engine);
Expand All @@ -47,6 +47,7 @@ class HcalElectronicsSim {

int theStartingCapId;
bool theStartingCapIdIsRandom;
bool PreMixDigis;
};


Expand Down
48 changes: 39 additions & 9 deletions SimCalorimetry/HcalSimAlgos/interface/HcalSignalGenerator.h
Expand Up @@ -13,6 +13,7 @@
#include "SimCalorimetry/HcalSimAlgos/interface/HcalElectronicsSim.h"
#include "SimCalorimetry/HcalSimAlgos/interface/HcalDigitizerTraits.h"
#include "DataFormats/Common/interface/Handle.h"
#include "DataFormats/HcalDetId/interface/HcalSubdetector.h"

/** Converts digis back into analog signals, to be used
* as noise
Expand Down Expand Up @@ -88,7 +89,6 @@ class HcalSignalGenerator : public HcalBaseSignalGenerator

if (digis)
{

// loop over digis, adding these to the existing maps
for(typename COLLECTION::const_iterator it = digis->begin();
it != digis->end(); ++it)
Expand Down Expand Up @@ -120,17 +120,47 @@ class HcalSignalGenerator : public HcalBaseSignalGenerator

CaloSamples samplesInPE(const DIGI & digi)
{

// calibration, for future reference: (same block for all Hcal types)
HcalDetId cell = digi.id();
// const HcalCalibrations& calibrations=conditions->getHcalCalibrations(cell);
const HcalQIECoder* channelCoder = theConditions->getHcalCoder (cell);
const HcalQIEShape* channelShape = theConditions->getHcalShape (cell);
HcalCoderDb coder (*channelCoder, *channelShape);
CaloSamples result;
coder.adc2fC(digi, result);
fC2pe(result);

// std::cout << " HcalSignalGenerator: noise input " << digi << std::endl;
CaloSamples result = CaloSamples(digi.id(),digi.size());

// first, check if there was an overflow in this fake digi:
bool overflow = false;
// find and list them

for(int isample=0; isample<digi.size(); ++isample) {
if(digi[isample].er()) overflow = true;
}

if(overflow) { // do full conversion, go back and overwrite fake entries

const HcalQIECoder* channelCoder = theConditions->getHcalCoder (cell);
const HcalQIEShape* channelShape = theConditions->getHcalShape (cell);
HcalCoderDb coder (*channelCoder, *channelShape);
coder.adc2fC(digi, result);

// overwrite with coded information
for(int isample=0; isample<digi.size(); ++isample) {
if(!digi[isample].er()) result[isample] = float(digi[isample].adc())/10.;
}
}
else { // saves creating the coder, etc., every time
// use coded information
for(int isample=0; isample<digi.size(); ++isample) {
result[isample] = float(digi[isample].adc())/10.;
}
result.setPresamples(digi.presamples());
}

std::cout << " HcalSignalGenerator: noise input ADC " << digi << std::endl;

std::cout << " HcalSignalGenerator: noise input in fC " << result << std::endl;

// translation done in fC, convert to pe:

fC2pe(result);

return result;
}
Expand Down
81 changes: 73 additions & 8 deletions SimCalorimetry/HcalSimAlgos/src/HcalElectronicsSim.cc
Expand Up @@ -7,15 +7,15 @@
#include "DataFormats/HcalDigi/interface/ZDCDataFrame.h"
#include "DataFormats/HcalDigi/interface/HcalUpgradeDataFrame.h"
#include "CLHEP/Random/RandFlat.h"
#include <math.h>



HcalElectronicsSim::HcalElectronicsSim(HcalAmplifier * amplifier, const HcalCoderFactory * coderFactory)
HcalElectronicsSim::HcalElectronicsSim(HcalAmplifier * amplifier, const HcalCoderFactory * coderFactory, bool PreMixing)
: theAmplifier(amplifier),
theCoderFactory(coderFactory),
theRandFlat(0),
theStartingCapId(0),
theStartingCapIdIsRandom(true)
theStartingCapIdIsRandom(false),
PreMixDigis(PreMixing)
{
}

Expand All @@ -41,26 +41,91 @@ template<class Digi>
void HcalElectronicsSim::convert(CaloSamples & frame, Digi & result) {
result.setSize(frame.size());
theAmplifier->amplify(frame);
theCoderFactory->coder(frame.id())->fC2adc(frame, result, theStartingCapId);

theCoderFactory->coder(frame.id())->fC2adc(frame, result, theStartingCapId);
}


void HcalElectronicsSim::analogToDigital(CaloSamples & lf, HBHEDataFrame & result) {
convert<HBHEDataFrame>(lf, result);
convert<HBHEDataFrame>(lf, result);
if(PreMixDigis) {
for(int isample = 0; isample !=lf.size(); ++isample) {
uint16_t theADC = round(10.0*lf[isample]);
unsigned capId = result[isample].capid();

if(theADC > 126) {
uint16_t keepADC = result[isample].adc();
HcalQIESample mysamp(keepADC, capId, 0, 0, true, true); // set error bit as a flag
result.setSample(isample, HcalQIESample(keepADC, capId, 0, 0, true, true) );
}
else {
result.setSample(isample, HcalQIESample(theADC, capId, 0, 0) ); // preserve fC, no noise
HcalQIESample mysamp(theADC, capId, 0, 0);
}
}
}
}


void HcalElectronicsSim::analogToDigital(CaloSamples & lf, HODataFrame & result) {
convert<HODataFrame>(lf, result);
if(PreMixDigis) {
for(int isample = 0; isample !=lf.size(); ++isample) {
uint16_t theADC = round(10.0*lf[isample]);
unsigned capId = result[isample].capid();

if(theADC > 126) {
uint16_t keepADC = result[isample].adc();
HcalQIESample mysamp(keepADC, capId, 0, 0, true, true);// set error bit as a flag
result.setSample(isample, HcalQIESample(keepADC, capId, 0, 0, true, true) );
}
else {
result.setSample(isample, HcalQIESample(theADC, capId, 0, 0) ); // preserve fC, no noise
HcalQIESample mysamp(theADC, capId, 0, 0);
}
}
}
}


void HcalElectronicsSim::analogToDigital(CaloSamples & lf, HFDataFrame & result) {
convert<HFDataFrame>(lf, result);
convert<HFDataFrame>(lf, result); //not enough dynamic range with simple multiplication
if(PreMixDigis) {
for(int isample = 0; isample !=lf.size(); ++isample) {
uint16_t theADC = round(10.0*lf[isample]);
unsigned capId = result[isample].capid();

if(theADC > 126) {
uint16_t keepADC = result[isample].adc();
HcalQIESample mysamp(keepADC, capId, 0, 0, true, true);// set error bit as a flag
result.setSample(isample, HcalQIESample(keepADC, capId, 0, 0, true, true) );
}
else {
result.setSample(isample, HcalQIESample(theADC, capId, 0, 0) ); // preserve fC, no noise
HcalQIESample mysamp(theADC, capId, 0, 0);
}
}
}
}

void HcalElectronicsSim::analogToDigital(CaloSamples & lf, ZDCDataFrame & result) {
convert<ZDCDataFrame>(lf, result);
convert<ZDCDataFrame>(lf, result);
if(PreMixDigis) {
for(int isample = 0; isample !=lf.size(); ++isample) {
uint16_t theADC = round(10.0*lf[isample]);
unsigned capId = result[isample].capid();

if(theADC > 126) {
uint16_t keepADC = result[isample].adc();
HcalQIESample mysamp(keepADC, capId, 0, 0, true, true);//set error bit as a flag
result.setSample(isample, HcalQIESample(keepADC, capId, 0, 0, true, true) );
}
else {
result.setSample(isample, HcalQIESample(theADC, capId, 0, 0) ); // preserve fC, no noise
HcalQIESample mysamp(theADC, capId, 0, 0);
}
}
}
}


Expand Down
2 changes: 1 addition & 1 deletion SimCalorimetry/HcalSimAlgos/test/HcalDigitizerTest.cpp
Expand Up @@ -223,7 +223,7 @@ int main() {
bool PM2 = false;
HcalAmplifier amplifier(&parameterMap, addNoise, PM1, PM2);
HcalCoderFactory coderFactory(HcalCoderFactory::NOMINAL);
HcalElectronicsSim electronicsSim(&amplifier, &coderFactory);
HcalElectronicsSim electronicsSim(&amplifier, &coderFactory,PM1);
amplifier.setDbService(&calibratorHandle);
amplifier.setRandomEngine(randomEngine);
electronicsSim.setRandomEngine(randomEngine);
Expand Down
12 changes: 6 additions & 6 deletions SimCalorimetry/HcalSimProducers/src/HcalDigitizer.cc
Expand Up @@ -184,12 +184,12 @@ HcalDigitizer::HcalDigitizer(const edm::ParameterSet& ps, edm::ConsumesCollector

// std::cout << "HcalDigitizer: theUpgradeCoderFactory created" << std::endl;

theHBHEElectronicsSim = new HcalElectronicsSim(theHBHEAmplifier, theCoderFactory);
theHFElectronicsSim = new HcalElectronicsSim(theHFAmplifier, theCoderFactory);
theHOElectronicsSim = new HcalElectronicsSim(theHOAmplifier, theCoderFactory);
theZDCElectronicsSim = new HcalElectronicsSim(theZDCAmplifier, theCoderFactory);
theUpgradeHBHEElectronicsSim = new HcalElectronicsSim(theHBHEAmplifier, theUpgradeCoderFactory);
theUpgradeHFElectronicsSim = new HcalElectronicsSim(theHFAmplifier, theUpgradeCoderFactory);
theHBHEElectronicsSim = new HcalElectronicsSim(theHBHEAmplifier, theCoderFactory, PreMix1);
theHFElectronicsSim = new HcalElectronicsSim(theHFAmplifier, theCoderFactory, PreMix1);
theHOElectronicsSim = new HcalElectronicsSim(theHOAmplifier, theCoderFactory, PreMix1);
theZDCElectronicsSim = new HcalElectronicsSim(theZDCAmplifier, theCoderFactory, PreMix1);
theUpgradeHBHEElectronicsSim = new HcalElectronicsSim(theHBHEAmplifier, theUpgradeCoderFactory, PreMix1);
theUpgradeHFElectronicsSim = new HcalElectronicsSim(theHFAmplifier, theUpgradeCoderFactory, PreMix1);

// std::cout << "HcalDigitizer: theUpgradeElectronicsSim created" << std::endl;

Expand Down
2 changes: 1 addition & 1 deletion SimCalorimetry/HcalTestBeam/src/HcalTBDigiProducer.cc
Expand Up @@ -48,7 +48,7 @@ HcalTBDigiProducer::HcalTBDigiProducer(const edm::ParameterSet& ps, edm::one::ED
bool dummy2 = false; // extra arguments for premixing
theAmplifier = new HcalAmplifier(theParameterMap, doNoise, dummy1, dummy2);
theCoderFactory = new HcalCoderFactory(HcalCoderFactory::DB);
theElectronicsSim = new HcalElectronicsSim(theAmplifier, theCoderFactory);
theElectronicsSim = new HcalElectronicsSim(theAmplifier, theCoderFactory, dummy1);

theHBHEDigitizer = new HBHEDigitizer(theHBHEResponse, theElectronicsSim, doNoise);
theHODigitizer = new HODigitizer(theHOResponse, theElectronicsSim, doNoise);
Expand Down
Expand Up @@ -54,6 +54,8 @@
#muonRPCDigis.InputLabel = 'rawDataCollector'
#castorDigis.InputLabel = 'rawDataCollector'

hcalDigis.FilterDataQuality = cms.bool(False)

hcalSimBlock.HcalPreMixStage2 = cms.bool(True)

mixData = cms.EDProducer("DataMixingModule",
Expand Down