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 {
16 }
17 
18 /*
19 CastorSimParameters::CastorSimParameters(double simHitToPhotoelectrons, double photoelectronsToAnalog,
20  double samplingFactor, double timePhase,
21  int readoutFrameSize, int binOfMaximum,
22  bool doPhotostatistics, bool syncPhase,
23  int firstRing, const std::vector<double> & samplingFactors)
24  : CaloSimParameters(simHitToPhotoelectrons, photoelectronsToAnalog, samplingFactor, timePhase,
25  readoutFrameSize, binOfMaximum, doPhotostatistics, syncPhase),
26  theDbService(0),
27  theFirstRing(firstRing),
28  theSamplingFactors(samplingFactors)
29 {
30 }
31 */
32 
35  theDbService(0),
36  theSamplingFactor( p.getParameter<double>("samplingFactor") )
37 {
38 }
39 /*
40 : CaloSimParameters(p),
41  theDbService(0),
42  theFirstRing( p.getParameter<int>("firstRing") ),
43  theSamplingFactors( p.getParameter<std::vector<double> >("samplingFactors") )
44 {
45 }
46 */
47 
48 /*
49 double CastorSimParameters::simHitToPhotoelectrons(const DetId & detId) const
50 {
51  // the gain is in units of GeV/fC. We want a constant with pe/dGeV
52  // pe/dGeV = (GeV/dGeV) / (GeV/fC) / (fC/pe)
53  double result = CaloSimParameters::simHitToPhotoelectrons(detId);
54  if(HcalGenericDetId(detId).genericSubdet() != HcalGenericDetId::HcalGenForward
55  || HcalGenericDetId(detId).genericSubdet() != HcalGenericDetId::HcalGenCastor)
56  {
57  result = samplingFactor(detId) / fCtoGeV(detId) / photoelectronsToAnalog();
58  }
59  return result;
60 }
61 */
62 
64 {
65  // pe/fC = pe/GeV * GeV/fC = (0.24 pe/GeV * 6 for photomult * 0.2146GeV/fC)
66  return 1./(theSamplingFactor * simHitToPhotoelectrons(detId) * fCtoGeV(detId));
67 }
68 
69 
70 
71 double CastorSimParameters::fCtoGeV(const DetId & detId) const
72 {
73  assert(theDbService != 0);
74  HcalGenericDetId hcalGenDetId(detId);
75  const CastorGain* gains = theDbService->getGain(hcalGenDetId);
76  const CastorGainWidth* gwidths = theDbService->getGainWidth(hcalGenDetId);
77  if (!gains || !gwidths )
78  {
79  edm::LogError("CastorAmplifier") << "Could not fetch HCAL conditions for channel " << hcalGenDetId;
80  }
81  // only one gain will be recorded per channel, so just use capID 0 for now
82 
83  double result = gains->getValue(0);
84 // if(doNoise_)
86 // result += CLHEP::RandGaussQ::shoot(0., gwidths->getValue(0));
87 // }
88  return result;
89 }
90 /*
91 double CastorSimParameters::samplingFactor(const DetId & detId) const {
92 HcalDetId hcalDetId(detId);
93 return theSamplingFactors[hcalDetId.ietaAbs()-theFirstRing];
94 }
95 */
const CastorDbService * theDbService
Main class for Parameters in different subdetectors.
CastorSimParameters(double simHitToPhotoelectrons, double photoelectronsToAnalog, double samplingFactor, double timePhase, bool syncPhase)
tuple result
Definition: query.py:137
double simHitToPhotoelectrons() const
const CastorGain * getGain(const HcalGenericDetId &fId) const
double fCtoGeV(const DetId &detId) const
Definition: DetId.h:20
const CastorGainWidth * getGainWidth(const HcalGenericDetId &fId) const
float getValue(int fCapId) const
get value for capId = 0..3
Definition: CastorGain.h:17
double photoelectronsToAnalog() const
the factor which goes from photoelectrons to whatever gets read by ADCs