CMS 3D CMS Logo

CastorSimParameters.cc
Go to the documentation of this file.
1 #include "CLHEP/Random/RandGaussQ.h"
9 #include <cassert>
10 
13  double samplingFactor,
14  double timePhase,
15  bool syncPhase)
18  theDbService(nullptr),
19  theSamplingFactor(samplingFactor),
20  nominalfCperPE(1) {}
21 
24  theDbService(nullptr),
25  theSamplingFactor(p.getParameter<double>("samplingFactor")),
26  nominalfCperPE(p.getParameter<double>("photoelectronsToAnalog")) {}
27 
29  // return the nominal PMT gain value of CASTOR from the config file.
30  return nominalfCperPE;
31 }
32 
34  // calculate factor (PMT gain) using sampling factor value & available
35  // electron gain
37 }
38 
40  assert(theDbService != nullptr);
41  HcalGenericDetId hcalGenDetId(detId);
42  const CastorGain *gains = theDbService->getGain(hcalGenDetId);
43  const CastorGainWidth *gwidths = theDbService->getGainWidth(hcalGenDetId);
44  double result = 0.0;
45  if (!gains || !gwidths) {
46  edm::LogError("CastorAmplifier") << "Could not fetch HCAL conditions for channel " << hcalGenDetId;
47  } else {
48  // only one gain will be recorded per channel, so just use capID 0 for now
49  result = gains->getValue(0);
50  }
51  return result;
52 }
double getNominalfCperPE() const
double fCtoGeV(const DetId &detId) const
const CastorDbService * theDbService
double photoelectronsToAnalog() const
the factor which goes from photoelectrons to whatever gets read by ADCs
const CastorGainWidth * getGainWidth(const HcalGenericDetId &fId) const
Log< level::Error, false > LogError
assert(be >=bs)
Main class for Parameters in different subdetectors.
CastorSimParameters(double simHitToPhotoelectrons, double photoelectronsToAnalog, double samplingFactor, double timePhase, bool syncPhase)
const CastorGain * getGain(const HcalGenericDetId &fId) const
constexpr float gains[NGAINS]
Definition: EcalConstants.h:20
Definition: DetId.h:17