CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HFSimParameters.cc
Go to the documentation of this file.
8 #include "CLHEP/Random/RandGaussQ.h"
9 
10 
11 HFSimParameters::HFSimParameters(double simHitToPhotoelectrons, double photoelectronsToAnalog,
12  double samplingFactor, double timePhase, bool syncPhase)
13 : CaloSimParameters(simHitToPhotoelectrons, photoelectronsToAnalog, samplingFactor, timePhase,
14  6, 4, false, syncPhase),
15  theDbService(0),
16  theSamplingFactor( samplingFactor )
17 {
18 }
19 
22  theDbService(0),
23  theSamplingFactor( p.getParameter<double>("samplingFactor") )
24 {
25 }
26 
27 
29 {
30  // pe/fC = pe/GeV * GeV/fC = (0.24 pe/GeV * 6 for photomult * 0.2146GeV/fC)
31  return 1./(theSamplingFactor * simHitToPhotoelectrons(detId) * fCtoGeV(detId));
32 }
33 
34 double HFSimParameters::fCtoGeV(const DetId & detId) const
35 {
36  assert(theDbService != 0);
37  HcalGenericDetId hcalGenDetId(detId);
38  const HcalGain* gains = theDbService->getGain(hcalGenDetId);
39  const HcalGainWidth* gwidths = theDbService->getGainWidth(hcalGenDetId);
40  if (!gains || !gwidths )
41  {
42  edm::LogError("HcalAmplifier") << "Could not fetch HCAL conditions for channel " << hcalGenDetId;
43  }
44  // only one gain will be recorded per channel, so just use capID 0 for now
45 
46  double result = gains->getValue(0);
47 // if(doNoise_)
49 // result += CLHEP::RandGaussQ::shoot(0., gwidths->getValue(0));
50 // }
51  return result;
52 }
53 
54 
const HcalGainWidth * getGainWidth(const HcalGenericDetId &fId) const
assert(m_qm.get())
double fCtoGeV(const DetId &detId) const
float getValue(int fCapId) const
get value for capId = 0..3
Definition: HcalGain.h:22
Main class for Parameters in different subdetectors.
tuple result
Definition: query.py:137
double simHitToPhotoelectrons() const
HFSimParameters(double simHitToPhotoelectrons, double photoelectronsToAnalog, double samplingFactor, double timePhase, bool syncPhase)
Definition: DetId.h:18
const HcalDbService * theDbService
double theSamplingFactor
const HcalGain * getGain(const HcalGenericDetId &fId) const
volatile std::atomic< bool > shutdown_flag false
double photoelectronsToAnalog() const
the factor which goes from photoelectrons to whatever gets read by ADCs