CMS 3D CMS Logo

CSCDbStripConditions.cc
Go to the documentation of this file.
10 
13  theConditions(pset),
14  theCapacitiveCrosstalk(pset.getParameter<double>("capacativeCrosstalk")),
15  theResistiveCrosstalkScaling(pset.getParameter<double>("resistiveCrosstalkScaling")),
16  theGainsConstant(pset.getParameter<double>("gainsConstant")),
17  doCorrelatedNoise_(pset.getParameter<bool>("doCorrelatedNoise")) {
18  // theCapacitiveCrosstalk = = 1/maxslope/maxsignal) = 1/ (0.00231/0.143);
19  // Howoever, need a bit more. Maybe the slope gets smeared?
20 }
21 
23  if (theNoisifier != nullptr)
24  delete theNoisifier;
25 }
26 
28 
29 float CSCDbStripConditions::gain(const CSCDetId &id, int channel) const {
30  return theConditions.gain(id, channel) * theGainsConstant;
31 }
32 
33 float CSCDbStripConditions::pedestal(const CSCDetId &id, int channel) const {
34  return theConditions.pedestal(id, channel);
35 }
36 
37 float CSCDbStripConditions::pedestalSigma(const CSCDetId &id, int channel) const {
38  return theConditions.pedestalSigma(id, channel);
39 }
40 
42  const CSCDetId &id, int channel, double stripLength, bool leftRight, float &capacitive, float &resistive) const {
43  resistive = theConditions.crosstalkIntercept(id, channel, leftRight) * theResistiveCrosstalkScaling;
44  float slope = theConditions.crosstalkSlope(id, channel, leftRight);
45  // ns before the peak where slope is max
46  float maxSlopeTime = 60.;
47  // some confusion about +/-
48  float capacitiveFraction = fabs(slope) * maxSlopeTime;
49  // theCapacitiveCrosstalk is the number needed for 100% xtalk, so
50  capacitive = theCapacitiveCrosstalk * capacitiveFraction;
51 }
52 
53 void CSCDbStripConditions::fetchNoisifier(const CSCDetId &id, int istrip) {
54  std::vector<float> me(12); // buffer for matrix elements
55  theConditions.noiseMatrixElements(id, istrip, me); // fill it
56 
58  // TODO get the pedestals right
59  matrix(2, 2) = me[0]; // item.elem33;
60  matrix(3, 3) = me[3]; // item.elem44;
61  matrix(4, 4) = me[6]; // item.elem55;
62  matrix(5, 5) = me[9]; // item.elem66;
63  matrix(6, 6) = me[11]; // item.elem77;
64 
65  if (doCorrelatedNoise_) {
66  matrix(2, 3) = me[1]; // item.elem34;
67  matrix(2, 4) = me[2]; // item.elem35;
68  matrix(3, 4) = me[4]; // item.elem45;
69  matrix(3, 5) = me[5]; // item.elem46;
70  matrix(4, 5) = me[7]; // item.elem56;
71  matrix(4, 6) = me[8]; // item.elem57;
72  matrix(5, 6) = me[10]; // item.elem67;
73  }
74 
75  // the other diagonal elements can just come from the pedestal sigma
76  float sigma = pedestalSigma(id, istrip);
77  //@@ float scaVariance = 2 * sigma * sigma;
78  //@@ The '2 *' IS strictly correct, but currently the value in the cond db is
79  // 2x too large since
80  //@@ it is the rms of the distribution of pedestals of all 8 time samples
81  // rather than the rms of
82  //@@ the average of the first two time samples
83  float scaVariance = sigma * sigma;
84  matrix(0, 0) = matrix(1, 1) = matrix(7, 7) = scaVariance;
85 
86  // unknown neighbors can be the average of the known neighbors
87  // float avgNeighbor = (matrix(2,3)+matrix(3,4)+matrix(4,5)+matrix(5,6))/4.;
88  // float avg2away = (matrix(2,4)+matrix(3,5)+matrix(4,6))/3.;
89  // matrix(0,1) = matrix(1,2) = matrix(6,7) = avgNeighbor;
90  // matrix(0,2) = matrix(1,3) = matrix(5,7) = avg2away;
91 
92  if (theNoisifier != nullptr)
93  delete theNoisifier;
95 }
96 
electrons_cff.bool
bool
Definition: electrons_cff.py:372
CSCDbStripConditions::initializeEvent
void initializeEvent(const edm::EventSetup &es) override
fetch the maps from the database
Definition: CSCDbStripConditions.cc:27
CSCDbStripConditions.h
makeMuonMisalignmentScenario.matrix
list matrix
Definition: makeMuonMisalignmentScenario.py:141
ESHandle.h
CSCDbStripConditions::isInBadChamber
bool isInBadChamber(const CSCDetId &id) const override
check list of bad chambers from db
Definition: CSCDbStripConditions.cc:97
CSCStripConditions::CSCCorrelatedNoiseMatrix
math::ErrorD< 8 >::type CSCCorrelatedNoiseMatrix
Definition: CSCStripConditions.h:15
CSCDbStripConditions::pedestal
float pedestal(const CSCDetId &detId, int channel) const override
in ADC counts
Definition: CSCDbStripConditions.cc:33
CSCConditions::pedestal
float pedestal(const CSCDetId &detId, int channel) const
static ped in ADC counts
Definition: CSCConditions.cc:237
CSCDbStripConditions::theResistiveCrosstalkScaling
float theResistiveCrosstalkScaling
Definition: CSCDbStripConditions.h:44
CSCDetId.h
CSCConditions::isInBadChamber
bool isInBadChamber(const CSCDetId &id) const
Is the gven chamber flagged as bad?
Definition: CSCConditions.cc:215
CSCDbStripConditions::doCorrelatedNoise_
bool doCorrelatedNoise_
Definition: CSCDbStripConditions.h:47
CSCConditions::gain
float gain(const CSCDetId &detId, int channel) const
gain per channel
Definition: CSCConditions.cc:229
CSCConditions::noiseMatrixElements
void noiseMatrixElements(const CSCDetId &id, int channel, std::vector< float > &me) const
Definition: CSCConditions.cc:285
CSCDbStripConditions::CSCDbStripConditions
CSCDbStripConditions(const edm::ParameterSet &pset)
Definition: CSCDbStripConditions.cc:11
CSCConditions::crosstalkIntercept
float crosstalkIntercept(const CSCDetId &detId, int channel, bool leftRight) const
crosstalk intercept for left and right
Definition: CSCConditions.cc:253
CSCDBNoiseMatrixRcd.h
CSCStripConditions
Definition: CSCStripConditions.h:13
CSCDbStripConditions::gain
float gain(const CSCDetId &detId, int channel) const override
channels count from 1
Definition: CSCDbStripConditions.cc:29
CSCDbStripConditions::fetchNoisifier
void fetchNoisifier(const CSCDetId &detId, int istrip) override
Definition: CSCDbStripConditions.cc:53
CSCStripConditions::CSCCorrelatedNoisifier
CorrelatedNoisifier< CSCCorrelatedNoiseMatrix > CSCCorrelatedNoisifier
Definition: CSCStripConditions.h:16
CSCStripConditions::theNoisifier
CSCCorrelatedNoisifier * theNoisifier
Definition: CSCStripConditions.h:52
edm::ParameterSet
Definition: ParameterSet.h:36
CSCDbStripConditions::~CSCDbStripConditions
~CSCDbStripConditions() override
Definition: CSCDbStripConditions.cc:22
CSCChannelTranslator.h
CSCDetId
Definition: CSCDetId.h:26
edm::EventSetup
Definition: EventSetup.h:57
CSCConditions::initializeEvent
void initializeEvent(const edm::EventSetup &es)
fetch database content via EventSetup
Definition: CSCConditions.cc:66
CSCDbStripConditions::pedestalSigma
float pedestalSigma(const CSCDetId &detId, int channel) const override
Definition: CSCDbStripConditions.cc:37
CSCConditions::pedestalSigma
float pedestalSigma(const CSCDetId &detId, int channel) const
static ped rms in ADC counts
Definition: CSCConditions.cc:245
CSCDBCrosstalkRcd.h
CSCConditions::crosstalkSlope
float crosstalkSlope(const CSCDetId &detId, int channel, bool leftRight) const
crosstalk slope for left and right
Definition: CSCConditions.cc:263
CSCDbStripConditions::theGainsConstant
float theGainsConstant
Definition: CSCDbStripConditions.h:46
CSCDBGainsRcd.h
CSCDbStripConditions::theConditions
CSCConditions theConditions
Definition: CSCDbStripConditions.h:38
Exception.h
hlt_dqm_clientPB-live_cfg.me
me
Definition: hlt_dqm_clientPB-live_cfg.py:56
slope
static const double slope[3]
Definition: CastorTimeSlew.cc:6
CSCDBPedestalsRcd.h
CSCDbStripConditions::crosstalk
void crosstalk(const CSCDetId &detId, int channel, double stripLength, bool leftRight, float &capacitive, float &resistive) const override
Definition: CSCDbStripConditions.cc:41
CSCDbStripConditions::theCapacitiveCrosstalk
float theCapacitiveCrosstalk
Definition: CSCDbStripConditions.h:41
muonDTDigis_cfi.pset
pset
Definition: muonDTDigis_cfi.py:27