00001 #include "SimCalorimetry/EcalSimAlgos/interface/EcalShapeBase.h" 00002 #include "FWCore/MessageLogger/interface/MessageLogger.h" 00003 00004 #include <cmath> 00005 #include <iostream> 00006 #include <fstream> 00007 #include <algorithm> 00008 00009 const double EcalShapeBase::qNSecPerBin = 1./(1.*kNBinsPerNSec) ; 00010 00011 EcalShapeBase::~EcalShapeBase() 00012 { 00013 delete m_derivPtr ; 00014 } 00015 00016 EcalShapeBase::EcalShapeBase( bool aSaveDerivative ) : 00017 m_firstIndexOverThreshold ( 0 ) , 00018 m_firstTimeOverThreshold ( 0.0 ) , 00019 m_indexOfMax ( 0 ) , 00020 m_timeOfMax ( 0.0 ) , 00021 m_shape ( DVec( kNBinsStored, 0.0 ) ) , 00022 m_derivPtr ( aSaveDerivative ? new DVec( kNBinsStored, 0.0 ) : 0 ) 00023 { 00024 } 00025 00026 double 00027 EcalShapeBase::timeOfThr() const 00028 { 00029 return m_firstTimeOverThreshold ; 00030 } 00031 00032 double 00033 EcalShapeBase::timeOfMax() const 00034 { 00035 return m_timeOfMax ; 00036 } 00037 00038 double 00039 EcalShapeBase::timeToRise() const 00040 { 00041 return timeOfMax() - timeOfThr() ; 00042 } 00043 00044 00045 void 00046 EcalShapeBase::buildMe() 00047 { 00048 DVec shapeArray( k1NSecBinsTotal , 0.0 ) ; 00049 00050 fillShape( shapeArray ) ; 00051 00052 const double maxel ( *max_element( shapeArray.begin(), shapeArray.end() ) ) ; 00053 00054 const double maxelt ( 1.e-5 < maxel ? maxel : 1 ) ; 00055 00056 for( unsigned int i ( 0 ) ; i != shapeArray.size(); ++i ) 00057 { 00058 shapeArray[i] = shapeArray[i]/maxelt ; 00059 } 00060 00061 const double thresh ( threshold()/maxelt ) ; 00062 00063 /* 00064 for( unsigned int i ( 0 ) ; i != k1NSecBinsTotal ; ++i ) 00065 { 00066 LogDebug("EcalShapeBase") << " time (ns) = " << (double)i << " tabulated ECAL pulse shape = " << shapeArray[i]; 00067 } 00068 */ 00069 00070 const double delta ( qNSecPerBin/2. ) ; 00071 00072 for( unsigned int j ( 0 ) ; j != kNBinsStored ; ++j ) 00073 { 00074 const double xb ( ( j + 0.5 )*qNSecPerBin ) ; 00075 00076 const unsigned int ibin ( j/kNBinsPerNSec ) ; 00077 00078 double value = 0.0 ; 00079 double deriv = 0.0 ; 00080 00081 if( 0 == ibin || 00082 shapeArray.size() == 1 + ibin ) // cannot do quadratic interpolation at ends 00083 { 00084 value = shapeArray[ibin]; 00085 deriv = 0.0 ; 00086 } 00087 else 00088 { 00089 const double x ( xb - ( ibin + 0.5 ) ) ; 00090 const double f1 ( shapeArray[ ibin - 1 ] ) ; 00091 const double f2 ( shapeArray[ ibin ] ) ; 00092 const double f3 ( shapeArray[ ibin + 1 ] ) ; 00093 const double a ( f2 ) ; 00094 const double b ( ( f3 - f1 )/2. ) ; 00095 const double c ( ( f1 + f3 )/2. - f2 ) ; 00096 value = a + b*x + c*x*x; 00097 deriv = ( b + 2*c*x )/delta ; 00098 } 00099 00100 m_shape[ j ] = value; 00101 if( 0 != m_derivPtr ) (*m_derivPtr)[ j ] = deriv; 00102 00103 if( 0 < j && 00104 thresh < value && 00105 0 == m_firstIndexOverThreshold ) 00106 { 00107 m_firstIndexOverThreshold = j - 1 ; 00108 m_firstTimeOverThreshold = m_firstIndexOverThreshold*qNSecPerBin ; 00109 } 00110 00111 if( m_shape[ m_indexOfMax ] < value ) 00112 { 00113 m_indexOfMax = j ; 00114 } 00115 00116 // LogDebug("EcalShapeBase") << " time (ns) = " << ( j + 1.0 )*qNSecPerBin - delta 00117 // << " interpolated ECAL pulse shape = " << m_shape[ j ] 00118 // << " derivative = " << ( 0 != m_derivPtr ? (*m_derivPtr)[ j ] : 0 ) ; 00119 } 00120 m_timeOfMax = m_indexOfMax*qNSecPerBin ; 00121 } 00122 00123 unsigned int 00124 EcalShapeBase::timeIndex( double aTime ) const 00125 { 00126 const int index ( m_firstIndexOverThreshold + 00127 (unsigned int) ( aTime*kNBinsPerNSec + 0.5 ) ) ; 00128 00129 const bool bad ( (int) m_firstIndexOverThreshold > index || 00130 (int) kNBinsStored <= index ) ; 00131 00132 if( (int) kNBinsStored <= index ) 00133 { 00134 LogDebug("EcalShapeBase") << " ECAL MGPA shape requested for out of range time " << aTime ; 00135 } 00136 return ( bad ? kNBinsStored : (unsigned int) index ) ; 00137 } 00138 00139 double 00140 EcalShapeBase::operator() ( double aTime ) const 00141 { 00142 // return pulse amplitude for request time in ns 00143 00144 const unsigned int index ( timeIndex( aTime ) ) ; 00145 return ( kNBinsStored == index ? 0 : m_shape[ index ] ) ; 00146 } 00147 00148 double 00149 EcalShapeBase::derivative( double aTime ) const 00150 { 00151 const unsigned int index ( timeIndex( aTime ) ) ; 00152 return ( 0 == m_derivPtr || 00153 kNBinsStored == index ? 0 : (*m_derivPtr)[ index ] ) ; 00154 }