CMS 3D CMS Logo

OptoScanAlgorithm Class Reference

#include <DQM/SiStripCommissioningAnalysis/interface/OptoScanAlgorithm.h>

Inheritance diagram for OptoScanAlgorithm:

CommissioningAlgorithm

List of all members.

Public Member Functions

Histo histo (const uint16_t &gain, const uint16_t &digital_level) const
 Histogram pointer and title.
 OptoScanAlgorithm (OptoScanAnalysis *const )
virtual ~OptoScanAlgorithm ()

Private Member Functions

void analyse ()
 Performs histogram anaysis.
void extract (const std::vector< TH1 * > &)
 Extracts and organises histograms.
 OptoScanAlgorithm ()

Private Attributes

std::vector< std::vector< Histo > > histos_
 Pointers and titles for histograms.


Detailed Description

Definition at line 17 of file OptoScanAlgorithm.h.


Constructor & Destructor Documentation

OptoScanAlgorithm::OptoScanAlgorithm ( OptoScanAnalysis * const  anal  ) 

Definition at line 17 of file OptoScanAlgorithm.cc.

00018   : CommissioningAlgorithm(anal),
00019     histos_( 4, std::vector<Histo>( 3, Histo(0,"") ) )
00020 {;}

virtual OptoScanAlgorithm::~OptoScanAlgorithm (  )  [inline, virtual]

Definition at line 23 of file OptoScanAlgorithm.h.

00023 {;}

OptoScanAlgorithm::OptoScanAlgorithm (  )  [inline, private]

Definition at line 31 of file OptoScanAlgorithm.h.

00031 {;}


Member Function Documentation

void OptoScanAlgorithm::analyse (  )  [private, virtual]

Performs histogram anaysis.

Implements CommissioningAlgorithm.

Definition at line 93 of file OptoScanAlgorithm.cc.

References sistrip::LinearFit::Params::a_, sistrip::LinearFit::add(), CommissioningAnalysis::addErrorCode(), CommissioningAlgorithm::anal(), sistrip::LinearFit::Params::b_, OptoScanAnalysis::baseSlope_, OptoScanAnalysis::bias_, OptoScanAnalysis::defaultBiasSetting_, OptoScanAnalysis::defaultGainSetting_, OptoScanAnalysis::fedAdcGain_, sistrip::LinearFit::fit(), OptoScanAnalysis::gain_, histos_, sistrip::invalid_, sistrip::invalidZeroLightLevel_, OptoScanAnalysis::liftOff_, OptoScanAnalysis::linkNoise_, max, sistrip::maximum_, OptoScanAnalysis::measGain_, min, sistrip::mlCommissioning_, sistrip::nullPtr_, sistrip::numberOfBins_, range, ss, OptoScanAnalysis::threshold_, OptoScanAnalysis::tickHeight_, tmp, sistrip::unexpectedBinNumber_, sistrip::valid_, width, and OptoScanAnalysis::zeroLight_.

