CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch2/src/DQM/SiStripCommissioningAnalysis/src/PedsOnlyAlgorithm.cc

Go to the documentation of this file.
00001 #include "DQM/SiStripCommissioningAnalysis/interface/PedsOnlyAlgorithm.h"
00002 #include "CondFormats/SiStripObjects/interface/PedsOnlyAnalysis.h" 
00003 #include "DataFormats/SiStripCommon/interface/SiStripHistoTitle.h"
00004 #include "DataFormats/SiStripCommon/interface/SiStripEnumsAndStrings.h"
00005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00006 #include "TProfile.h"
00007 #include "TH1.h"
00008 #include <iostream>
00009 #include <iomanip>
00010 #include <cmath>
00011 
00012 using namespace sistrip;
00013 
00014 // ----------------------------------------------------------------------------
00015 // 
00016 PedsOnlyAlgorithm::PedsOnlyAlgorithm( const edm::ParameterSet & pset, PedsOnlyAnalysis* const anal ) 
00017   : CommissioningAlgorithm(anal),
00018     hPeds_(0,""),
00019     hNoise_(0,"")
00020 {}
00021 
00022 // ----------------------------------------------------------------------------
00023 // 
00024 void PedsOnlyAlgorithm::extract( const std::vector<TH1*>& histos ) { 
00025 
00026   if ( !anal() ) {
00027     edm::LogWarning(mlCommissioning_)
00028       << "[PedsOnlyAlgorithm::" << __func__ << "]"
00029       << " NULL pointer to Analysis object!";
00030     return; 
00031   }
00032 
00033   // Check number of histograms
00034   if ( histos.size() != 2 ) {
00035     anal()->addErrorCode(sistrip::numberOfHistos_);
00036   }
00037   
00038   // Extract FED key from histo title
00039   if ( !histos.empty() ) { anal()->fedKey( extractFedKey( histos.front() ) ); }
00040   
00041   // Extract histograms
00042   std::vector<TH1*>::const_iterator ihis = histos.begin();
00043   for ( ; ihis != histos.end(); ihis++ ) {
00044     
00045     // Check for NULL pointer
00046     if ( !(*ihis) ) { continue; }
00047 
00048     // Check run type
00049     SiStripHistoTitle title( (*ihis)->GetName() );
00050     if ( title.runType() != sistrip::PEDS_ONLY ) {
00051       anal()->addErrorCode(sistrip::unexpectedTask_);
00052       continue;
00053     }
00054     
00055     // Extract peds and raw noise histos (check for legacy names first!)
00056     if ( title.extraInfo().find(sistrip::extrainfo::pedsAndRawNoise_) != std::string::npos ) {
00057       hPeds_.first = *ihis;
00058       hPeds_.second = (*ihis)->GetName();
00059       hNoise_.first = *ihis;
00060       hNoise_.second = (*ihis)->GetName();
00061       PedsOnlyAnalysis* a = dynamic_cast<PedsOnlyAnalysis*>( const_cast<CommissioningAnalysis*>( anal() ) );
00062       if ( a ) { a->legacy_ = true; }
00063     } else if ( title.extraInfo().find(sistrip::extrainfo::pedestals_) != std::string::npos ) {
00064       hPeds_.first = *ihis;
00065       hPeds_.second = (*ihis)->GetName();
00066     } else if ( title.extraInfo().find(sistrip::extrainfo::rawNoise_) != std::string::npos ) {
00067       hNoise_.first = *ihis;
00068       hNoise_.second = (*ihis)->GetName();
00069     } else { 
00070       anal()->addErrorCode(sistrip::unexpectedExtraInfo_);
00071     }
00072     
00073   }
00074 
00075 }
00076 
00077 // -----------------------------------------------------------------------------
00078 // 
00079 void PedsOnlyAlgorithm::analyse() {
00080   
00081   if ( !anal() ) {
00082     edm::LogWarning(mlCommissioning_)
00083       << "[PedsOnlyAlgorithm::" << __func__ << "]"
00084       << " NULL pointer to base Analysis object!";
00085     return; 
00086   }
00087   
00088   CommissioningAnalysis* tmp = const_cast<CommissioningAnalysis*>( anal() );
00089   PedsOnlyAnalysis* anal = dynamic_cast<PedsOnlyAnalysis*>( tmp );
00090   if ( !anal ) {
00091     edm::LogWarning(mlCommissioning_)
00092       << "[PedsOnlyAlgorithm::" << __func__ << "]"
00093       << " NULL pointer to derived Analysis object!";
00094     return; 
00095   }
00096   
00097   if ( !hPeds_.first ) {
00098     anal->addErrorCode(sistrip::nullPtr_);
00099     return;
00100   }
00101   
00102   if ( !hNoise_.first ) {
00103     anal->addErrorCode(sistrip::nullPtr_);
00104     return;
00105   }
00106   
00107   TProfile* peds_histo = dynamic_cast<TProfile*>(hPeds_.first);
00108   TProfile* raw_histo = dynamic_cast<TProfile*>(hNoise_.first);
00109   
00110   if ( !peds_histo ) {
00111     anal->addErrorCode(sistrip::nullPtr_);
00112     return;
00113   }
00114 
00115   if ( !raw_histo ) {
00116     anal->addErrorCode(sistrip::nullPtr_);
00117     return;
00118   }
00119 
00120   if ( peds_histo->GetNbinsX() != 256 ) {
00121     anal->addErrorCode(sistrip::numberOfBins_);
00122     return;
00123   }
00124 
00125   if ( raw_histo->GetNbinsX() != 256 ) {
00126     anal->addErrorCode(sistrip::numberOfBins_);
00127     return;
00128   }
00129   
00130   // Iterate through APVs 
00131   for ( uint16_t iapv = 0; iapv < 2; iapv++ ) {
00132 
00133     // Used to calc mean and rms for peds and noise
00134     float p_sum = 0., p_sum2 = 0., p_max = -1.*sistrip::invalid_, p_min = sistrip::invalid_;
00135     float r_sum = 0., r_sum2 = 0., r_max = -1.*sistrip::invalid_, r_min = sistrip::invalid_;
00136 
00137     // Iterate through strips of APV
00138     for ( uint16_t istr = 0; istr < 128; istr++ ) {
00139 
00140       static uint16_t strip;
00141       strip = iapv*128 + istr;
00142 
00143       // Pedestals 
00144       if ( peds_histo ) {
00145         if ( peds_histo->GetBinEntries(strip+1) ) {
00146           anal->peds_[iapv][istr] = peds_histo->GetBinContent(strip+1);
00147           p_sum += anal->peds_[iapv][istr];
00148           p_sum2 += (anal->peds_[iapv][istr] * anal->peds_[iapv][istr]);
00149           if ( anal->peds_[iapv][istr] > p_max ) { p_max = anal->peds_[iapv][istr]; }
00150           if ( anal->peds_[iapv][istr] < p_min ) { p_min = anal->peds_[iapv][istr]; }
00151         }
00152       } 
00153       
00154       // Raw noise
00155       if ( !anal->legacy_ ) {
00156         if ( raw_histo ) {
00157           if ( raw_histo->GetBinEntries(strip+1) ) {
00158             anal->raw_[iapv][istr] = raw_histo->GetBinContent(strip+1);
00159             r_sum += anal->raw_[iapv][istr];
00160             r_sum2 += ( anal->raw_[iapv][istr] * anal->raw_[iapv][istr] );
00161             if ( anal->raw_[iapv][istr] > r_max ) { r_max = anal->raw_[iapv][istr]; }
00162             if ( anal->raw_[iapv][istr] < r_min ) { r_min = anal->raw_[iapv][istr]; }
00163           }
00164         }
00165       } else {
00166         if ( peds_histo ) {
00167           if ( peds_histo->GetBinEntries(strip+1) ) {
00168             anal->raw_[iapv][istr] = raw_histo->GetBinError(strip+1);
00169             r_sum += anal->raw_[iapv][istr];
00170             r_sum2 += ( anal->raw_[iapv][istr] * anal->raw_[iapv][istr] );
00171             if ( anal->raw_[iapv][istr] > r_max ) { r_max = anal->raw_[iapv][istr]; }
00172             if ( anal->raw_[iapv][istr] < r_min ) { r_min = anal->raw_[iapv][istr]; }
00173           }
00174         }
00175       }
00176       
00177     } // strip loop
00178     
00179     // Calc mean and rms for peds
00180     if ( !anal->peds_[iapv].empty() ) { 
00181       p_sum /= static_cast<float>( anal->peds_[iapv].size() );
00182       p_sum2 /= static_cast<float>( anal->peds_[iapv].size() );
00183       anal->pedsMean_[iapv] = p_sum;
00184       anal->pedsSpread_[iapv] = sqrt( fabs(p_sum2 - p_sum*p_sum) );
00185     }
00186     
00187     // Calc mean and rms for raw noise
00188     if ( !anal->raw_[iapv].empty() ) { 
00189       r_sum /= static_cast<float>( anal->raw_[iapv].size() );
00190       r_sum2 /= static_cast<float>( anal->raw_[iapv].size() );
00191       anal->rawMean_[iapv] = r_sum;
00192       anal->rawSpread_[iapv] = sqrt( fabs(r_sum2 - r_sum*r_sum) );
00193     }
00194     
00195     // Set max and min values for peds and raw noise
00196     if ( p_max > -1.*sistrip::maximum_ ) { anal->pedsMax_[iapv] = p_max; }
00197     if ( p_min < 1.*sistrip::maximum_ )  { anal->pedsMin_[iapv] = p_min; }
00198     if ( r_max > -1.*sistrip::maximum_ ) { anal->rawMax_[iapv] = r_max; }
00199     if ( r_min < 1.*sistrip::maximum_ )  { anal->rawMin_[iapv] = r_min; }
00200     
00201   } // apv loop
00202 
00203 }