File indexing completed on 2024-04-06 12:29:21
0001 #include "CalibFormats/CaloObjects/interface/CaloSamples.h"
0002 #include "CalibFormats/CastorObjects/interface/CastorDbService.h"
0003 #include "CondFormats/CastorObjects/interface/CastorGain.h"
0004 #include "CondFormats/CastorObjects/interface/CastorGainWidth.h"
0005 #include "CondFormats/CastorObjects/interface/CastorPedestal.h"
0006 #include "CondFormats/CastorObjects/interface/CastorPedestalWidth.h"
0007 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
0008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0009 #include "SimCalorimetry/CastorSim/interface/CastorAmplifier.h"
0010 #include "SimCalorimetry/CastorSim/interface/CastorSimParameters.h"
0011
0012 #include "CLHEP/Random/RandGaussQ.h"
0013
0014 #include <iostream>
0015 #include <cassert>
0016
0017 CastorAmplifier::CastorAmplifier(const CastorSimParameterMap *parameters, bool addNoise)
0018 : theDbService(nullptr), theParameterMap(parameters), theStartingCapId(0), addNoise_(addNoise) {}
0019
0020 void CastorAmplifier::amplify(CaloSamples &frame, CLHEP::HepRandomEngine *engine) const {
0021 const CastorSimParameters ¶meters = theParameterMap->castorParameters();
0022 assert(theDbService != nullptr);
0023 HcalGenericDetId hcalGenDetId(frame.id());
0024 const CastorPedestal *peds = theDbService->getPedestal(hcalGenDetId);
0025 const CastorPedestalWidth *pwidths = theDbService->getPedestalWidth(hcalGenDetId);
0026 if (!peds || !pwidths) {
0027 edm::LogError("CastorAmplifier") << "Could not fetch HCAL/CASTOR conditions for channel " << hcalGenDetId;
0028 } else {
0029 double gauss[32];
0030 double noise[32];
0031 double fCperPE = parameters.photoelectronsToAnalog(frame.id());
0032 double nominalfCperPE = parameters.getNominalfCperPE();
0033
0034 for (int i = 0; i < frame.size(); i++) {
0035 gauss[i] = CLHEP::RandGaussQ::shoot(engine, 0., 1.);
0036 }
0037 if (addNoise_) {
0038 pwidths->makeNoise(frame.size(), gauss, noise);
0039 }
0040 for (int tbin = 0; tbin < frame.size(); ++tbin) {
0041 int capId = (theStartingCapId + tbin) % 4;
0042 double pedestal = peds->getValue(capId);
0043 if (addNoise_) {
0044 pedestal += noise[tbin] * (fCperPE / nominalfCperPE);
0045 }
0046 frame[tbin] *= fCperPE;
0047 frame[tbin] += pedestal;
0048 }
0049 }
0050 LogDebug("CastorAmplifier") << frame;
0051 }