00093                                 { 
00094   
00095   if ( !anal() ) {
00096     edm::LogWarning(mlCommissioning_)
00097       << "[OptoScanAlgorithm::" << __func__ << "]"
00098       << " NULL pointer to base Analysis object!";
00099     return; 
00100   }
00101   
00102   CommissioningAnalysis* tmp = const_cast<CommissioningAnalysis*>( anal() );
00103   OptoScanAnalysis* anal = dynamic_cast<OptoScanAnalysis*>( tmp );
00104   if ( !anal ) {
00105     edm::LogWarning(mlCommissioning_)
00106       << "[OptoScanAlgorithm::" << __func__ << "]"
00107       << " NULL pointer to derived Analysis object!";
00108     return; 
00109   }
00110   
00111   // Iterate through four gain settings
00112   for ( uint16_t igain = 0; igain < 4; igain++ ) {
00113     
00114     // Select histos appropriate for gain setting
00115     TH1* base_his = histos_[igain][0].first; 
00116     TH1* peak_his = histos_[igain][1].first;
00117     TH1* noise_his = histos_[igain][2].first;
00118 
00119     if ( !base_his ) {
00120       anal->addErrorCode(sistrip::nullPtr_);
00121       return;
00122     }
00123     
00124     if ( !peak_his ) {
00125       anal->addErrorCode(sistrip::nullPtr_);
00126       return;
00127     }
00128 
00129     if ( !noise_his ) {
00130       anal->addErrorCode(sistrip::nullPtr_);
00131       return;
00132     }
00133     
00134     TProfile* base_histo = dynamic_cast<TProfile*>(base_his);
00135     if ( !base_histo ) {
00136       anal->addErrorCode(sistrip::nullPtr_);
00137       return;
00138     }
00139     
00140     TProfile* peak_histo = dynamic_cast<TProfile*>(peak_his);
00141     if ( !peak_histo ) {
00142       anal->addErrorCode(sistrip::nullPtr_);
00143       return;
00144     }
00145     
00146     TProfile* noise_histo = dynamic_cast<TProfile*>(noise_his);
00147     if ( !noise_histo ) {
00148       anal->addErrorCode(sistrip::nullPtr_);
00149       return;
00150     }
00151 
00152     // Check histogram binning
00153     uint16_t nbins = static_cast<uint16_t>( peak_histo->GetNbinsX() );
00154     if ( static_cast<uint16_t>( base_histo->GetNbinsX() ) != nbins ) {
00155       anal->addErrorCode(sistrip::numberOfBins_);
00156       if ( static_cast<uint16_t>( base_histo->GetNbinsX() ) < nbins ) {
00157         nbins = static_cast<uint16_t>( base_histo->GetNbinsX() );
00158       }
00159     }
00160 
00161     // Some containers
00162     std::vector<float> peak_contents(0);
00163     std::vector<float> peak_errors(0);
00164     std::vector<float> peak_entries(0);
00165     std::vector<float> base_contents(0);
00166     std::vector<float> base_errors(0);
00167     std::vector<float> base_entries(0);
00168     std::vector<float> noise_contents(0);
00169     std::vector<float> noise_errors(0);
00170     std::vector<float> noise_entries(0);
00171     float peak_max = -1.*sistrip::invalid_;
00172     float peak_min =  1.*sistrip::invalid_;
00173     float base_max = -1.*sistrip::invalid_;
00174     float base_min =  1.*sistrip::invalid_;
00175     float noise_max = -1.*sistrip::invalid_;
00176     float noise_min =  1.*sistrip::invalid_;
00177 
00178     // Transfer histogram contents/errors/stats to containers
00179     for ( uint16_t ibin = 0; ibin < nbins; ibin++ ) {
00180 
00181       // Peak histogram
00182       peak_contents.push_back( peak_histo->GetBinContent(ibin+1) );
00183       peak_errors.push_back( peak_histo->GetBinError(ibin+1) );
00184       peak_entries.push_back( peak_histo->GetBinEntries(ibin+1) );
00185       if ( peak_entries[ibin] ) { 
00186         if ( peak_contents[ibin] > peak_max ) { peak_max = peak_contents[ibin]; }
00187         if ( peak_contents[ibin] < peak_min && ibin ) { peak_min = peak_contents[ibin]; }
00188       }
00189 
00190       // Base histogram
00191       base_contents.push_back( base_histo->GetBinContent(ibin+1) );
00192       base_errors.push_back( base_histo->GetBinError(ibin+1) );
00193       base_entries.push_back( base_histo->GetBinEntries(ibin+1) );
00194       if ( base_entries[ibin] ) { 
00195         if ( base_contents[ibin] > base_max ) { base_max = base_contents[ibin]; }
00196         if ( base_contents[ibin] < base_min && ibin ) { base_min = base_contents[ibin]; }
00197       }
00198 
00199       // Noise histogram
00200       noise_contents.push_back( noise_histo->GetBinContent(ibin+1) );
00201       noise_errors.push_back( noise_histo->GetBinError(ibin+1) );
00202       noise_entries.push_back( noise_histo->GetBinEntries(ibin+1) );
00203       if ( noise_entries[ibin] ) { 
00204         if ( noise_contents[ibin] > noise_max ) { noise_max = noise_contents[ibin]; }
00205         if ( noise_contents[ibin] < noise_min && ibin ) { noise_min = noise_contents[ibin]; }
00206       }
00207 
00208     }
00209     
00210     // Find "zero light" level and error
00211     //@@ record bias setting used for zero light level
00212     //@@ zero light error changes wrt gain setting ???
00213     float zero_light_level = sistrip::invalid_;
00214     float zero_light_error = sistrip::invalid_;
00215     for ( uint16_t ibin = 0; ibin < nbins; ibin++ ) {
00216       if ( base_entries[ibin] ) { 
00217         zero_light_level = base_contents[ibin];
00218         zero_light_error = base_errors[ibin];
00219         break;
00220       }
00221     }
00222 
00223     float zero_light_thres = sistrip::invalid_;
00224     if ( zero_light_level <= sistrip::maximum_ && 
00225          zero_light_error <= sistrip::maximum_ ) { 
00226       zero_light_thres = zero_light_level + 5. * zero_light_error;
00227     } else {
00228       std::stringstream ss;
00229       ss << sistrip::invalidZeroLightLevel_ << "ForGain" << igain;
00230       anal->addErrorCode( ss.str() );
00231       continue;
00232     }
00233     
00234     // Find range of base histogram
00235     float base_range = base_max - base_min;
00236 
00237     // Find overlapping max/min region that constrains range of linear fit
00238     float max = peak_max < base_max ? peak_max : base_max;
00239     float min = peak_min > base_min ? peak_min : base_min;
00240     float range = max - min;
00241 
00242     // Container identifying whether samples from 'base' histo are above "zero light" 
00243     std::vector<bool> above_zero_light;
00244     above_zero_light.resize(3,true);
00245     
00246     // Linear fits to top of peak and base curves and one to bottom of base curve
00247     sistrip::LinearFit peak_high;
00248     sistrip::LinearFit base_high;
00249     sistrip::LinearFit base_low;
00250 
00251     // Iterate through histogram bins
00252     uint16_t peak_bin = 0;
00253     uint16_t base_bin = 0;
00254     uint16_t low_bin = 0;
00255     for ( uint16_t ibin = 0; ibin < nbins; ibin++ ) {
00256       
00257       // Record whether consecutive samples from 'base' histo are above the "zero light" level
00258       if ( base_entries[ibin] ) {
00259         above_zero_light.erase( above_zero_light.begin() );
00260         if ( base_contents[ibin] > zero_light_thres ) { above_zero_light.push_back( true ); }
00261         else { above_zero_light.push_back( false ); }
00262         if ( above_zero_light.size() != 3 ) { above_zero_light.resize(3,false); } 
00263       }
00264       
00265       // Linear fit to peak histogram
00266       if ( peak_entries[ibin] &&
00267            peak_contents[ibin] > ( min + 0.2*range ) &&
00268            peak_contents[ibin] < ( min + 0.8*range ) ) {
00269         if ( !peak_bin ) { peak_bin = ibin; }
00270         if ( ( ibin - peak_bin ) < 10 ) { 
00271           peak_high.add( ibin, peak_contents[ibin] ); //@@ should weight using bin error or bin contents (sqrt(N)/N)
00272         }
00273       }
00274       // Linear fit to base histogram
00275       if ( base_entries[ibin] &&
00276            base_contents[ibin] > ( min + 0.2*range ) &&
00277            base_contents[ibin] < ( min + 0.8*range ) ) {
00278         if ( !base_bin ) { base_bin = ibin; }
00279         if ( ( ibin - base_bin ) < 10 ) { 
00280           base_high.add( ibin, base_contents[ibin] ); //@@ should weight using bin error or bin contents (sqrt(N)/N)
00281         }
00282       }
00283       // Low linear fit to base histogram
00284       if ( base_entries[ibin] &&
00285            //@@ above_zero_light[0] && above_zero_light[1] && above_zero_light[2] && 
00286            base_contents[ibin] > ( base_min + 0.2*base_range ) &&
00287            base_contents[ibin] < ( base_min + 0.6*base_range ) ) { 
00288         if ( !low_bin ) { low_bin = ibin; }
00289         if ( ( ibin - low_bin ) < 10 ) { 
00290           base_low.add( ibin, base_contents[ibin] ); //@@ should weight using bin error or bin contents (sqrt(N)/N)
00291         }
00292       }
00293       
00294     }
00295     
00296     // Extract width between two curves at midpoint within range
00297     float mid = min + 0.5*range;
00298     sistrip::LinearFit::Params peak_params;
00299     sistrip::LinearFit::Params base_params;
00300     peak_high.fit( peak_params );
00301     base_high.fit( base_params );
00302 
00303     float peak_pos = sistrip::invalid_;
00304     float base_pos = sistrip::invalid_;
00305     float width    = sistrip::invalid_;
00306     if ( peak_params.b_ > 0. ) {
00307       peak_pos = ( mid - peak_params.a_ ) / peak_params.b_;
00308     }
00309     if ( base_params.b_ > 0. ) {
00310       base_pos = ( mid - base_params.a_ ) / base_params.b_;
00311     }
00312     if ( base_pos < sistrip::valid_ &&
00313          peak_pos < sistrip::valid_ ) {
00314       width = base_pos - peak_pos;
00315     }
00316 
00317     // Extrapolate to zero light to find "lift off"
00318     sistrip::LinearFit::Params low_params;
00319     base_low.fit( low_params );
00320     float lift_off = sistrip::invalid_;
00321     if ( low_params.b_ > 0. ) {
00322       lift_off = ( zero_light_level - low_params.a_ ) / low_params.b_;
00323     }
00324     
00325     // ---------- Set all parameters ----------
00326     
00327     // Slope of baseline
00328     if ( low_params.b_ > 0. ) {
00329       anal->baseSlope_[igain] = low_params.b_;
00330     } 
00331     
00332     // Check "lift off" value and set bias setting accordingly
00333     if ( lift_off <= sistrip::maximum_ ) {
00334       anal->bias_[igain] = static_cast<uint16_t>( lift_off ) + 2;
00335     } else { anal->bias_[igain] = OptoScanAnalysis::defaultBiasSetting_; } 
00336     
00337     // Calculate "lift off" (in mA)
00338     if ( lift_off <= sistrip::maximum_ ) {
00339       anal->liftOff_[igain] = 0.45 * lift_off;
00340     }
00341     
00342     // Calculate laser threshold (in mA)
00343     if ( width < sistrip::invalid_ ) {
00344       anal->threshold_[igain] = 0.45 * ( lift_off - width/2. );
00345     }
00346 
00347     // Set "zero light" level
00348     anal->zeroLight_[igain] = zero_light_level;
00349     
00350     // Set link noise
00351     uint16_t bin_number = sistrip::invalid_;
00352     if ( anal->threshold_[igain] < sistrip::valid_ ) {
00353       bin_number = static_cast<uint16_t>( anal->threshold_[igain] / 0.45 ); 
00354     }
00355     if ( bin_number < sistrip::valid_ ) {
00356       if ( bin_number < noise_contents.size() ) { 
00357         anal->linkNoise_[igain] = noise_contents[bin_number]; 
00358       } else { anal->addErrorCode(sistrip::unexpectedBinNumber_); }
00359     }
00360       
00361     // Calculate tick mark height
00362     if ( low_params.b_ <= sistrip::maximum_ &&
00363          width <= sistrip::maximum_ ) {
00364       anal->tickHeight_[igain] = width * low_params.b_;
00365     }
00366       
00367     // Set measured gain 
00368     if ( anal->tickHeight_[igain] < sistrip::invalid_-1. ) {
00369       anal->measGain_[igain] = anal->tickHeight_[igain] * OptoScanAnalysis::fedAdcGain_ / 0.800;
00370     } else { anal->measGain_[igain] = sistrip::invalid_; }
00371       
00372   } // gain loop
00373 
00374   // Iterate through four gain settings and identify optimum gain setting
00375   const float target_gain = 0.8;
00376   float diff_in_gain = sistrip::invalid_;
00377   for ( uint16_t igain = 0; igain < 4; igain++ ) {
00378     
00379     // Check for sensible gain value
00380     if ( anal->measGain_[igain] > sistrip::maximum_ ) { continue; }
00381     
00382     // Find optimum gain setting
00383     if ( fabs( anal->measGain_[igain] - target_gain ) < diff_in_gain ) {
00384       anal->gain_ = igain;
00385       diff_in_gain = fabs( anal->measGain_[igain] - target_gain );
00386     }
00387     
00388   } 
00389   
00390   // Check optimum gain setting
00391   if ( anal->gain_ > sistrip::maximum_ ) { anal->gain_ = OptoScanAnalysis::defaultGainSetting_; }
00392   
00393 }

