CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoLocalTracker/SiStripClusterizer/src/ThreeThresholdAlgorithm.cc

Go to the documentation of this file.
00001 #include "RecoLocalTracker/SiStripClusterizer/interface/ThreeThresholdAlgorithm.h"
00002 #include "DataFormats/SiStripDigi/interface/SiStripDigi.h"
00003 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
00004 #include <cmath>
00005 #include <numeric>
00006 
00007 ThreeThresholdAlgorithm::
00008 ThreeThresholdAlgorithm(float chan, float seed, float cluster, unsigned holes, unsigned bad, unsigned adj, std::string qL) 
00009   : ChannelThreshold( chan ), SeedThreshold( seed ), ClusterThresholdSquared( cluster*cluster ),
00010     MaxSequentialHoles( holes ), MaxSequentialBad( bad ), MaxAdjacentBad( adj ) {
00011   qualityLabel = (qL);
00012   ADCs.reserve(128);
00013 }
00014 
00015 template<class digiDetSet>
00016 inline
00017 void ThreeThresholdAlgorithm::
00018 clusterizeDetUnit_(const digiDetSet& digis, output_t::FastFiller& output) {
00019   if( !isModuleUsable( digis.detId() )) return;
00020   setDetId( digis.detId() );
00021   
00022   typename digiDetSet::const_iterator  
00023     scan( digis.begin() ), 
00024     end(  digis.end() );
00025 
00026   clearCandidate();
00027   while( scan != end ) {
00028     while( scan != end  && !candidateEnded( scan->strip() ) ) 
00029       addToCandidate(*scan++);
00030     endCandidate(output);
00031   }
00032 }
00033 
00034 inline 
00035 bool ThreeThresholdAlgorithm::
00036 candidateEnded(const uint16_t& testStrip) const {
00037   uint16_t holes = testStrip - lastStrip - 1;
00038   return ( !ADCs.empty() &&                    // a candidate exists, and
00039            holes > MaxSequentialHoles &&       // too many holes if not all are bad strips, and
00040            ( holes > MaxSequentialBad ||       // (too many bad strips anyway, or 
00041              !allBadBetween( lastStrip, testStrip ))); // not all holes are bad strips)
00042 }
00043 
00044 inline 
00045 void ThreeThresholdAlgorithm::
00046 addToCandidate(const SiStripDigi& digi) { 
00047   float Noise = noise( digi.strip() );
00048   if( bad(digi.strip()) || digi.adc() < static_cast<uint16_t>( Noise * ChannelThreshold))
00049     return;
00050 
00051   if(candidateLacksSeed) candidateLacksSeed  =  digi.adc() < static_cast<uint16_t>( Noise * SeedThreshold);
00052   if(ADCs.empty()) lastStrip = digi.strip() - 1; // begin candidate
00053   while( ++lastStrip < digi.strip() ) ADCs.push_back(0); // pad holes
00054 
00055   ADCs.push_back( digi.adc() );
00056   noiseSquared += Noise*Noise;
00057 }
00058 
00059 template <class T>
00060 inline
00061 void ThreeThresholdAlgorithm::
00062 endCandidate(T& out) {
00063   if(candidateAccepted()) {
00064     applyGains();
00065     appendBadNeighbors();
00066     out.push_back(SiStripCluster(currentId(), firstStrip(), ADCs.begin(), ADCs.end()));
00067   }
00068   clearCandidate();  
00069 }
00070 
00071 inline 
00072 bool ThreeThresholdAlgorithm::
00073 candidateAccepted() const {
00074   return ( !candidateLacksSeed &&
00075            noiseSquared * ClusterThresholdSquared
00076            <=  std::pow( std::accumulate(ADCs.begin(),ADCs.end(),float(0)), 2));
00077 }
00078 
00079 inline
00080 void ThreeThresholdAlgorithm::
00081 applyGains() {
00082   uint16_t strip = firstStrip();
00083   for( std::vector<uint16_t>::iterator adc = ADCs.begin();  adc != ADCs.end();  adc++) {
00084     if(*adc > 255) throw InvalidChargeException( SiStripDigi(strip,*adc) );
00085     if(*adc > 253) continue; //saturated, do not scale
00086     uint16_t charge = static_cast<uint16_t>( *adc/gain(strip++) + 0.5 ); //adding 0.5 turns truncation into rounding
00087     *adc = ( charge > 1022 ? 255 : 
00088            ( charge >  253 ? 254 : charge ));
00089   }
00090 }
00091 
00092 inline 
00093 void ThreeThresholdAlgorithm::
00094 appendBadNeighbors() {
00095   uint8_t max = MaxAdjacentBad;
00096   while(0 < max--) {
00097     if( bad( firstStrip()-1) ) { ADCs.insert( ADCs.begin(), 0);  }
00098     if( bad(  lastStrip + 1) ) { ADCs.push_back(0); lastStrip++; }
00099   }
00100 }
00101 
00102 
00103 void ThreeThresholdAlgorithm::clusterizeDetUnit(const    edm::DetSet<SiStripDigi>& digis, output_t::FastFiller& output) {clusterizeDetUnit_(digis,output);}
00104 void ThreeThresholdAlgorithm::clusterizeDetUnit(const edmNew::DetSet<SiStripDigi>& digis, output_t::FastFiller& output) {clusterizeDetUnit_(digis,output);}
00105 
00106 inline
00107 bool ThreeThresholdAlgorithm::
00108 stripByStripBegin(uint32_t id) {
00109   if( !isModuleUsable( id )) return false;
00110   setDetId( id );
00111   clearCandidate();
00112   return true;
00113 }
00114 
00115 inline
00116 void ThreeThresholdAlgorithm::
00117 stripByStripAdd(uint16_t strip, uint16_t adc, std::vector<SiStripCluster>& out) {
00118   if(candidateEnded(strip))
00119     endCandidate(out);
00120   addToCandidate(SiStripDigi(strip,adc));
00121 }
00122 
00123 inline
00124 void ThreeThresholdAlgorithm::
00125 stripByStripEnd(std::vector<SiStripCluster>& out) { 
00126   endCandidate(out);
00127 }