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

Reduced radom number calls in the runSimple Method of HGCDigitizerBase for Potential Speedup. #27980

Merged
merged 7 commits into from Sep 24, 2019
Merged
Show file tree
Hide file tree
Changes from 3 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
10 changes: 10 additions & 0 deletions SimCalorimetry/HGCalSimProducers/interface/HGCDigitizerBase.h
Expand Up @@ -68,6 +68,10 @@ class HGCDigitizerBase {
@short CTOR
*/
HGCDigitizerBase(const edm::ParameterSet& ps);
/**
@short Gaussian Noise Generation Member Function
*/
void GenerateGaussianNoise(CLHEP::HepRandomEngine* engine, const double NoiseMean, const double NoiseStd);
/**
@short steer digitization mode
*/
Expand Down Expand Up @@ -150,6 +154,12 @@ class HGCDigitizerBase {

//if set to true, threshold will be computed based on the expected meap peak/2
bool thresholdFollowsMIP_;

// New NoiseArray Parameters
const double NoiseMean_, NoiseStd_;
static const size_t NoiseArrayLength_ = 20000;
static const size_t samplesize_ = 15;
std::array<std::array<double, samplesize_>, NoiseArrayLength_> GaussianNoiseArray_;
};

#endif
Expand Up @@ -6,7 +6,7 @@

#include "CLHEP/Random/RandGauss.h"
#include "CLHEP/Random/RandGaussQ.h"

#include "CLHEP/Random/RandFlat.h"
#include "SimCalorimetry/HGCalSimProducers/interface/HGCDigitizerTypes.h"

/**
Expand Down
4 changes: 4 additions & 0 deletions SimCalorimetry/HGCalSimProducers/python/hgcalDigitizer_cfi.py
Expand Up @@ -66,6 +66,7 @@
accumulatorType = cms.string("HGCDigiProducer"),
hitCollection = cms.string("HGCHitsEE"),
digiCollection = cms.string("HGCDigisEE"),
seedOffset = cms.uint32(9182),
maxSimHitsAccTime = cms.uint32(100),
bxTime = cms.double(25),
eVPerEleHolePair = cms.double(eV_per_eh_pair),
Expand Down Expand Up @@ -135,6 +136,7 @@
accumulatorType = cms.string("HGCDigiProducer"),
hitCollection = cms.string("HGCHitsHEfront"),
digiCollection = cms.string("HGCDigisHEfront"),
seedOffset = cms.uint32(8273),
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

since these parameters are no longer used, they can be deleted

maxSimHitsAccTime = cms.uint32(100),
bxTime = cms.double(25),
tofDelay = cms.double(5),
Expand Down Expand Up @@ -203,6 +205,7 @@
accumulatorType = cms.string("HGCDigiProducer"),
hitCollection = cms.string("HcalHits"),
digiCollection = cms.string("HGCDigisHEback"),
seedOffset = cms.uint32(7364),
maxSimHitsAccTime = cms.uint32(100),
bxTime = cms.double(25),
tofDelay = cms.double(1),
Expand Down Expand Up @@ -246,6 +249,7 @@
accumulatorType = cms.string("HGCDigiProducer"),
hitCollection = cms.string("HFNoseHits"),
digiCollection = cms.string("HFNoseDigis"),
seedOffset = cms.uint32(5647),
maxSimHitsAccTime = cms.uint32(100),
bxTime = cms.double(25),
eVPerEleHolePair = cms.double(eV_per_eh_pair),
Expand Down
24 changes: 20 additions & 4 deletions SimCalorimetry/HGCalSimProducers/src/HGCDigitizerBase.cc
Expand Up @@ -6,12 +6,12 @@ using namespace hgc_digi;
using namespace hgc_digi_utils;

template <class DFr>
HGCDigitizerBase<DFr>::HGCDigitizerBase(const edm::ParameterSet& ps) : scaleByDose_(false) {
HGCDigitizerBase<DFr>::HGCDigitizerBase(const edm::ParameterSet& ps)
: scaleByDose_(false), NoiseMean_(0.0), NoiseStd_(1.0) {
bxTime_ = ps.getParameter<double>("bxTime");
myCfg_ = ps.getParameter<edm::ParameterSet>("digiCfg");
doTimeSamples_ = myCfg_.getParameter<bool>("doTimeSamples");
thresholdFollowsMIP_ = myCfg_.getParameter<bool>("thresholdFollowsMIP");

if (myCfg_.exists("keV2fC"))
keV2fC_ = myCfg_.getParameter<double>("keV2fC");
else
Expand Down Expand Up @@ -53,6 +53,17 @@ HGCDigitizerBase<DFr>::HGCDigitizerBase(const edm::ParameterSet& ps) : scaleByDo
myFEelectronics_->SetNoiseValues(noise_fC_);
}

template <class DFr>
void HGCDigitizerBase<DFr>::GenerateGaussianNoise(CLHEP::HepRandomEngine* engine,
const double NoiseMean,
const double NoiseStd) {
for (size_t i = 0; i < NoiseArrayLength_; i++) {
for (size_t j = 0; j < samplesize_; j++) {
GaussianNoiseArray_[i][j] = CLHEP::RandGaussQ::shoot(NoiseMean, NoiseStd);
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

here you need to use CLHEP::RandGaussQ::shoot(engine, NoiseMean, NoiseStd);, otherwise the engine is not actually used.

}
}
}

template <class DFr>
void HGCDigitizerBase<DFr>::run(std::unique_ptr<HGCDigitizerBase::DColl>& digiColl,
HGCSimHitDataAccumulator& simData,
Expand Down Expand Up @@ -80,13 +91,17 @@ void HGCDigitizerBase<DFr>::runSimple(std::unique_ptr<HGCDigitizerBase::DColl>&
HGCCellInfo zeroData;
zeroData.hit_info[0].fill(0.f); //accumulated energy
zeroData.hit_info[1].fill(0.f); //time-of-flight
GenerateGaussianNoise(engine, NoiseMean_, NoiseStd_);

for (const auto& id : validIds) {
chargeColl.fill(0.f);
toa.fill(0.f);
HGCSimHitDataAccumulator::iterator it = simData.find(id);
HGCCellInfo& cell = (simData.end() == it ? zeroData : it->second);
addCellMetadata(cell, theGeom, id);
size_t hash_index = (CLHEP::RandFlat::shootInt(engine, (NoiseArrayLength_ - 1)) + id) % NoiseArrayLength_;

auto cellNoiseArray = GaussianNoiseArray_[hash_index];

//set the noise,cce, LSB and threshold to be used
float cce(1.f), noiseWidth(0.f), lsbADC(-1.f), maxADC(-1.f);
Expand Down Expand Up @@ -122,7 +137,7 @@ void HGCDigitizerBase<DFr>::runSimple(std::unique_ptr<HGCDigitizerBase::DColl>&
toa[i] = cell.hit_info[1][i] / rawCharge;

//final charge estimation
float noise = CLHEP::RandGaussQ::shoot(engine, 0.0, noiseWidth);
float noise = (float)cellNoiseArray[i] * noiseWidth;
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

it might be worth putting a switch here to go back to the old behavior, in case someone wants to run a test in the future.

float totalCharge(rawCharge * cce + noise);
if (totalCharge < 0.f)
totalCharge = 0.f;
Expand All @@ -131,7 +146,8 @@ void HGCDigitizerBase<DFr>::runSimple(std::unique_ptr<HGCDigitizerBase::DColl>&

//run the shaper to create a new data frame
DFr rawDataFrame(id);
myFEelectronics_->runShaper(rawDataFrame, chargeColl, toa, engine, thrADC, lsbADC, maxADC, cell.thickness);
int thickness = cell.thickness > 0 ? cell.thickness : 1;
myFEelectronics_->runShaper(rawDataFrame, chargeColl, toa, engine, thrADC, lsbADC, maxADC, thickness);

//update the output according to the final shape
updateOutput(coll, rawDataFrame);
Expand Down