CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CSCDbStripConditions.cc
Go to the documentation of this file.
10 
11 
14  theConditions( pset ),
15  theCapacitiveCrosstalk(pset.getParameter<double>("capacativeCrosstalk")),
16  theResistiveCrosstalkScaling(pset.getParameter<double>("resistiveCrosstalkScaling")),
17  theGainsConstant(pset.getParameter<double>("gainsConstant")),
18  doCorrelatedNoise_(pset.getParameter<bool>("doCorrelatedNoise"))
19 {
20 // theCapacitiveCrosstalk = = 1/maxslope/maxsignal) = 1/ (0.00231/0.143);
21 // Howoever, need a bit more. Maybe the slope gets smeared?
22 }
23 
24 
26 {
27  if(theNoisifier != 0) delete theNoisifier;
28 }
29 
30 
32 {
34 }
35 
36 
37 float CSCDbStripConditions::gain(const CSCDetId & detId, int channel) const
38 {
39  CSCChannelTranslator translate;
40  CSCDetId idraw = translate.rawCSCDetId( detId );
41  int iraw = translate.rawStripChannel( detId, channel );
42  return theConditions.gain(idraw, iraw) * theGainsConstant;
43 }
44 
45 
46 
47 float CSCDbStripConditions::pedestal(const CSCDetId & detId, int channel) const
48 {
49  CSCChannelTranslator translate;
50  CSCDetId idraw = translate.rawCSCDetId( detId );
51  int iraw = translate.rawStripChannel( detId, channel );
52  return theConditions.pedestal(idraw, iraw);
53 }
54 
55 
56 float CSCDbStripConditions::pedestalSigma(const CSCDetId&detId, int channel) const
57 {
58  CSCChannelTranslator translate;
59  CSCDetId idraw = translate.rawCSCDetId( detId );
60  int iraw = translate.rawStripChannel( detId, channel );
61  return theConditions.pedestalSigma(idraw, iraw);
62 }
63 
64 
65 void CSCDbStripConditions::crosstalk(const CSCDetId&detId, int channel,
66  double stripLength, bool leftRight,
67  float & capacitive, float & resistive) const
68 {
69  CSCChannelTranslator translate;
70  CSCDetId idraw = translate.rawCSCDetId( detId );
71  int iraw = translate.rawStripChannel( detId, channel );
72 
73  resistive = theConditions.crosstalkIntercept(idraw, iraw, leftRight)
75  float slope = theConditions.crosstalkSlope(idraw, iraw, leftRight);
76  // ns before the peak where slope is max
77  float maxSlopeTime = 60.;
78  // some confusion about +/-
79  float capacitiveFraction = fabs(slope)*maxSlopeTime;
80  // theCapacitiveCrosstalk is the number needed for 100% xtalk, so
81  capacitive = theCapacitiveCrosstalk * capacitiveFraction;
82 }
83 
84 
85 
86 void CSCDbStripConditions::fetchNoisifier(const CSCDetId & detId, int istrip)
87 {
88  CSCChannelTranslator translate;
89  CSCDetId idraw = translate.rawCSCDetId( detId );
90  int iraw = translate.rawStripChannel( detId, istrip );
91 
92  std::vector<float> me(12); // buffer for matrix elements
93  theConditions.noiseMatrixElements( idraw, iraw, me ); // fill it
94 
96  //TODO get the pedestals right
97  matrix(2,2) = me[0]; // item.elem33;
98  matrix(3,3) = me[3]; // item.elem44;
99  matrix(4,4) = me[6]; // item.elem55;
100  matrix(5,5) = me[9]; // item.elem66;
101  matrix(6,6) = me[11]; // item.elem77;
102 
104  {
105  matrix(2,3) = me[1]; // item.elem34;
106  matrix(2,4) = me[2]; // item.elem35;
107  matrix(3,4) = me[4]; // item.elem45;
108  matrix(3,5) = me[5]; // item.elem46;
109  matrix(4,5) = me[7]; // item.elem56;
110  matrix(4,6) = me[8]; // item.elem57;
111  matrix(5,6) = me[10]; // item.elem67;
112  }
113 
114  // the other diagonal elements can just come from the pedestal sigma
115  float sigma = pedestalSigma(detId, istrip);
116  //@@ float scaVariance = 2 * sigma * sigma;
117  //@@ The '2 *' IS strictly correct, but currently the value in the cond db is 2x too large since
118  //@@ it is the rms of the distribution of pedestals of all 8 time samples rather than the rms of
119  //@@ the average of the first two time samples
120  float scaVariance = sigma * sigma;
121  matrix(0,0) = matrix(1,1) = matrix(7,7) = scaVariance;
122 
123  // unknown neighbors can be the average of the known neighbors
124  //float avgNeighbor = (matrix(2,3)+matrix(3,4)+matrix(4,5)+matrix(5,6))/4.;
125  //float avg2away = (matrix(2,4)+matrix(3,5)+matrix(4,6))/3.;
126  //matrix(0,1) = matrix(1,2) = matrix(6,7) = avgNeighbor;
127  //matrix(0,2) = matrix(1,3) = matrix(5,7) = avg2away;
128 
129  if(theNoisifier != 0) delete theNoisifier;
130  theNoisifier = new CSCCorrelatedNoisifier(matrix);
131 }
132 
134 {
135  return theConditions.isInBadChamber( detId );
136 }
float pedestalSigma(const CSCDetId &detId, int channel) const
void noiseMatrixElements(const CSCDetId &id, int channel, std::vector< float > &me) const
fill vector (dim 12, must be allocated by caller) with noise matrix elements (scaled to float) ...
static const double slope[3]
math::ErrorD< 8 >::type CSCCorrelatedNoiseMatrix
float pedestal(const CSCDetId &detId, int channel) const
in ADC counts
CSCDbStripConditions(const edm::ParameterSet &pset)
virtual float pedestal(const CSCDetId &detId, int channel) const
in ADC counts
tuple pset
Definition: CrabTask.py:85
void initializeEvent(const edm::EventSetup &es)
fetch the maps from the database
float crosstalkSlope(const CSCDetId &detId, int channel, bool leftRight) const
float crosstalkIntercept(const CSCDetId &detId, int channel, bool leftRight) const
bool isInBadChamber(const CSCDetId &id) const
Is the gven chamber flagged as bad?
float gain(const CSCDetId &detId, int channel) const
channels count from 1
virtual float pedestalSigma(const CSCDetId &detId, int channel) const
virtual void fetchNoisifier(const CSCDetId &detId, int istrip)
virtual void initializeEvent(const edm::EventSetup &es)
fetch the maps from the database
virtual bool isInBadChamber(const CSCDetId &id) const
check list of bad chambers from db
virtual float gain(const CSCDetId &detId, int channel) const
channels count from 1
CSCCorrelatedNoisifier * theNoisifier
int rawStripChannel(const CSCDetId &id, int igeom) const
Return raw strip channel number for input geometrical channel number.
CorrelatedNoisifier< CSCCorrelatedNoiseMatrix > CSCCorrelatedNoisifier
CSCDetId rawCSCDetId(const CSCDetId &id) const
virtual void crosstalk(const CSCDetId &detId, int channel, double stripLength, bool leftRight, float &capacitive, float &resistive) const