CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HPDIonFeedbackSim.cc
Go to the documentation of this file.
1 
2 // --------------------------------------------------------
3 // A class to simulated HPD ion feedback noise.
4 // The deliverable of the class is the ion feedback noise
5 // for an HcalDetId units of fC or GeV
6 //
7 // Project: HPD ion feedback
8 // Author: T.Yetkin University of Iowa, Feb. 16, 2010
9 // $Id: HPDIonFeedbackSim.cc,v 1.6 2011/11/14 11:08:16 abdullin Exp $
10 // --------------------------------------------------------
11 
24 
25 
26 using namespace edm;
27 using namespace std;
28 
29 // constants for simulation/parameterization
30 double pe2Charge = 0.333333; // fC/p.e.
31 
33 : theDbService(0), theShapes(shapes),
34  theRandBinomial(0), theRandFlat(0), theRandGauss(0), theRandPoissonQ(0)
35 {
36 }
37 
39 {
41  if (theRandFlat) delete theRandFlat;
42  if (theRandGauss) delete theRandGauss;
43  if (theRandPoissonQ) delete theRandPoissonQ;
44 }
45 
46 
47 void HPDIonFeedbackSim::setRandomEngine(CLHEP::HepRandomEngine & engine)
48 {
50  if (theRandFlat) delete theRandFlat;
51  if (theRandGauss) delete theRandGauss;
52  if (theRandPoissonQ) delete theRandPoissonQ;
53 
54  theRandBinomial = new CLHEP::RandBinomial(engine);
55  theRandFlat = new CLHEP::RandFlat(engine);
56  theRandGauss = new CLHEP::RandGaussQ(engine);
57  theRandPoissonQ = new CLHEP::RandPoissonQ(engine);
58 }
59 
60 double HPDIonFeedbackSim::getIonFeedback(DetId detId, double signal, double pedWidth, bool doThermal, bool isInGeV)
61 {
62 
63  // HcalDetId id = detId;
64 
65  double GeVperfC = 1.;
66  if(isInGeV) GeVperfC = 1./fCtoGeV(detId);
67 
68  double charge = signal / GeVperfC;
69 
70  double noise = 0.; // fC
71  if (charge > 3. * pedWidth) { // 3 sigma away from pedestal mean
72  int npe = int (charge / pe2Charge);
73  if(doThermal) {
74  double electronEmission = 0.08;
75  npe += theRandPoissonQ->fire(electronEmission);
76  }
77 
78  noise = correctPE(detId, npe) - npe;
79  }
80  return (noise * GeVperfC);
81 
82 }
83 
84 
85 double HPDIonFeedbackSim::correctPE(const DetId & detId, double npe) const
86 {
87  double rateInTail = 0.000211988;//read this from XML file
88  double rateInSecondTail = 4.61579e-06;//read this from XML file
89 
90  // three gauss fit is applied to data to get ion feedback distribution
91  // parameters (in fC)
92  // first gaussian
93  // double p0 = 9.53192e+05;
94  // double p1 = -3.13653e-01;
95  // double p2 = 2.78350e+00;
96 
97  // second gaussian
98  // double p3 = 2.41611e+03;
99  double p4 = 2.06117e+01;
100  double p5 = 1.09239e+01;
101 
102  // third gaussian
103  // double p6 = 3.42793e+01;
104  double p7 = 5.45548e+01;
105  double p8 = 1.59696e+01;
106 
107  double noise = 0.; // fC
108  int nFirst = (int)(theRandBinomial->fire(npe, rateInTail));
109  int nSecond = (int)(theRandBinomial->fire(npe, rateInSecondTail));
110 
111  for (int j = 0; j < nFirst; ++j) {
112  noise += theRandGauss->fire(p4, p5);
113  }
114  for (int j = 0; j < nSecond; ++j) {
115  noise += theRandGauss->fire(p7, p8);
116  }
117 
118  return npe + std::max(noise/pe2Charge, 0.);
119 }
120 
121 
123 {
124  // make some chance to add a PE (with a chance of feedback)
125  // for each time sample
126  double meanPE = 0.02;
127  DetId detId(samples.id());
128  int nSamples = samples.size();
129  const CaloVShape * shape = theShapes->shape(detId);
130  for(int i = 0; i < nSamples; ++i)
131  {
132  double npe = theRandPoissonQ->fire(meanPE);
133  // TODOprobably should time-smear these
134  if(npe > 0.)
135  {
136  // chance of feedback
137  npe = correctPE(detId, npe);
138  for(int j = i; j < nSamples; ++j)
139  {
140  double timeFromPE = (j-i) * 25.;
141  samples[j] += (*shape)(timeFromPE) * npe;
142  }
143  }
144  }
145 }
146 
147 
148 double HPDIonFeedbackSim::fCtoGeV(const DetId & detId) const
149 {
150  assert(theDbService != 0);
151  HcalGenericDetId hcalGenDetId(detId);
152  const HcalGain* gains = theDbService->getGain(hcalGenDetId);
153  const HcalGainWidth* gwidths = theDbService->getGainWidth(hcalGenDetId);
154  if (!gains || !gwidths )
155  {
156  edm::LogError("HcalAmplifier") << "Could not fetch HCAL conditions for channel " << hcalGenDetId;
157  }
158  // only one gain will be recorded per channel, so just use capID 0 for now
159  double result = gains->getValue(0);
160  return result;
161 }
HPDIonFeedbackSim(const edm::ParameterSet &, const CaloShapes *shapes)
need a shaper in order to set thermal noise
const HcalGainWidth * getGainWidth(const HcalGenericDetId &fId) const
int i
Definition: DBlmapReader.cc:9
void addThermalNoise(CaloSamples &samples)
CLHEP::RandGaussQ * theRandGauss
double fCtoGeV(const DetId &detId) const
CLHEP::RandFlat * theRandFlat
Electronic response of the preamp.
Definition: CaloVShape.h:11
float getValue(int fCapId) const
get value for capId = 0..3
Definition: HcalGain.h:20
double pe2Charge
double charge(const std::vector< uint8_t > &Ampls)
const HcalDbService * theDbService
virtual double correctPE(const DetId &detId, double npe) const
const T & max(const T &a, const T &b)
void setRandomEngine(CLHEP::HepRandomEngine &engine)
need a shaper in order to set thermal noise
double p4[4]
Definition: TauolaWrapper.h:92
tuple result
Definition: query.py:137
int j
Definition: DBlmapReader.cc:9
CLHEP::RandBinomial * theRandBinomial
Definition: DetId.h:20
CLHEP::RandPoissonQ * theRandPoissonQ
int size() const
get the size
Definition: CaloSamples.h:24
const HcalGain * getGain(const HcalGenericDetId &fId) const
const CaloShapes * theShapes
DetId id() const
get the (generic) id
Definition: CaloSamples.h:21
double getIonFeedback(DetId detId, double signal, double pedWidth, bool doThermal, bool isInGeV)
virtual const CaloVShape * shape(const DetId &detId) const
Definition: CaloShapes.h:15