CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HcalSimParameters.cc
Go to the documentation of this file.
8 #include "CLHEP/Random/RandGaussQ.h"
9 using namespace std;
10 
11 HcalSimParameters::HcalSimParameters(double simHitToPhotoelectrons, const std::vector<double> & photoelectronsToAnalog,
12  double samplingFactor, double timePhase,
13  int readoutFrameSize, int binOfMaximum,
14  bool doPhotostatistics, bool syncPhase,
15  int firstRing, const std::vector<double> & samplingFactors)
16 : CaloSimParameters(simHitToPhotoelectrons, photoelectronsToAnalog[0], samplingFactor, timePhase,
17  readoutFrameSize, binOfMaximum, doPhotostatistics, syncPhase),
18  theDbService(0),
19  theFirstRing(firstRing),
20  theSamplingFactors(samplingFactors),
21  thePE2fCByRing(photoelectronsToAnalog),
22  thePixels(0),
23  theSiPMSmearing(false),
24  doTimeSmear_(true)
25 {
27 }
28 
31  theDbService(0),
32  theFirstRing( p.getParameter<int>("firstRing") ),
33  theSamplingFactors( p.getParameter<std::vector<double> >("samplingFactors") ),
34  thePE2fCByRing( p.getParameter<std::vector<double> >("photoelectronsToAnalog") ),
35  thePixels(0), theSiPMSmearing(false),
36  doTimeSmear_( p.getParameter<bool>("timeSmearing"))
37 {
38  if(p.exists("pixels"))
39  {
40  thePixels = p.getParameter<int>("pixels");
41  }
42  if (p.exists("doSiPMSmearing"))
43  theSiPMSmearing = p.getParameter<bool>("doSiPMSmearing");
45 }
46 
48 {
49  // the gain is in units of GeV/fC. We want a constant with pe/dGeV
50  // pe/dGeV = (GeV/dGeV) / (GeV/fC) / (fC/pe)
52  if(HcalGenericDetId(detId).genericSubdet() != HcalGenericDetId::HcalGenForward
53  || HcalGenericDetId(detId).genericSubdet() != HcalGenericDetId::HcalGenZDC)
54  {
55  result = samplingFactor(detId) / fCtoGeV(detId) / photoelectronsToAnalog(detId);
56  }
57  return result;
58 }
59 
60 
61 double HcalSimParameters::fCtoGeV(const DetId & detId) const
62 {
63  assert(theDbService != 0);
64  HcalGenericDetId hcalGenDetId(detId);
65  const HcalGain* gains = theDbService->getGain(hcalGenDetId);
66  const HcalGainWidth* gwidths = theDbService->getGainWidth(hcalGenDetId);
67  if (!gains || !gwidths )
68  {
69  edm::LogError("HcalAmplifier") << "Could not fetch HCAL conditions for channel " << hcalGenDetId;
70  }
71  // only one gain will be recorded per channel, so just use capID 0 for now
72  double result = gains->getValue(0);
73 // if(doNoise_)
75 // result += CLHEP::RandGaussQ::shoot(0., gwidths->getValue(0));
76 // }
77  return result;
78 }
79 
80 double HcalSimParameters::samplingFactor(const DetId & detId) const
81 {
82  HcalDetId hcalDetId(detId);
83  return theSamplingFactors.at(hcalDetId.ietaAbs()-theFirstRing);
84 }
85 
86 
88 {
89  HcalDetId hcalDetId(detId);
90  return thePE2fCByRing.at(hcalDetId.ietaAbs()-theFirstRing);
91 }
92 
93 
94 //static const double GeV2fC = 1.0/0.145;
95 static const double GeV2fC = 1.0/0.4;
96 
98  // GeV->ampl (fC), time (ns)
99  theSmearSettings.emplace_back( 4.00*GeV2fC, 4.050);
100  theSmearSettings.emplace_back( 20.00*GeV2fC, 3.300);
101  theSmearSettings.emplace_back( 25.00*GeV2fC, 2.925);
102  theSmearSettings.emplace_back( 30.00*GeV2fC, 2.714);
103  theSmearSettings.emplace_back( 37.00*GeV2fC, 2.496);
104  theSmearSettings.emplace_back( 44.50*GeV2fC, 2.278);
105  theSmearSettings.emplace_back( 56.00*GeV2fC, 2.138);
106  theSmearSettings.emplace_back( 63.50*GeV2fC, 2.022);
107  theSmearSettings.emplace_back( 81.00*GeV2fC, 1.788);
108  theSmearSettings.emplace_back( 88.50*GeV2fC, 1.695);
109  theSmearSettings.emplace_back(114.50*GeV2fC, 1.716);
110  theSmearSettings.emplace_back(175.50*GeV2fC, 1.070);
111  theSmearSettings.emplace_back(350.00*GeV2fC, 1.564);
112  theSmearSettings.emplace_back(99999.00*GeV2fC, 1.564);
113 }
114 
115 double HcalSimParameters::timeSmearRMS(double ampl) const {
117  double smearsigma=0;
118 
119  for (i=0; i<theSmearSettings.size(); i++)
120  if (theSmearSettings[i].first > ampl)
121  break;
122 
123  // Smearing occurs only within the envelope definitions.
124  if (i!=0 && (i < theSmearSettings.size())) {
125  double energy1 = theSmearSettings[i-1].first;
126  double sigma1 = theSmearSettings[i-1].second;
127  double energy2 = theSmearSettings[i].first;
128  double sigma2 = theSmearSettings[i].second;
129 
130  if (energy2 != energy1)
131  smearsigma = sigma1 + ((sigma2-sigma1)*(ampl-energy1)/(energy2-energy1));
132  else
133  smearsigma = (sigma2+sigma1)/2.;
134  }
135 
136  return smearsigma;
137 
138 }
T getParameter(std::string const &) const
const HcalGainWidth * getGainWidth(const HcalGenericDetId &fId) const
int i
Definition: DBlmapReader.cc:9
assert(m_qm.get())
double fCtoGeV(const DetId &detId) const
bool exists(std::string const &parameterName) const
checks if a parameter exists
float getValue(int fCapId) const
get value for capId = 0..3
Definition: HcalGain.h:22
uint16_t size_type
Main class for Parameters in different subdetectors.
tuple result
Definition: mps_fire.py:84
static const double GeV2fC
HcalTimeSmearSettings theSmearSettings
std::vector< double > thePE2fCByRing
double timeSmearRMS(double ampl) const
double simHitToPhotoelectrons() const
std::vector< double > theSamplingFactors
virtual double samplingFactor(const DetId &detId) const
int ietaAbs() const
get the absolute value of the cell ieta
Definition: HcalDetId.cc:119
HcalSimParameters(double simHitToPhotoelectrons, const std::vector< double > &photoelectronsToAnalog, double samplingFactor, double timePhase, int readoutFrameSize, int binOfMaximum, bool doPhotostatistics, bool syncPhase, int firstRing, const std::vector< double > &samplingFactors)
Definition: DetId.h:18
const HcalGain * getGain(const HcalGenericDetId &fId) const
const HcalDbService * theDbService
volatile std::atomic< bool > shutdown_flag false
double photoelectronsToAnalog() const
the factor which goes from photoelectrons to whatever gets read by ADCs