CMS 3D CMS Logo

HFSimParameters.cc
Go to the documentation of this file.
8 #include "CLHEP/Random/RandGaussQ.h"
9 
12  double samplingFactor,
13  double timePhase,
14  bool syncPhase)
17  theDbService(nullptr),
18  theSamplingFactor(samplingFactor) {}
19 
22  theDbService(nullptr),
23  theSamplingFactor(p.getParameter<double>("samplingFactor")),
24  threshold_currentTDC_(p.getParameter<double>("threshold_currentTDC")) {}
25 
27  // pe/fC = pe/GeV * GeV/fC = (0.24 pe/GeV * 6 for photomult * 0.2146GeV/fC)
29 }
30 
31 double HFSimParameters::fCtoGeV(const DetId& detId) const {
32  assert(theDbService != nullptr);
33  HcalGenericDetId hcalGenDetId(detId);
34  const HcalGain* gains = theDbService->getGain(hcalGenDetId);
35  const HcalGainWidth* gwidths = theDbService->getGainWidth(hcalGenDetId);
36  double result = 0.0;
37  if (!gains || !gwidths) {
38  edm::LogError("HcalAmplifier") << "Could not fetch HCAL conditions for channel " << hcalGenDetId;
39  } else {
40  // only one gain will be recorded per channel, so just use capID 0 for now
41  result = gains->getValue(0);
42  // if(doNoise_)
44  // result += CLHEP::RandGaussQ::shoot(0., gwidths->getValue(0));
45  // }
46  }
47  return result;
48 }
49 
const HcalGainWidth * getGainWidth(const HcalGenericDetId &fId) const
double photoelectronsToAnalog() const
the factor which goes from photoelectrons to whatever gets read by ADCs
Log< level::Error, false > LogError
assert(be >=bs)
Main class for Parameters in different subdetectors.
double simHitToPhotoelectrons() const
double samplingFactor() const
constexpr float gains[NGAINS]
Definition: EcalConstants.h:20
HFSimParameters(double simHitToPhotoelectrons, double photoelectronsToAnalog, double samplingFactor, double timePhase, bool syncPhase)
Definition: DetId.h:17
const HcalDbService * theDbService
double theSamplingFactor
const HcalGain * getGain(const HcalGenericDetId &fId) const
double fCtoGeV(const DetId &detId) const