00001 #include "DQM/SiStripCommissioningAnalysis/interface/SamplingAlgorithm.h"
00002 #include "CondFormats/SiStripObjects/interface/SamplingAnalysis.h"
00003 #include "DataFormats/SiStripCommon/interface/SiStripHistoTitle.h"
00004 #include "DataFormats/SiStripCommon/interface/SiStripEnumsAndStrings.h"
00005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00006 #include "DQM/SiStripCommissioningAnalysis/interface/SiStripPulseShape.h"
00007 #include "TProfile.h"
00008 #include "TF1.h"
00009 #include <iostream>
00010 #include <sstream>
00011 #include <iomanip>
00012 #include <cmath>
00013
00014 using namespace sistrip;
00015
00016
00017
00018 SamplingAlgorithm::SamplingAlgorithm( const edm::ParameterSet & pset, SamplingAnalysis* const anal, uint32_t latencyCode )
00019 : CommissioningAlgorithm(anal),
00020 histo_(0,""),
00021 deconv_fitter_(0),
00022 peak_fitterA_(0),
00023 peak_fitterB_(0),
00024 latencyCode_(latencyCode),
00025 samp_(0)
00026 {
00027 peak_fitterA_ = new TF1("peak_fitterA",fpeak_convoluted,-4800,0,5);
00028 peak_fitterA_->SetNpx(2000);
00029 peak_fitterA_->FixParameter(0,0);
00030 peak_fitterA_->SetParLimits(1,0,4800);
00031 peak_fitterA_->SetParLimits(2,0,20);
00032 peak_fitterA_->FixParameter(3,50);
00033 peak_fitterA_->SetParLimits(4,0,75);
00034 peak_fitterA_->SetParameters(0.,1250,10,50,10);
00035
00036 peak_fitterB_ = new TF1("peak_fitterB",fpeak_convoluted,-100,100,5);
00037 peak_fitterB_->SetNpx(200);
00038 peak_fitterB_->FixParameter(0,0);
00039 peak_fitterB_->SetParLimits(1,-100,100);
00040 peak_fitterB_->SetParLimits(2,0,20);
00041 peak_fitterB_->FixParameter(3,50);
00042 peak_fitterB_->SetParLimits(4,0,75);
00043 peak_fitterB_->SetParameters(0.,-50,10,50,10);
00044
00045 deconv_fitter_ = new TF1("deconv_fitter",fdeconv_convoluted,-50,50,5);
00046 deconv_fitter_->SetNpx(1000);
00047 deconv_fitter_->FixParameter(0,0);
00048 deconv_fitter_->SetParLimits(1,-50,50);
00049 deconv_fitter_->SetParLimits(2,0,200);
00050 deconv_fitter_->SetParLimits(3,5,100);
00051 deconv_fitter_->FixParameter(3,50);
00052 deconv_fitter_->SetParLimits(4,0,75);
00053 deconv_fitter_->SetParameters(0.,-2.82,0.96,50,20);
00054 }
00055
00056
00057
00058 void SamplingAlgorithm::extract( const std::vector<TH1*>& histos) {
00059
00060 if ( !anal() ) {
00061 edm::LogWarning(mlCommissioning_)
00062 << "[SamplingAlgorithm::" << __func__ << "]"
00063 << " NULL pointer to Analysis object!";
00064 return;
00065 }
00066
00067 CommissioningAnalysis* tmp = const_cast<CommissioningAnalysis*>( anal() );
00068 samp_ = dynamic_cast<SamplingAnalysis*>( tmp );
00069 if ( !samp_ ) {
00070 edm::LogWarning(mlCommissioning_)
00071 << "[SamplingAlgorithm::" << __func__ << "]"
00072 << " NULL pointer to derived Analysis object!";
00073 return;
00074 }
00075
00076
00077 if ( histos.size() != 1 && histos.size() != 2 ) {
00078 samp_->addErrorCode(sistrip::numberOfHistos_);
00079 }
00080
00081
00082 if ( !histos.empty() ) { samp_->fedKey( extractFedKey( histos.front() ) ); }
00083
00084
00085 std::vector<TH1*>::const_iterator ihis = histos.begin();
00086 for ( ; ihis != histos.end(); ihis++ ) {
00087
00088
00089 if ( !(*ihis) ) {
00090 edm::LogWarning(mlCommissioning_) << " NULL pointer to histogram!";
00091 continue;
00092 }
00093
00094
00095 SiStripHistoTitle title( (*ihis)->GetName() );
00096 if ( title.runType() != sistrip::APV_LATENCY && title.runType() != sistrip::FINE_DELAY) {
00097 samp_->addErrorCode(sistrip::unexpectedTask_);
00098 continue;
00099 }
00100
00101 samp_->runType_ = title.runType();
00102
00103
00104 samp_->granularity_ = title.granularity();
00105
00106
00107 if ( title.extraInfo().find(sistrip::extrainfo::clusterCharge_) != std::string::npos ) {
00108 histo_.first = *ihis;
00109 histo_.second = (*ihis)->GetName();
00110 }
00111
00112 }
00113
00114 }
00115
00116
00117
00118 void SamplingAlgorithm::analyse() {
00119
00120 if ( !samp_ ) {
00121 edm::LogWarning(mlCommissioning_)
00122 << "[SamplingAlgorithm::" << __func__ << "]"
00123 << " NULL pointer to derived Analysis object!";
00124 return;
00125 }
00126
00127 TProfile* prof = (TProfile*)(histo_.first);
00128 if ( !prof ) {
00129 edm::LogWarning(mlCommissioning_)
00130 << " NULL pointer to histogram!" ;
00131 return;
00132 }
00133
00134
00135 prof->SetErrorOption(" ");
00136
00137
00138
00139 for(int i=0;i<prof->GetNbinsX();++i) {
00140 if(prof->GetBinEntries(i)>0)
00141 prof->SetBinError(i,prof->GetBinError(i)/sqrt(prof->GetBinEntries(i)));
00142 }
00143
00144
00145 pruneProfile(prof);
00146
00147
00148 correctBinning(prof);
00149
00150
00151 correctProfile(prof,samp_->sOnCut_);
00152
00153
00154 if(samp_->runType_==sistrip::APV_LATENCY) {
00155
00156
00157 float max = prof->GetBinCenter(prof->GetMaximumBin());
00158 float ampl = prof->GetMaximum();
00159 peak_fitterA_->SetParameters(0.,50-max,ampl/20.,50,10);
00160
00161
00162 if(prof->Fit(peak_fitterA_,"Q")==0)
00163 prof->Fit(peak_fitterA_,"QEM");
00164
00165
00166 samp_->max_ = peak_fitterA_->GetMaximumX();
00167 samp_->error_ = peak_fitterA_->GetParError(1);
00168
00169 } else {
00170
00171
00172 float max = prof->GetBinCenter(prof->GetMaximumBin());
00173 float ampl = prof->GetMaximum();
00174 deconv_fitter_->SetParameters(0.,-max,ampl/10.,50,20);
00175 peak_fitterB_->SetParameters(0.,50-max,ampl/20.,50,10);
00176 if(latencyCode_&0x80) {
00177
00178 if(prof->Fit(deconv_fitter_,"Q")==0)
00179 prof->Fit(deconv_fitter_,"QEM");
00180
00181 samp_->max_ = deconv_fitter_->GetMaximumX();
00182 samp_->error_ = deconv_fitter_->GetParError(1);
00183 } else {
00184
00185 if(prof->Fit(peak_fitterB_,"Q")==0)
00186 prof->Fit(peak_fitterB_,"QEM");
00187
00188 samp_->max_ = peak_fitterB_->GetMaximumX();
00189 samp_->error_ = peak_fitterB_->GetParError(1);
00190 }
00191
00192 }
00193
00194 }
00195
00196
00197
00198 void SamplingAlgorithm::pruneProfile(TProfile* profile) const
00199 {
00200
00201 uint32_t nbins=profile->GetNbinsX();
00202 uint32_t max = 0;
00203 for(uint32_t bin=1;bin<=nbins;++bin) max = max < profile->GetBinEntries(bin) ? uint32_t(profile->GetBinEntries(bin)) : max;
00204
00205 uint32_t min = max/10;
00206 for(uint32_t bin=1;bin<=nbins;++bin)
00207 if(profile->GetBinEntries(bin)<min) {
00208 profile->SetBinContent(bin,0.);
00209 profile->SetBinError(bin,0.);
00210 }
00211 }
00212
00213
00214
00215 void SamplingAlgorithm::correctBinning(TProfile* prof) const
00216 {
00217 prof->GetXaxis()->SetLimits(prof->GetXaxis()->GetXmin()-prof->GetBinWidth(1)/2.,
00218 prof->GetXaxis()->GetXmax()-prof->GetBinWidth(1)/2.);
00219 }
00220
00221
00222
00223 void SamplingAlgorithm::correctProfile(TProfile* profile, float SoNcut) const
00224 {
00225 if ( !samp_ ) { return; }
00226 uint32_t nbins=profile->GetNbinsX();
00227 float min = samp_->limit(SoNcut);
00228 for(uint32_t bin=1;bin<=nbins;++bin)
00229 if(profile->GetBinContent(bin)<min) {
00230 profile->SetBinContent(bin,0.);
00231 profile->SetBinError(bin,0.);
00232 profile->SetBinEntries(bin,0);
00233 }
00234 else {
00235 profile->SetBinContent(bin,
00236 profile->GetBinEntries(bin)*
00237 samp_->correctMeasurement(profile->GetBinContent(bin),SoNcut));
00238 }
00239 }