CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/SimMuon/CSCDigitizer/src/CSCDbStripConditions.cc

Go to the documentation of this file.
00001 #include "SimMuon/CSCDigitizer/src/CSCDbStripConditions.h"
00002 #include "DataFormats/MuonDetId/interface/CSCDetId.h"
00003 #include "FWCore/Utilities/interface/Exception.h"
00004 #include "FWCore/Framework/interface/ESHandle.h"
00005 #include "CondFormats/CSCObjects/interface/CSCChannelTranslator.h"
00006 #include "CondFormats/DataRecord/interface/CSCDBGainsRcd.h"
00007 #include "CondFormats/DataRecord/interface/CSCDBPedestalsRcd.h"
00008 #include "CondFormats/DataRecord/interface/CSCDBNoiseMatrixRcd.h"
00009 #include "CondFormats/DataRecord/interface/CSCDBCrosstalkRcd.h"
00010 
00011 
00012 CSCDbStripConditions::CSCDbStripConditions(const edm::ParameterSet & pset) 
00013 : CSCStripConditions(),
00014   theConditions( pset ),
00015   theCapacitiveCrosstalk(pset.getParameter<double>("capacativeCrosstalk")),
00016   theResistiveCrosstalkScaling(pset.getParameter<double>("resistiveCrosstalkScaling")),
00017   theGainsConstant(pset.getParameter<double>("gainsConstant")),
00018   doCorrelatedNoise_(pset.getParameter<bool>("doCorrelatedNoise"))
00019 {
00020 //  theCapacitiveCrosstalk = = 1/maxslope/maxsignal) = 1/ (0.00231/0.143);
00021 // Howoever, need a bit more.  Maybe the slope gets smeared?
00022 }
00023 
00024 
00025 CSCDbStripConditions::~CSCDbStripConditions()
00026 {
00027   if(theNoisifier != 0) delete theNoisifier;
00028 }
00029 
00030 
00031 void CSCDbStripConditions::initializeEvent(const edm::EventSetup & es)
00032 {
00033   theConditions.initializeEvent(es);
00034 }
00035 
00036 
00037 float CSCDbStripConditions::gain(const CSCDetId & detId, int channel) const
00038 {
00039   CSCChannelTranslator translate;
00040   CSCDetId idraw = translate.rawCSCDetId( detId );
00041   int iraw = translate.rawStripChannel( detId, channel );
00042   return theConditions.gain(idraw, iraw)  * theGainsConstant;
00043 }
00044 
00045 
00046 
00047 float CSCDbStripConditions::pedestal(const CSCDetId & detId, int channel) const
00048 {
00049   CSCChannelTranslator translate;
00050   CSCDetId idraw = translate.rawCSCDetId( detId );
00051   int iraw = translate.rawStripChannel( detId, channel );
00052   return theConditions.pedestal(idraw, iraw);
00053 }
00054 
00055 
00056 float CSCDbStripConditions::pedestalSigma(const CSCDetId&detId, int channel) const
00057 {
00058   CSCChannelTranslator translate;
00059   CSCDetId idraw = translate.rawCSCDetId( detId );
00060   int iraw = translate.rawStripChannel( detId, channel );
00061   return theConditions.pedestalSigma(idraw, iraw);
00062 }
00063 
00064 
00065 void CSCDbStripConditions::crosstalk(const CSCDetId&detId, int channel, 
00066                         double stripLength, bool leftRight, 
00067                         float & capacitive, float & resistive) const
00068 {
00069   CSCChannelTranslator translate;
00070   CSCDetId idraw = translate.rawCSCDetId( detId );
00071   int iraw = translate.rawStripChannel( detId, channel );
00072         
00073   resistive = theConditions.crosstalkIntercept(idraw, iraw, leftRight)
00074              * theResistiveCrosstalkScaling;
00075   float slope = theConditions.crosstalkSlope(idraw, iraw, leftRight);
00076   // ns before the peak where slope is max
00077   float maxSlopeTime = 60.; 
00078   // some confusion about +/-
00079   float capacitiveFraction = fabs(slope)*maxSlopeTime;
00080   // theCapacitiveCrosstalk is the number needed for 100% xtalk, so
00081   capacitive = theCapacitiveCrosstalk * capacitiveFraction;
00082 }
00083 
00084 
00085 
00086 void CSCDbStripConditions::fetchNoisifier(const CSCDetId & detId, int istrip)
00087 {
00088   CSCChannelTranslator translate;
00089   CSCDetId idraw = translate.rawCSCDetId( detId );
00090   int iraw = translate.rawStripChannel( detId, istrip );
00091         
00092   std::vector<float> me(12); // buffer for matrix elements
00093   theConditions.noiseMatrixElements( idraw, iraw, me ); // fill it
00094 
00095   CSCCorrelatedNoiseMatrix matrix;
00096   //TODO get the pedestals right
00097   matrix(2,2) = me[0]; // item.elem33;
00098   matrix(3,3) = me[3]; // item.elem44;
00099   matrix(4,4) = me[6]; // item.elem55;
00100   matrix(5,5) = me[9]; // item.elem66;
00101   matrix(6,6) = me[11]; // item.elem77;
00102   
00103   if(doCorrelatedNoise_)
00104   {
00105     matrix(2,3) = me[1]; // item.elem34;
00106     matrix(2,4) = me[2]; // item.elem35;
00107     matrix(3,4) = me[4]; // item.elem45;
00108     matrix(3,5) = me[5]; // item.elem46;
00109     matrix(4,5) = me[7]; // item.elem56;
00110     matrix(4,6) = me[8]; // item.elem57;
00111     matrix(5,6) = me[10]; // item.elem67;
00112   }
00113 
00114   // the other diagonal elements can just come from the pedestal sigma
00115   float sigma = pedestalSigma(detId, istrip);
00116   //@@  float scaVariance = 2 * sigma * sigma;
00117   //@@ The '2 *' IS strictly correct, but currently the value in the cond db is 2x too large since
00118   //@@ it is the rms of the distribution of pedestals of all 8 time samples rather than the rms of
00119   //@@ the average of the first two time samples
00120   float scaVariance = sigma * sigma;
00121   matrix(0,0) = matrix(1,1) = matrix(7,7) = scaVariance;
00122 
00123   // unknown neighbors can be the average of the known neighbors
00124   //float avgNeighbor = (matrix(2,3)+matrix(3,4)+matrix(4,5)+matrix(5,6))/4.;
00125   //float avg2away = (matrix(2,4)+matrix(3,5)+matrix(4,6))/3.;
00126   //matrix(0,1) = matrix(1,2) = matrix(6,7) = avgNeighbor;
00127   //matrix(0,2) = matrix(1,3) = matrix(5,7) = avg2away;
00128 
00129   if(theNoisifier != 0) delete theNoisifier;
00130   theNoisifier = new CSCCorrelatedNoisifier(matrix);
00131 }
00132 
00133 bool CSCDbStripConditions::isInBadChamber( const CSCDetId& detId ) const
00134 {
00135   return theConditions.isInBadChamber( detId );
00136 }