test
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 & id, int channel) const
38 {
39  return theConditions.gain(id, channel) * theGainsConstant;
40 }
41 
42 
43 float CSCDbStripConditions::pedestal(const CSCDetId & id, int channel) const
44 {
45  return theConditions.pedestal(id, channel);
46 }
47 
48 
49 float CSCDbStripConditions::pedestalSigma(const CSCDetId& id, int channel) const
50 {
51  return theConditions.pedestalSigma(id, channel);
52 }
53 
54 
55 void CSCDbStripConditions::crosstalk(const CSCDetId& id, int channel,
56  double stripLength, bool leftRight,
57  float & capacitive, float & resistive) const
58 {
59  resistive = theConditions.crosstalkIntercept(id, channel, leftRight)
61  float slope = theConditions.crosstalkSlope(id, channel, leftRight);
62  // ns before the peak where slope is max
63  float maxSlopeTime = 60.;
64  // some confusion about +/-
65  float capacitiveFraction = fabs(slope)*maxSlopeTime;
66  // theCapacitiveCrosstalk is the number needed for 100% xtalk, so
67  capacitive = theCapacitiveCrosstalk * capacitiveFraction;
68 }
69 
70 
71 void CSCDbStripConditions::fetchNoisifier(const CSCDetId & id, int istrip)
72 {
73  std::vector<float> me(12); // buffer for matrix elements
74  theConditions.noiseMatrixElements( id, istrip, me ); // fill it
75 
77  //TODO get the pedestals right
78  matrix(2,2) = me[0]; // item.elem33;
79  matrix(3,3) = me[3]; // item.elem44;
80  matrix(4,4) = me[6]; // item.elem55;
81  matrix(5,5) = me[9]; // item.elem66;
82  matrix(6,6) = me[11]; // item.elem77;
83 
85  {
86  matrix(2,3) = me[1]; // item.elem34;
87  matrix(2,4) = me[2]; // item.elem35;
88  matrix(3,4) = me[4]; // item.elem45;
89  matrix(3,5) = me[5]; // item.elem46;
90  matrix(4,5) = me[7]; // item.elem56;
91  matrix(4,6) = me[8]; // item.elem57;
92  matrix(5,6) = me[10]; // item.elem67;
93  }
94 
95  // the other diagonal elements can just come from the pedestal sigma
96  float sigma = pedestalSigma(id, istrip);
97  //@@ float scaVariance = 2 * sigma * sigma;
98  //@@ The '2 *' IS strictly correct, but currently the value in the cond db is 2x too large since
99  //@@ it is the rms of the distribution of pedestals of all 8 time samples rather than the rms of
100  //@@ the average of the first two time samples
101  float scaVariance = sigma * sigma;
102  matrix(0,0) = matrix(1,1) = matrix(7,7) = scaVariance;
103 
104  // unknown neighbors can be the average of the known neighbors
105  //float avgNeighbor = (matrix(2,3)+matrix(3,4)+matrix(4,5)+matrix(5,6))/4.;
106  //float avg2away = (matrix(2,4)+matrix(3,5)+matrix(4,6))/3.;
107  //matrix(0,1) = matrix(1,2) = matrix(6,7) = avgNeighbor;
108  //matrix(0,2) = matrix(1,3) = matrix(5,7) = avg2away;
109 
110  if(theNoisifier != 0) delete theNoisifier;
111  theNoisifier = new CSCCorrelatedNoisifier(matrix);
112 }
113 
115 {
116  return theConditions.isInBadChamber( id );
117 }
float pedestalSigma(const CSCDetId &detId, int channel) const
static ped rms in ADC counts
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
static ped in ADC counts
CSCDbStripConditions(const edm::ParameterSet &pset)
virtual float pedestal(const CSCDetId &detId, int channel) const
in ADC counts
void initializeEvent(const edm::EventSetup &es)
fetch database content via EventSetup
float crosstalkSlope(const CSCDetId &detId, int channel, bool leftRight) const
crosstalk slope for left and right
float crosstalkIntercept(const CSCDetId &detId, int channel, bool leftRight) const
crosstalk intercept for left and right
bool isInBadChamber(const CSCDetId &id) const
Is the gven chamber flagged as bad?
float gain(const CSCDetId &detId, int channel) const
gain per channel
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
CorrelatedNoisifier< CSCCorrelatedNoiseMatrix > CSCCorrelatedNoisifier
virtual void crosstalk(const CSCDetId &detId, int channel, double stripLength, bool leftRight, float &capacitive, float &resistive) const