CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CastorSimParameters.cc
Go to the documentation of this file.
8 #include "CLHEP/Random/RandGaussQ.h"
9 
10 
11 CastorSimParameters::CastorSimParameters(double simHitToPhotoelectrons, double photoelectronsToAnalog,double samplingFactor, double timePhase, bool syncPhase)
12 : CaloSimParameters(simHitToPhotoelectrons, photoelectronsToAnalog, samplingFactor, timePhase, 6, 4, false, syncPhase),
13  theDbService(0),
14  theSamplingFactor( samplingFactor ),
15  nominalfCperPE( 1),
16  dynamicNoise(false)
17 {
18 }
19 
20 
23  theDbService(0),
24  theSamplingFactor( p.getParameter<double>("samplingFactor") ),
25  nominalfCperPE( p.getParameter<double>("photoelectronsToAnalog") ),
26  dynamicNoise(p.getParameter<bool>("doDynamicNoise") )
27 {
28 }
29 
31 {
32  // return the nominal PMT gain value of CASTOR from the config file.
33  return nominalfCperPE;
34 }
35 
37 {
38  // activate proper noise treatment depending on used gain values.
39  return dynamicNoise;
40 }
41 
43 {
44  // calculate factor (PMT gain) using sampling factor value & available electron gain
45  return theSamplingFactor/fCtoGeV(detId);
46 }
47 
48 
49 
50 double CastorSimParameters::fCtoGeV(const DetId & detId) const
51 {
52  assert(theDbService != nullptr);
53  HcalGenericDetId hcalGenDetId(detId);
54  const CastorGain* gains = theDbService->getGain(hcalGenDetId);
55  const CastorGainWidth* gwidths = theDbService->getGainWidth(hcalGenDetId);
56  double result = 0.0;
57  if (!gains || !gwidths )
58  {
59  edm::LogError("CastorAmplifier") << "Could not fetch HCAL conditions for channel " << hcalGenDetId;
60  }
61  else
62  {
63  // only one gain will be recorded per channel, so just use capID 0 for now
64  result = gains->getValue(0);
65  }
66  return result;
67 }
bool doDynamicNoise() const
double getNominalfCperPE() const
const CastorDbService * theDbService
assert(m_qm.get())
Main class for Parameters in different subdetectors.
CastorSimParameters(double simHitToPhotoelectrons, double photoelectronsToAnalog, double samplingFactor, double timePhase, bool syncPhase)
tuple result
Definition: mps_fire.py:84
const CastorGain * getGain(const HcalGenericDetId &fId) const
double fCtoGeV(const DetId &detId) const
Definition: DetId.h:18
const CastorGainWidth * getGainWidth(const HcalGenericDetId &fId) const
volatile std::atomic< bool > shutdown_flag false
float getValue(int fCapId) const
get value for capId = 0..3
Definition: CastorGain.h:19
double photoelectronsToAnalog() const
the factor which goes from photoelectrons to whatever gets read by ADCs