CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
HPDIonFeedbackSim Class Reference

#include <HPDIonFeedbackSim.h>

Inheritance diagram for HPDIonFeedbackSim:
CaloVPECorrection

Public Member Functions

void addThermalNoise (CaloSamples &samples, CLHEP::HepRandomEngine *)
 
double correctPE (const DetId &detId, double npe, CLHEP::HepRandomEngine *) const override
 
double getIonFeedback (DetId detId, double signal, double pedWidth, bool doThermal, bool isInGeV, CLHEP::HepRandomEngine *)
 
 HPDIonFeedbackSim (const edm::ParameterSet &, const CaloShapes *shapes)
 need a shaper in order to set thermal noise More...
 
void setDbService (const HcalDbService *service)
 
 ~HPDIonFeedbackSim () override
 
- Public Member Functions inherited from CaloVPECorrection
virtual ~CaloVPECorrection ()
 

Private Member Functions

double fCtoGeV (const DetId &detId) const
 

Private Attributes

const HcalDbServicetheDbService
 
const CaloShapestheShapes
 

Detailed Description

Definition at line 27 of file HPDIonFeedbackSim.h.

Constructor & Destructor Documentation

HPDIonFeedbackSim::HPDIonFeedbackSim ( const edm::ParameterSet iConfig,
const CaloShapes shapes 
)

need a shaper in order to set thermal noise

Definition at line 32 of file HPDIonFeedbackSim.cc.

33 : theDbService(nullptr), theShapes(shapes)
34 {
35 }
const HcalDbService * theDbService
const CaloShapes * theShapes
HPDIonFeedbackSim::~HPDIonFeedbackSim ( )
override

Definition at line 37 of file HPDIonFeedbackSim.cc.

38 {
39 }

Member Function Documentation

void HPDIonFeedbackSim::addThermalNoise ( CaloSamples samples,
CLHEP::HepRandomEngine *  engine 
)

Definition at line 104 of file HPDIonFeedbackSim.cc.

References correctPE(), mps_fire::i, CaloSamples::id(), hgc_digi::nSamples, CaloShapes::shape(), CaloSamples::size(), and theShapes.

Referenced by HcalAmplifier::amplify().

105 {
106  // make some chance to add a PE (with a chance of feedback)
107  // for each time sample
108  double meanPE = 0.02;
109  DetId detId(samples.id());
110  int nSamples = samples.size();
111  const CaloVShape * shape = theShapes->shape(detId);
112  for(int i = 0; i < nSamples; ++i)
113  {
114  CLHEP::RandPoissonQ theRandPoissonQ(*engine, meanPE);
115  double npe = theRandPoissonQ.fire();
116  // TODOprobably should time-smear these
117  if(npe > 0.)
118  {
119  // chance of feedback
120  npe = correctPE(detId, npe, engine);
121  for(int j = i; j < nSamples; ++j)
122  {
123  double timeFromPE = (j-i) * 25.;
124  samples[j] += (*shape)(timeFromPE) * npe;
125  }
126  }
127  }
128 }
virtual const CaloVShape * shape(const DetId &detId, bool precise=false) const
Definition: CaloShapes.h:15
Electronic response of the preamp.
Definition: CaloVShape.h:11
constexpr size_t nSamples
Definition: DetId.h:18
int size() const
get the size
Definition: CaloSamples.h:24
const CaloShapes * theShapes
DetId id() const
get the (generic) id
Definition: CaloSamples.h:21
double correctPE(const DetId &detId, double npe, CLHEP::HepRandomEngine *) const override
double HPDIonFeedbackSim::correctPE ( const DetId detId,
double  npe,
CLHEP::HepRandomEngine *  engine 
) const
overridevirtual

Implements CaloVPECorrection.

Definition at line 67 of file HPDIonFeedbackSim.cc.

References createfilelist::int, SiStripPI::max, p4, and pe2Charge.

Referenced by addThermalNoise(), and getIonFeedback().

