#include <CSCDbStripConditions.h>
Public Member Functions | |
virtual void | crosstalk (const CSCDetId &detId, int channel, double stripLength, bool leftRight, float &capacitive, float &resistive) const |
CSCDbStripConditions (const edm::ParameterSet &pset) | |
virtual float | gain (const CSCDetId &detId, int channel) const |
channels count from 1 | |
virtual float | gainSigma (const CSCDetId &detId, int channel) const |
total calibration precision | |
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 | pedestal (const CSCDetId &detId, int channel) const |
in ADC counts | |
virtual float | pedestalSigma (const CSCDetId &detId, int channel) const |
virtual | ~CSCDbStripConditions () |
Private Member Functions | |
virtual void | fetchNoisifier (const CSCDetId &detId, int istrip) |
Private Attributes | |
bool | doCorrelatedNoise_ |
float | theCapacitiveCrosstalk |
CSCConditions | theConditions |
float | theGainsConstant |
float | theResistiveCrosstalkScaling |
Definition at line 8 of file CSCDbStripConditions.h.
CSCDbStripConditions::CSCDbStripConditions | ( | const edm::ParameterSet & | pset | ) | [explicit] |
Definition at line 12 of file CSCDbStripConditions.cc.
: CSCStripConditions(), theConditions( pset ), theCapacitiveCrosstalk(pset.getParameter<double>("capacativeCrosstalk")), theResistiveCrosstalkScaling(pset.getParameter<double>("resistiveCrosstalkScaling")), theGainsConstant(pset.getParameter<double>("gainsConstant")), doCorrelatedNoise_(pset.getParameter<bool>("doCorrelatedNoise")) { // theCapacitiveCrosstalk = = 1/maxslope/maxsignal) = 1/ (0.00231/0.143); // Howoever, need a bit more. Maybe the slope gets smeared? }
CSCDbStripConditions::~CSCDbStripConditions | ( | ) | [virtual] |
Definition at line 25 of file CSCDbStripConditions.cc.
References CSCStripConditions::theNoisifier.
{ if(theNoisifier != 0) delete theNoisifier; }
void CSCDbStripConditions::crosstalk | ( | const CSCDetId & | detId, |
int | channel, | ||
double | stripLength, | ||
bool | leftRight, | ||
float & | capacitive, | ||
float & | resistive | ||
) | const [virtual] |
Implements CSCStripConditions.
Definition at line 55 of file CSCDbStripConditions.cc.
References CSCConditions::crosstalkIntercept(), CSCConditions::crosstalkSlope(), slope, theCapacitiveCrosstalk, theConditions, and theResistiveCrosstalkScaling.
{ resistive = theConditions.crosstalkIntercept(id, channel, leftRight) * theResistiveCrosstalkScaling; float slope = theConditions.crosstalkSlope(id, channel, leftRight); // ns before the peak where slope is max float maxSlopeTime = 60.; // some confusion about +/- float capacitiveFraction = fabs(slope)*maxSlopeTime; // theCapacitiveCrosstalk is the number needed for 100% xtalk, so capacitive = theCapacitiveCrosstalk * capacitiveFraction; }
void CSCDbStripConditions::fetchNoisifier | ( | const CSCDetId & | detId, |
int | istrip | ||
) | [private, virtual] |
Implements CSCStripConditions.
Definition at line 71 of file CSCDbStripConditions.cc.
References doCorrelatedNoise_, makeMuonMisalignmentScenario::matrix, CSCConditions::noiseMatrixElements(), pedestalSigma(), theConditions, and CSCStripConditions::theNoisifier.
{ std::vector<float> me(12); // buffer for matrix elements theConditions.noiseMatrixElements( id, istrip, me ); // fill it CSCCorrelatedNoiseMatrix matrix; //TODO get the pedestals right matrix(2,2) = me[0]; // item.elem33; matrix(3,3) = me[3]; // item.elem44; matrix(4,4) = me[6]; // item.elem55; matrix(5,5) = me[9]; // item.elem66; matrix(6,6) = me[11]; // item.elem77; if(doCorrelatedNoise_) { matrix(2,3) = me[1]; // item.elem34; matrix(2,4) = me[2]; // item.elem35; matrix(3,4) = me[4]; // item.elem45; matrix(3,5) = me[5]; // item.elem46; matrix(4,5) = me[7]; // item.elem56; matrix(4,6) = me[8]; // item.elem57; matrix(5,6) = me[10]; // item.elem67; } // the other diagonal elements can just come from the pedestal sigma float sigma = pedestalSigma(id, istrip); //@@ float scaVariance = 2 * sigma * sigma; //@@ The '2 *' IS strictly correct, but currently the value in the cond db is 2x too large since //@@ it is the rms of the distribution of pedestals of all 8 time samples rather than the rms of //@@ the average of the first two time samples float scaVariance = sigma * sigma; matrix(0,0) = matrix(1,1) = matrix(7,7) = scaVariance; // unknown neighbors can be the average of the known neighbors //float avgNeighbor = (matrix(2,3)+matrix(3,4)+matrix(4,5)+matrix(5,6))/4.; //float avg2away = (matrix(2,4)+matrix(3,5)+matrix(4,6))/3.; //matrix(0,1) = matrix(1,2) = matrix(6,7) = avgNeighbor; //matrix(0,2) = matrix(1,3) = matrix(5,7) = avg2away; if(theNoisifier != 0) delete theNoisifier; theNoisifier = new CSCCorrelatedNoisifier(matrix); }
float CSCDbStripConditions::gain | ( | const CSCDetId & | detId, |
int | channel | ||
) | const [virtual] |
channels count from 1
Implements CSCStripConditions.
Definition at line 37 of file CSCDbStripConditions.cc.
References CSCConditions::gain(), theConditions, and theGainsConstant.
{ return theConditions.gain(id, channel) * theGainsConstant; }
virtual float CSCDbStripConditions::gainSigma | ( | const CSCDetId & | detId, |
int | channel | ||
) | const [inline, virtual] |
total calibration precision
Implements CSCStripConditions.
Definition at line 20 of file CSCDbStripConditions.h.
{return 0.005;}
void CSCDbStripConditions::initializeEvent | ( | const edm::EventSetup & | es | ) | [virtual] |
fetch the maps from the database
Reimplemented from CSCStripConditions.
Definition at line 31 of file CSCDbStripConditions.cc.
References CSCConditions::initializeEvent(), and theConditions.
{ theConditions.initializeEvent(es); }
bool CSCDbStripConditions::isInBadChamber | ( | const CSCDetId & | id | ) | const [virtual] |
check list of bad chambers from db
Reimplemented from CSCStripConditions.
Definition at line 114 of file CSCDbStripConditions.cc.
References CSCConditions::isInBadChamber(), and theConditions.
{ return theConditions.isInBadChamber( id ); }
float CSCDbStripConditions::pedestal | ( | const CSCDetId & | detId, |
int | channel | ||
) | const [virtual] |
in ADC counts
Implements CSCStripConditions.
Definition at line 43 of file CSCDbStripConditions.cc.
References CSCConditions::pedestal(), and theConditions.
{ return theConditions.pedestal(id, channel); }
float CSCDbStripConditions::pedestalSigma | ( | const CSCDetId & | detId, |
int | channel | ||
) | const [virtual] |
Implements CSCStripConditions.
Definition at line 49 of file CSCDbStripConditions.cc.
References CSCConditions::pedestalSigma(), and theConditions.
Referenced by fetchNoisifier().
{ return theConditions.pedestalSigma(id, channel); }
bool CSCDbStripConditions::doCorrelatedNoise_ [private] |
Definition at line 45 of file CSCDbStripConditions.h.
Referenced by fetchNoisifier().
float CSCDbStripConditions::theCapacitiveCrosstalk [private] |
Definition at line 39 of file CSCDbStripConditions.h.
Referenced by crosstalk().
Definition at line 36 of file CSCDbStripConditions.h.
Referenced by crosstalk(), fetchNoisifier(), gain(), initializeEvent(), isInBadChamber(), pedestal(), and pedestalSigma().
float CSCDbStripConditions::theGainsConstant [private] |
Definition at line 44 of file CSCDbStripConditions.h.
Referenced by gain().
float CSCDbStripConditions::theResistiveCrosstalkScaling [private] |
Definition at line 42 of file CSCDbStripConditions.h.
Referenced by crosstalk().