CMS 3D CMS Logo

VpspScanAlgorithm.cc

Go to the documentation of this file.
00001 #include "DQM/SiStripCommissioningAnalysis/interface/VpspScanAlgorithm.h"
00002 #include "CondFormats/SiStripObjects/interface/VpspScanAnalysis.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 <sstream>
00010 #include <iomanip>
00011 #include <cmath>
00012 
00013 using namespace sistrip;
00014 
00015 // -----------------------------------------------------------------------------
00016 //
00017 VpspScanAlgorithm::VpspScanAlgorithm( VpspScanAnalysis* const anal )
00018   : CommissioningAlgorithm(anal),
00019     histos_( 2, Histo(0,"") )
00020 {;}
00021 
00022 // ----------------------------------------------------------------------------
00023 // 
00024 void VpspScanAlgorithm::extract( const std::vector<TH1*>& histos ) { 
00025 
00026   if ( !anal() ) {
00027     edm::LogWarning(mlCommissioning_)
00028       << "[VpspScanAlgorithm::" << __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 pointer
00046     if ( !(*ihis) ) { continue; }
00047     
00048     // Check name
00049     SiStripHistoTitle title( (*ihis)->GetName() );
00050     if ( title.runType() != sistrip::VPSP_SCAN ) {
00051       anal()->addErrorCode(sistrip::unexpectedTask_);
00052       continue;
00053     }
00054     
00055     // Extract APV number
00056     uint16_t apv = sistrip::invalid_; 
00057     if ( title.extraInfo().find(sistrip::apv_) != std::string::npos ) {
00058       std::stringstream ss;
00059       ss << title.extraInfo().substr( title.extraInfo().find(sistrip::apv_) + sistrip::apv_.size(), 1 );
00060       ss >> std::dec >> apv;
00061     }
00062 
00063     if ( apv <= 1 ) {
00064       histos_[apv].first = *ihis; 
00065       histos_[apv].second = (*ihis)->GetName();
00066     } else {
00067       anal()->addErrorCode(sistrip::unexpectedExtraInfo_);
00068     }
00069     
00070   }
00071 
00072 }
00073 
00074 // -----------------------------------------------------------------------------
00075 //
00076 void VpspScanAlgorithm::analyse() {
00077 
00078   if ( !anal() ) {
00079     edm::LogWarning(mlCommissioning_)
00080       << "[VpspScanAlgorithm::" << __func__ << "]"
00081       << " NULL pointer to base Analysis object!";
00082     return; 
00083   }
00084 
00085   CommissioningAnalysis* tmp = const_cast<CommissioningAnalysis*>( anal() );
00086   VpspScanAnalysis* anal = dynamic_cast<VpspScanAnalysis*>( tmp );
00087   if ( !anal ) {
00088     edm::LogWarning(mlCommissioning_)
00089       << "[VpspScanAlgorithm::" << __func__ << "]"
00090       << " NULL pointer to derived Analysis object!";
00091     return; 
00092   }
00093 
00094   // from deprecated()
00095   
00096   std::vector<const TProfile*> histos; 
00097   std::vector<uint16_t> monitorables;
00098   for ( uint16_t iapv = 0; iapv < 2; iapv++ ) {
00099     
00100     monitorables.clear();
00101     monitorables.resize( 7, sistrip::invalid_ );
00102     
00103     histos.clear();
00104     histos.push_back( const_cast<const TProfile*>( dynamic_cast<TProfile*>(histos_[iapv].first) ) );
00105     
00106     if ( !histos[0] ) {
00107       anal->addErrorCode(sistrip::nullPtr_);
00108       continue;
00109     }
00110 
00111     // Find "top" plateau
00112     int first = sistrip::invalid_;
00113     float top = -1. * sistrip::invalid_;;
00114     for ( int k = 5; k < 55; k++ ) {
00115       if ( !histos[0]->GetBinEntries(k) ) { continue; }
00116       if ( histos[0]->GetBinContent(k) >= top ) { 
00117         first = k; 
00118         top = histos[0]->GetBinContent(k); 
00119       }
00120     }
00121     if ( top < -1. * sistrip::valid_ ) { top = sistrip::invalid_; } //@@ just want +ve invalid number here
00122     if ( top > 1. * sistrip::valid_ ) { 
00123       anal->addErrorCode(sistrip::noTopPlateau_);
00124       continue;
00125     } 
00126     monitorables[5] = static_cast<uint16_t>(top);
00127     monitorables[3] = first;
00128         
00129     // Find "bottom" plateau
00130     int last = sistrip::invalid_;
00131     float bottom = 1. * sistrip::invalid_;
00132     for ( int k = 55; k > 5; k-- ) {
00133       if ( !histos[0]->GetBinEntries(k) ) { continue; }
00134       if ( histos[0]->GetBinContent(k) <= bottom ) { 
00135         last = k; 
00136         bottom = histos[0]->GetBinContent(k); 
00137       }
00138     }
00139     if ( bottom > 1. * sistrip::valid_ ) {
00140       anal->addErrorCode(sistrip::noBottomPlateau_);
00141       continue;
00142     } 
00143     monitorables[6] = static_cast<uint16_t>(bottom);
00144     monitorables[4] = last;
00145       
00146     // Set optimum baseline level
00147     float opt = bottom + ( top - bottom ) * 1./3.; 
00148     monitorables[1] = static_cast<uint16_t>(opt);
00149       
00150     // Find optimum VPSP setting 
00151     uint16_t vpsp = sistrip::invalid_;
00152     if ( opt < 1. * sistrip::valid_ ) {
00153       uint16_t ivpsp = 5; 
00154       for ( ; ivpsp < 55; ivpsp++ ) { 
00155         if ( histos[0]->GetBinContent(ivpsp) < opt ) { break; }
00156       }
00157       if ( ivpsp != 54 ) { 
00158         vpsp = ivpsp; 
00159         monitorables[0] = vpsp;
00160       }
00161       else { 
00162         anal->addErrorCode(sistrip::noVpspSetting_); 
00163         continue;
00164       }
00165         
00166     } else { 
00167       anal->addErrorCode(sistrip::noBaselineLevel_); 
00168       continue;
00169     }
00170     
00171     // Set analysis results for both APVs
00172     anal->vpsp_[iapv]        = monitorables[0];
00173     anal->adcLevel_[iapv]    = monitorables[1];
00174     anal->fraction_[iapv]    = monitorables[2];
00175     anal->topEdge_[iapv]     = monitorables[3];
00176     anal->bottomEdge_[iapv]  = monitorables[4];
00177     anal->topLevel_[iapv]    = monitorables[5];
00178     anal->bottomLevel_[iapv] = monitorables[6];
00179     
00180   }
00181   
00182 }

Generated on Tue Jun 9 17:33:27 2009 for CMSSW by  doxygen 1.5.4