68 {
69  double rateInTail = 0.000211988;//read this from XML file
70  double rateInSecondTail = 4.61579e-06;//read this from XML file
71 
72  // three gauss fit is applied to data to get ion feedback distribution
73  // parameters (in fC)
74  // first gaussian
75  // double p0 = 9.53192e+05;
76  // double p1 = -3.13653e-01;
77  // double p2 = 2.78350e+00;
78 
79  // second gaussian
80  // double p3 = 2.41611e+03;
81  double p4 = 2.06117e+01;
82  double p5 = 1.09239e+01;
83 
84  // third gaussian
85  // double p6 = 3.42793e+01;
86  double p7 = 5.45548e+01;
87  double p8 = 1.59696e+01;
88 
89  double noise = 0.; // fC
90  int nFirst = (int)(CLHEP::RandBinomial::shoot(engine, npe, rateInTail));
91  int nSecond = (int)(CLHEP::RandBinomial::shoot(engine, npe, rateInSecondTail));
92 
93  for (int j = 0; j < nFirst; ++j) {
94  noise += CLHEP::RandGaussQ::shoot(engine, p4, p5);
95  }
96  for (int j = 0; j < nSecond; ++j) {
97  noise += CLHEP::RandGaussQ::shoot(engine, p7, p8);
98  }
99 
100  return npe + std::max(noise/pe2Charge, 0.);
101 }
static const double pe2Charge
double p4[4]
Definition: TauolaWrapper.h:92
double HPDIonFeedbackSim::fCtoGeV ( const DetId detId) const
private

Definition at line 131 of file HPDIonFeedbackSim.cc.

References HcalDbService::getGain(), HcalDbService::getGainWidth(), HcalGain::getValue(), mps_fire::result, and theDbService.

Referenced by getIonFeedback().

132 {
133  assert(theDbService != nullptr);
134  HcalGenericDetId hcalGenDetId(detId);
135  const HcalGain* gains = theDbService->getGain(hcalGenDetId);
136  const HcalGainWidth* gwidths = theDbService->getGainWidth(hcalGenDetId);
137  double result = 0.0;
138  if (!gains || !gwidths )
139  {
140  edm::LogError("HcalAmplifier") << "Could not fetch HCAL conditions for channel " << hcalGenDetId;
141  }
142  else
143  {
144  // only one gain will be recorded per channel, so just use capID 0 for now
145  result = gains->getValue(0);
146  }
147  return result;
148 }
const HcalGainWidth * getGainWidth(const HcalGenericDetId &fId) const
float getValue(int fCapId) const
get value for capId = 0..3
Definition: HcalGain.h:22
const HcalDbService * theDbService
const HcalGain * getGain(const HcalGenericDetId &fId) const
double HPDIonFeedbackSim::getIonFeedback ( DetId  detId,
double  signal,
double  pedWidth,
bool  doThermal,
bool  isInGeV,
CLHEP::HepRandomEngine *  engine 
)

Definition at line 41 of file HPDIonFeedbackSim.cc.

References ALCARECOTkAlJpsiMuMu_cff::charge, correctPE(), fCtoGeV(), createfilelist::int, and pe2Charge.

42 {
43 
44  // HcalDetId id = detId;
45 
46  double GeVperfC = 1.;
47  if(isInGeV) GeVperfC = 1./fCtoGeV(detId);
48 
49  double charge = signal / GeVperfC;
50 
51  double noise = 0.; // fC
52  if (charge > 3. * pedWidth) { // 3 sigma away from pedestal mean
53  int npe = int (charge / pe2Charge);
54  if(doThermal) {
55  double electronEmission = 0.08;
56  CLHEP::RandPoissonQ theRandPoissonQ(*engine, electronEmission);
57  npe += theRandPoissonQ.fire();
58  }
59 
60  noise = correctPE(detId, npe, engine) - npe;
61  }
62  return (noise * GeVperfC);
63 
64 }
static const double pe2Charge
double fCtoGeV(const DetId &detId) const
double correctPE(const DetId &detId, double npe, CLHEP::HepRandomEngine *) const override
void HPDIonFeedbackSim::setDbService ( const HcalDbService service)
inline

Definition at line 35 of file HPDIonFeedbackSim.h.

References hcalTTPDigis_cfi::samples.

Referenced by HcalAmplifier::setDbService().

35 {theDbService = service;}
const HcalDbService * theDbService

Member Data Documentation

const HcalDbService* HPDIonFeedbackSim::theDbService
private

Definition at line 45 of file HPDIonFeedbackSim.h.

Referenced by fCtoGeV().

const CaloShapes* HPDIonFeedbackSim::theShapes
private

Definition at line 46 of file HPDIonFeedbackSim.h.

Referenced by addThermalNoise().