CMS 3D CMS Logo

HcalSimParameters.cc
Go to the documentation of this file.
9 #include "CLHEP/Random/RandGaussQ.h"
10 using namespace std;
11 
13  double samplingFactor, double timePhase,
15  bool doPhotostatistics, bool syncPhase,
16  int firstRing, const std::vector<double> & samplingFactors,
17  double sipmTau
18  )
19 : CaloSimParameters(simHitToPhotoelectrons, 0.0, samplingFactor, timePhase,
20  readoutFrameSize, binOfMaximum, doPhotostatistics, syncPhase),
21  theDbService(nullptr),
22  theSiPMcharacteristics(nullptr),
23  theFirstRing(firstRing),
24  theSamplingFactors(samplingFactors),
25  theSiPMSmearing(false),
26  doTimeSmear_(true),
27  theSiPMTau(sipmTau)
28 {
30 
31  edm::LogInfo("HcalSimParameters:") << " doSiPMsmearing = " << theSiPMSmearing;
32 }
33 
37  theFirstRing( p.getParameter<int>("firstRing") ),
38  theSamplingFactors( p.getParameter<std::vector<double> >("samplingFactors") ),
39  theSiPMSmearing( p.getParameter<bool>("doSiPMSmearing") ),
40  doTimeSmear_( p.getParameter<bool>("timeSmearing") ),
41  theSiPMTau( p.getParameter<double>("sipmTau") )
42 {
44 
45  edm::LogInfo("HcalSimParameters:") << " doSiPMsmearing = " << theSiPMSmearing;
46 }
47 
49 {
50  assert(service);
51  theDbService = service;
53  assert(theSiPMcharacteristics);
54 }
55 
57 {
58  // the gain is in units of GeV/fC. We want a constant with pe/dGeV
59  // pe/dGeV = (GeV/dGeV) / (GeV/fC) / (fC/pe)
61  if(HcalGenericDetId(detId).genericSubdet() != HcalGenericDetId::HcalGenForward
62  || HcalGenericDetId(detId).genericSubdet() != HcalGenericDetId::HcalGenZDC)
63  {
64  result = samplingFactor(detId) / fCtoGeV(detId) / photoelectronsToAnalog(detId);
65  }
66  return result;
67 }
68 
69 
70 double HcalSimParameters::fCtoGeV(const DetId & detId) const
71 {
72  assert(theDbService != nullptr);
73  HcalGenericDetId hcalGenDetId(detId);
74  const HcalGain* gains = theDbService->getGain(hcalGenDetId);
75  const HcalGainWidth* gwidths = theDbService->getGainWidth(hcalGenDetId);
76  double result = 0.0;
77  if (!gains || !gwidths )
78  {
79  edm::LogError("HcalAmplifier") << "Could not fetch HCAL conditions for channel " << hcalGenDetId;
80  }
81  else
82  {
83  // only one gain will be recorded per channel, so just use capID 0 for now
84  result = gains->getValue(0);
85 // if(doNoise_)
87 // result += CLHEP::RandGaussQ::shoot(0., gwidths->getValue(0));
88 // }
89  }
90  return result;
91 }
92 
93 double HcalSimParameters::samplingFactor(const DetId & detId) const
94 {
95  HcalDetId hcalDetId(detId);
96  return theSamplingFactors.at(hcalDetId.ietaAbs()-theFirstRing);
97 }
98 
99 
101 {
102  //now always taken from database for HPDs or SiPMs (HB, HE, HO)
103  assert(theDbService);
104  return theDbService->getHcalSiPMParameter(detId)->getFCByPE();
105 }
106 
107 
108 //static const double GeV2fC = 1.0/0.145;
109 static const double GeV2fC = 1.0/0.4;
110 
112  // GeV->ampl (fC), time (ns)
113  theSmearSettings.emplace_back( 4.00*GeV2fC, 4.050);
114  theSmearSettings.emplace_back( 20.00*GeV2fC, 3.300);
115  theSmearSettings.emplace_back( 25.00*GeV2fC, 2.925);
116  theSmearSettings.emplace_back( 30.00*GeV2fC, 2.714);
117  theSmearSettings.emplace_back( 37.00*GeV2fC, 2.496);
118  theSmearSettings.emplace_back( 44.50*GeV2fC, 2.278);
119  theSmearSettings.emplace_back( 56.00*GeV2fC, 2.138);
120  theSmearSettings.emplace_back( 63.50*GeV2fC, 2.022);
121  theSmearSettings.emplace_back( 81.00*GeV2fC, 1.788);
122  theSmearSettings.emplace_back( 88.50*GeV2fC, 1.695);
123  theSmearSettings.emplace_back(114.50*GeV2fC, 1.716);
124  theSmearSettings.emplace_back(175.50*GeV2fC, 1.070);
125  theSmearSettings.emplace_back(350.00*GeV2fC, 1.564);
126  theSmearSettings.emplace_back(99999.00*GeV2fC, 1.564);
127 }
128 
129 double HcalSimParameters::timeSmearRMS(double ampl) const {
131  double smearsigma=0;
132 
133  for (i=0; i<theSmearSettings.size(); i++)
134  if (theSmearSettings[i].first > ampl)
135  break;
136 
137  // Smearing occurs only within the envelope definitions.
138  if (i!=0 && (i < theSmearSettings.size())) {
139  double energy1 = theSmearSettings[i-1].first;
140  double sigma1 = theSmearSettings[i-1].second;
141  double energy2 = theSmearSettings[i].first;
142  double sigma2 = theSmearSettings[i].second;
143 
144  if (energy2 != energy1)
145  smearsigma = sigma1 + ((sigma2-sigma1)*(ampl-energy1)/(energy2-energy1));
146  else
147  smearsigma = (sigma2+sigma1)/2.;
148  }
149 
150  return smearsigma;
151 
152 }
153 
154 int HcalSimParameters::pixels(const DetId & detId) const {
155  assert(theDbService);
157  return theSiPMcharacteristics->getPixels(type);
158 }
159 
160 double HcalSimParameters::sipmDarkCurrentuA(const DetId & detId) const
161 {
162  assert(theDbService);
164 }
165 
166 double HcalSimParameters::sipmCrossTalk(const DetId & detId) const
167 {
168  assert(theDbService);
170  return theSiPMcharacteristics->getCrossTalk(type);
171 }
172 std::vector<float> HcalSimParameters::sipmNonlinearity(const DetId & detId) const
173 {
174  assert(theDbService);
177 }
178 
179 unsigned int HcalSimParameters::signalShape(const DetId & detId) const {
180  assert(theDbService);
181  return theDbService->getHcalMCParam(detId)->signalShape();
182 }
std::vector< float > sipmNonlinearity(const DetId &detId) const
const HcalSiPMCharacteristics * theSiPMcharacteristics
type
Definition: HCALResponse.h:21
const HcalGainWidth * getGainWidth(const HcalGenericDetId &fId) const
double fCtoGeV(const DetId &detId) const
int getType() const
get SiPM type
int pixels(const DetId &detId) const
#define nullptr
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.
float getCrossTalk(int type) const
get cross talk
static const double GeV2fC
HcalTimeSmearSettings theSmearSettings
unsigned int signalShape(const DetId &detId) const
double timeSmearRMS(double ampl) const
double simHitToPhotoelectrons() const
const HcalMCParam * getHcalMCParam(const HcalGenericDetId &fId) const
std::vector< double > theSamplingFactors
virtual double samplingFactor(const DetId &detId) const
float getDarkCurrent() const
get dark current
int ietaAbs() const
get the absolute value of the cell ieta
Definition: HcalDetId.cc:119
Definition: DetId.h:18
HcalSimParameters(double simHitToPhotoelectrons, double samplingFactor, double timePhase, int readoutFrameSize, int binOfMaximum, bool doPhotostatistics, bool syncPhase, int firstRing, const std::vector< double > &samplingFactors, double sipmTau)
double sipmDarkCurrentuA(const DetId &detId) const
const HcalGain * getGain(const HcalGenericDetId &fId) const
const HcalSiPMCharacteristics * getHcalSiPMCharacteristics() const
std::vector< float > getNonLinearities(int type) const
get nonlinearity constants
unsigned int signalShape() const
Definition: HcalMCParam.h:40
void setDbService(const HcalDbService *service)
double sipmCrossTalk(const DetId &detId) const
const HcalDbService * theDbService
int getPixels(int type) const
get # of pixels
const HcalSiPMParameter * getHcalSiPMParameter(const HcalGenericDetId &fId) const
float getFCByPE() const
get fcByPE
double photoelectronsToAnalog() const
the factor which goes from photoelectrons to whatever gets read by ADCs