void OptoScanAlgorithm::extract ( const std::vector< TH1 * > &  histos  )  [private, virtual]

Extracts and organises histograms.

Implements CommissioningAlgorithm.

Definition at line 24 of file OptoScanAlgorithm.cc.

References CommissioningAnalysis::addErrorCode(), CommissioningAlgorithm::anal(), sistrip::extrainfo::baselineRms_, sistrip::extrainfo::digital_, CommissioningAlgorithm::extractFedKey(), CommissioningAnalysis::fedKey(), sistrip::extrainfo::gain_, histos_, sistrip::invalid_, sistrip::mlCommissioning_, sistrip::numberOfHistos_, sistrip::OPTO_SCAN, ss, indexGen::title, sistrip::unexpectedExtraInfo_, and sistrip::unexpectedTask_.

00024                                                                { 
00025 
00026   if ( !anal() ) {
00027     edm::LogWarning(mlCommissioning_)
00028       << "[ApvTimingAlgorithm::" << __func__ << "]"
00029       << " NULL pointer to Analysis object!";
00030     return; 
00031   }
00032 
00033   // Check number of histograms
00034   if ( histos.size() != 12 ) {
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 name
00049     SiStripHistoTitle title( (*ihis)->GetName() );
00050     if ( title.runType() != sistrip::OPTO_SCAN ) {
00051       anal()->addErrorCode(sistrip::unexpectedTask_);
00052       continue;
00053     }
00054 
00055     // Extract gain setting and digital high/low info
00056     uint16_t gain = sistrip::invalid_; 
00057     if ( title.extraInfo().find(sistrip::extrainfo::gain_) != std::string::npos ) {
00058       std::stringstream ss;
00059       ss << title.extraInfo().substr( title.extraInfo().find(sistrip::extrainfo::gain_) + sistrip::extrainfo::gain_.size(), 1 );
00060       ss >> std::dec >> gain;
00061     }
00062     uint16_t digital = sistrip::invalid_; 
00063     if ( title.extraInfo().find(sistrip::extrainfo::digital_) != std::string::npos ) {
00064       std::stringstream ss;
00065       ss << title.extraInfo().substr( title.extraInfo().find(sistrip::extrainfo::digital_) + sistrip::extrainfo::digital_.size(), 1 );
00066       ss >> std::dec >> digital;
00067     }
00068     bool baseline_rms = false;
00069     if ( title.extraInfo().find(sistrip::extrainfo::baselineRms_) != std::string::npos ) {
00070       baseline_rms = true;
00071     }
00072     
00073     if ( gain <= 3 ) { 
00074       if ( digital <= 1 ) {
00075         histos_[gain][digital].first = *ihis; 
00076         histos_[gain][digital].second = (*ihis)->GetName();
00077       } else if ( baseline_rms ) {
00078         histos_[gain][2].first = *ihis; 
00079         histos_[gain][2].second = (*ihis)->GetName();
00080       } else {
00081         anal()->addErrorCode(sistrip::unexpectedExtraInfo_);
00082       }
00083     } else {
00084       anal()->addErrorCode(sistrip::unexpectedExtraInfo_);
00085     }
00086     
00087   }
00088 
00089 }

CommissioningAlgorithm::Histo OptoScanAlgorithm::histo ( const uint16_t &  gain,
const uint16_t &  digital_level 
) const

Histogram pointer and title.

Definition at line 397 of file OptoScanAlgorithm.cc.

References histos_.

00398                                                                                               {
00399   if ( gain <= 3 && digital_level <= 1 ) { return histos_[gain][digital_level]; }
00400   else { return Histo(0,""); }
00401 }


Member Data Documentation

std::vector< std::vector<Histo> > OptoScanAlgorithm::histos_ [private]

Pointers and titles for histograms.

Definition at line 42 of file OptoScanAlgorithm.h.

Referenced by analyse(), extract(), and histo().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:29:19 2009 for CMSSW by  doxygen 1.5.4