CMS 3D CMS Logo

CSCConditions.cc

Go to the documentation of this file.
00001 #include "CondFormats/CSCObjects/interface/CSCConditions.h"
00002 #include "FWCore/Framework/interface/ESHandle.h"
00003 #include "CondFormats/DataRecord/interface/CSCDBPedestalsRcd.h"
00004 #include "CondFormats/DataRecord/interface/CSCDBNoiseMatrixRcd.h"
00005 #include "CondFormats/DataRecord/interface/CSCDBCrosstalkRcd.h"
00006 #include "CondFormats/CSCObjects/interface/CSCDBGains.h"
00007 #include "CondFormats/CSCObjects/interface/CSCDBPedestals.h"
00008 #include "CondFormats/CSCObjects/interface/CSCDBCrosstalk.h"
00009 #include "CondFormats/CSCObjects/interface/CSCBadStrips.h"
00010 #include "CondFormats/CSCObjects/interface/CSCBadWires.h"
00011 #include "CondFormats/CSCObjects/interface/CSCBadChambers.h"
00012 #include "DataFormats/MuonDetId/interface/CSCIndexer.h"
00013 
00014 CSCConditions::CSCConditions( const edm::ParameterSet& ps ) 
00015  :theNoiseMatrix(0),
00016   theGains(0),
00017   thePedestals(0),
00018   theCrosstalk(0),
00019   theBadStrips(0),
00020   theBadWires(0),
00021   theBadChambers(0),
00022   readBadChannels_(false), readBadChambers_(false),
00023   theAverageGain( -1.0 )
00024 {
00025   readBadChannels_ = ps.getParameter<bool>("readBadChannels");
00026   readBadChambers_ = ps.getParameter<bool>("readBadChambers");
00027   // initialize #layers = 2808
00028   badStripWords.resize( 2808, 0 );
00029   badWireWords.resize( 2808, 0 );
00030 }
00031 
00032 
00033 CSCConditions::~CSCConditions()
00034 {
00035 }
00036 
00037 
00038 void CSCConditions::initializeEvent(const edm::EventSetup & es)
00039 {
00040   // Strip gains
00041   edm::ESHandle<CSCDBGains> hGains;
00042   es.get<CSCDBGainsRcd>().get( hGains );
00043   theGains = hGains.product();
00044   // Strip X-talk
00045   edm::ESHandle<CSCDBCrosstalk> hCrosstalk;
00046   es.get<CSCDBCrosstalkRcd>().get( hCrosstalk );
00047   theCrosstalk = hCrosstalk.product();
00048   // Strip pedestals
00049   edm::ESHandle<CSCDBPedestals> hPedestals;
00050   es.get<CSCDBPedestalsRcd>().get( hPedestals );
00051   thePedestals = hPedestals.product();
00052   // Strip autocorrelation noise matrix
00053   edm::ESHandle<CSCDBNoiseMatrix> hNoiseMatrix;
00054   es.get<CSCDBNoiseMatrixRcd>().get(hNoiseMatrix);
00055   theNoiseMatrix = hNoiseMatrix.product();
00056 
00057   if ( readBadChannels() ) {
00058   // Bad strip channels
00059     edm::ESHandle<CSCBadStrips> hBadS;
00060     es.get<CSCBadStripsRcd>().get( hBadS );
00061     theBadStrips = hBadS.product();
00062   // Bad wiregroup channels
00063     edm::ESHandle<CSCBadWires> hBadW;
00064     es.get<CSCBadWiresRcd>().get( hBadW );
00065     theBadWires = hBadW.product();
00066 
00067     //@@    if( badStripsWatcher_.check( es ) ) { 
00068       fillBadStripWords();
00069     //@@    }
00070     //@@    if( badWiresWatcher_.check( es ) ) { 
00071       fillBadWireWords();
00072     //@    }
00073 
00074   }
00075 
00076   // Has GainsRcd changed?
00077   if( gainsWatcher_.check( es ) ) { // Yes...
00078     theAverageGain = -1.0; // ...reset, so next access will recalculate it
00079   }
00080  
00081   if ( readBadChambers() ) {
00082   // Entire bad chambers
00083     edm::ESHandle<CSCBadChambers> hBadC;
00084     es.get<CSCBadChambersRcd>().get( hBadC );
00085     theBadChambers = hBadC.product();
00086   }
00087 
00088 //  print();
00089 }
00090 
00091 void CSCConditions::fillBadStripWords(){
00092   // reset existing values
00093   badStripWords.assign( 2808, 0 );
00094   if ( readBadChannels() ) {
00095     // unpack what we've read from theBadStrips
00096 
00097     // chambers is a vector<BadChamber>
00098     // channels is a vector<BadChannel>
00099     // Each BadChamber contains its index (1-468), the no. of bad channels, 
00100     // and the index within vector<BadChannel> where this chamber's bad channels start.
00101 
00102     CSCIndexer indexer;
00103 
00104     for ( size_t i=0; i<theBadStrips->chambers.size(); ++i ) { // loop over bad chambers
00105       int indexc = theBadStrips->chambers[i].chamber_index;
00106       int start =  theBadStrips->chambers[i].pointer;  // where this chamber's bad channels start in vector<BadChannel>
00107       int nbad  =  theBadStrips->chambers[i].bad_channels;
00108 
00109       CSCDetId id = indexer.detIdFromChamberIndex( indexc ); // We need this to build layer index (1-2808)
00110 
00111       for ( int j=start-1; j<start-1+nbad; ++j ) { // bad channels in this chamber
00112         short lay  = theBadStrips->channels[j].layer;    // value 1-6
00113         short chan = theBadStrips->channels[j].channel;  // value 1-80
00114     //    short f1 = theBadStrips->channels[j].flag1;
00115     //    short f2 = theBadStrips->channels[j].flag2;
00116     //    short f3 = theBadStrips->channels[j].flag3;
00117         int indexl = indexer.layerIndex( id.endcap(), id.station(), id.ring(), id.chamber(), lay );
00118         badStripWords[indexl-1].set( chan-1, 1 ); // set bit 0-79 in 80-bit bitset representing this layer
00119       } // j
00120     } // i
00121 
00122   } 
00123 }
00124 
00125 void CSCConditions::fillBadWireWords(){
00126   // reset existing values
00127   badWireWords.assign( 2808, 0 );
00128   if ( readBadChannels() ) {
00129     // unpack what we've read from theBadWires
00130     CSCIndexer indexer;
00131 
00132     for ( size_t i=0; i<theBadWires->chambers.size(); ++i ) { // loop over bad chambers
00133       int indexc = theBadWires->chambers[i].chamber_index;
00134       int start =  theBadWires->chambers[i].pointer;  // where this chamber's bad channels start in vector<BadChannel>
00135       int nbad  =  theBadWires->chambers[i].bad_channels;
00136 
00137       CSCDetId id = indexer.detIdFromChamberIndex( indexc ); // We need this to build layer index (1-2808)
00138 
00139       for ( int j=start-1; j<start-1+nbad; ++j ) { // bad channels in this chamber
00140         short lay  = theBadWires->channels[j].layer;    // value 1-6
00141         short chan = theBadWires->channels[j].channel;  // value 1-80
00142     //    short f1 = theBadWires->channels[j].flag1;
00143     //    short f2 = theBadWires->channels[j].flag2;
00144     //    short f3 = theBadWires->channels[j].flag3;
00145         int indexl = indexer.layerIndex( id.endcap(), id.station(), id.ring(), id.chamber(), lay );
00146         badWireWords[indexl-1].set( chan-1, 1 ); // set bit 0-111 in 112-bit bitset representing this layer
00147       } // j
00148     } // i
00149 
00150   } 
00151 }
00152 
00153 bool CSCConditions::isInBadChamber( const CSCDetId& id ) const {
00154 
00155   if ( theBadChambers->numberOfBadChambers == 0 ) return false;
00156 
00157   short int iri = id.ring();
00158   if ( iri == 4 ) iri = 1; // reset ME1A to ME11
00159   CSCIndexer indexer;
00160   int ilin = indexer.chamberIndex( id.endcap(), id.station(), iri, id.chamber() );
00161   std::vector<int>::const_iterator badbegin = theBadChambers->chambers.begin();
00162   std::vector<int>::const_iterator badend = theBadChambers->chambers.end();
00163   std::vector<int>::const_iterator it = std::find( badbegin, badend, ilin );
00164   if ( it != badend ) return true; // id is in the list of bad chambers
00165   else return false;
00166 }
00167 
00168 void CSCConditions::print() const
00169   //@@ NEEDS THOROUGH UPDATING
00170 {
00171 /*
00172   std::cout << "SIZES: GAINS: " << theGains->gains.size()
00173             << "   PEDESTALS: " << thePedestals->pedestals.size()
00174             << "   NOISES "  << theNoiseMatrix->matrix.size() << std::endl;;
00175 
00176   std::map< int,std::vector<CSCDBGains::Item> >::const_iterator layerGainsItr = theGains->gains.begin(), 
00177       lastGain = theGains->gains.end();
00178   for( ; layerGainsItr != lastGain; ++layerGainsItr)
00179   {
00180     std::cout << "GAIN " << layerGainsItr->first 
00181               << " STRIPS " << layerGainsItr->second.size() << " "
00182               << layerGainsItr->second[0].gain_slope 
00183               << " " << layerGainsItr->second[0].gain_intercept << std::endl;
00184   }
00185 
00186   std::map< int,std::vector<CSCDBPedestals::Item> >::const_iterator pedestalItr = thePedestals->pedestals.begin(), 
00187                                                                   lastPedestal = thePedestals->pedestals.end();
00188   for( ; pedestalItr != lastPedestal; ++pedestalItr)
00189   {
00190     std::cout << "PEDS " << pedestalItr->first << " " 
00191               << " STRIPS " << pedestalItr->second.size() << " ";
00192     for(int i = 1; i < 80; ++i)
00193     {
00194        std::cout << pedestalItr->second[i-1].rms << " " ;
00195      }
00196      std::cout << std::endl;
00197   }
00198 
00199   std::map< int,std::vector<CSCDBCrosstalk::Item> >::const_iterator crosstalkItr = theCrosstalk->crosstalk.begin(),
00200                                                                   lastCrosstalk = theCrosstalk->crosstalk.end();
00201   for( ; crosstalkItr != lastCrosstalk; ++crosstalkItr)
00202   {
00203     std::cout << "XTALKS " << crosstalkItr->first 
00204       << " STRIPS " << crosstalkItr->second.size() << " "  
00205      << crosstalkItr->second[5].xtalk_slope_left << " " 
00206      << crosstalkItr->second[5].xtalk_slope_right << " " 
00207      << crosstalkItr->second[5].xtalk_intercept_left << " " 
00208      << crosstalkItr->second[5].xtalk_intercept_right << std::endl;
00209   }
00210 */
00211 }
00212 
00213 float CSCConditions::gain(const CSCDetId & detId, int channel) const
00214 {
00215   assert(theGains != 0);
00216   return float( theGains->item(detId, channel).gain_slope )/theGains->factor_gain;
00217 }
00218 
00219 
00220 float CSCConditions::pedestal(const CSCDetId & detId, int channel) const
00221 {
00222   assert(thePedestals != 0);
00223   return float ( thePedestals->item(detId, channel).ped )/thePedestals->factor_ped;
00224 }
00225 
00226 
00227 float CSCConditions::pedestalSigma(const CSCDetId&detId, int channel) const
00228 {
00229   assert(thePedestals != 0);
00230   return float ( thePedestals->item(detId, channel).rms )/thePedestals->factor_rms;
00231 }
00232 
00233 
00234 float CSCConditions::crosstalkIntercept(const CSCDetId&detId, int channel, bool leftRight) const
00235 {
00236   assert(theCrosstalk != 0);
00237   const CSCDBCrosstalk::Item & item = theCrosstalk->item(detId, channel);
00238 
00239   // resistive fraction is at the peak, where t=0
00240   return leftRight ? float ( item.xtalk_intercept_right )/theCrosstalk->factor_intercept 
00241                    : float ( item.xtalk_intercept_left )/theCrosstalk->factor_intercept ;
00242 }
00243 
00244 
00245 float CSCConditions::crosstalkSlope(const CSCDetId&detId, int channel, bool leftRight) const
00246 {
00247   assert(theCrosstalk != 0);
00248   const CSCDBCrosstalk::Item & item = theCrosstalk->item(detId, channel);
00249 
00250   // resistive fraction is at the peak, where t=0
00251   return leftRight ? float ( item.xtalk_slope_right )/theCrosstalk->factor_slope
00252                    : float ( item.xtalk_slope_left )/theCrosstalk->factor_slope ;
00253 }
00254 
00255 const CSCDBNoiseMatrix::Item & CSCConditions::noiseMatrix(const CSCDetId& detId, int channel) const
00256 {
00257   assert(theNoiseMatrix != 0);
00258   return theNoiseMatrix->item(detId, channel);
00259 }
00260 
00261 void CSCConditions::noiseMatrixElements( const CSCDetId& id, int channel, std::vector<float>& me ) const {
00262   assert(me.size()>11);
00263   const CSCDBNoiseMatrix::Item& item = noiseMatrix(id, channel);
00264   me[0] = float ( item.elem33 )/theNoiseMatrix->factor_noise;
00265   me[1] = float ( item.elem34 )/theNoiseMatrix->factor_noise;
00266   me[2] = float ( item.elem35 )/theNoiseMatrix->factor_noise;
00267   me[3] = float ( item.elem44 )/theNoiseMatrix->factor_noise;
00268   me[4] = float ( item.elem45 )/theNoiseMatrix->factor_noise;
00269   me[5] = float ( item.elem46 )/theNoiseMatrix->factor_noise;
00270   me[6] = float ( item.elem55 )/theNoiseMatrix->factor_noise;
00271   me[7] = float ( item.elem56 )/theNoiseMatrix->factor_noise;
00272   me[8] = float ( item.elem57 )/theNoiseMatrix->factor_noise;
00273   me[9] = float ( item.elem66 )/theNoiseMatrix->factor_noise;
00274   me[10] = float ( item.elem67 )/theNoiseMatrix->factor_noise;
00275   me[11] = float ( item.elem77 )/theNoiseMatrix->factor_noise;
00276 }
00277 
00278 void CSCConditions::crossTalk( const CSCDetId& id, int channel, std::vector<float>& ct ) const {
00279   assert(theCrosstalk != 0);
00280   const CSCDBCrosstalk::Item & item = theCrosstalk->item(id, channel);
00281   ct[0] = float ( item.xtalk_slope_left )/theCrosstalk->factor_slope;
00282   ct[1] = float ( item.xtalk_intercept_left )/theCrosstalk->factor_intercept;
00283   ct[2] = float ( item.xtalk_slope_right )/theCrosstalk->factor_slope;
00284   ct[3] = float ( item.xtalk_intercept_right )/theCrosstalk->factor_intercept;
00285 }
00286 
00287 const std::bitset<80>& CSCConditions::badStripWord( const CSCDetId& id ) const {
00288   CSCIndexer indexer;
00289   return badStripWords[indexer.layerIndex(id) - 1];
00290 }
00291 
00292 const std::bitset<112>& CSCConditions::badWireWord( const CSCDetId& id ) const {
00293   CSCIndexer indexer;
00294   return badWireWords[indexer.layerIndex(id) - 1];
00295 }
00296 
00300 float CSCConditions::averageGain() const {
00301 
00302   const float loEdge = 5.0; // consider gains above this
00303   const float hiEdge = 10.0; // consider gains below this
00304   const float loLimit = 6.0; // lowest acceptable average gain
00305   const float hiLimit = 9.0; // highest acceptable average gain
00306   const float expectedAverage = 7.5; // default average gain
00307 
00308   if ( theAverageGain > 0. ) return theAverageGain; // only recalculate if necessary
00309 
00310   int  n_strip   = 0;
00311   float gain_tot = 0.;
00312 
00313   CSCDBGains::GainContainer::const_iterator it;
00314   for ( it=theGains->gains.begin(); it!=theGains->gains.end(); ++it ) {
00315     float the_gain = float( it->gain_slope )/theGains->factor_gain;
00316     if (the_gain > loEdge && the_gain < hiEdge ) {
00317       gain_tot += the_gain;
00318       ++n_strip;
00319     }
00320   }
00321 
00322   // Average gain
00323   if ( n_strip > 0 ) {
00324     theAverageGain = gain_tot / n_strip;
00325   }
00326 
00327   // Average gain has been around 7.5 in real data
00328   if ( theAverageGain < loLimit || theAverageGain > hiLimit ) {
00329     //    LogTrace("CSC") << "Average CSC strip gain = "
00330     //                    << theAverageGain << "  is reset to expected value " << expectedAverage;
00331     theAverageGain = expectedAverage;
00332   }
00333 
00334   return theAverageGain;
00335 }

Generated on Tue Jun 9 17:26:25 2009 for CMSSW by  doxygen 1.5.4