00001 #include "DQM/SiStripCommissioningAnalysis/interface/CalibrationAlgorithm.h"
00002 #include "CondFormats/SiStripObjects/interface/CalibrationAnalysis.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 "TH1.h"
00010 #include "TVirtualFitter.h"
00011 #include <iostream>
00012 #include <sstream>
00013 #include <iomanip>
00014 #include <cmath>
00015
00016 using namespace sistrip;
00017
00018
00019
00020 CalibrationAlgorithm::CalibrationAlgorithm( CalibrationAnalysis* const anal )
00021 : CommissioningAlgorithm(anal),
00022 deconv_fitter_(0),
00023 peak_fitter_(0),
00024 cal_(0)
00025 {
00026 deconv_fitter_ = new TF1("deconv_fitter",fdeconv_convoluted,-50,50,5);
00027 deconv_fitter_->FixParameter(0,0);
00028 deconv_fitter_->SetParLimits(1,-100,0);
00029 deconv_fitter_->SetParLimits(2,0,200);
00030 deconv_fitter_->SetParLimits(3,5,100);
00031 deconv_fitter_->FixParameter(3,50);
00032 deconv_fitter_->SetParLimits(4,0,50);
00033 deconv_fitter_->SetParameters(0.,-10,0.96,50,20);
00034 peak_fitter_ = new TF1("peak_fitter",fpeak_convoluted,-50,50,5);
00035 peak_fitter_->FixParameter(0,0);
00036 peak_fitter_->SetParLimits(1,-100,0);
00037 peak_fitter_->SetParLimits(2,0,400);
00038 peak_fitter_->SetParLimits(3,5,100);
00039 peak_fitter_->FixParameter(3,50);
00040 peak_fitter_->SetParLimits(4,0,50);
00041 peak_fitter_->SetParameters(0.,-10,0.96,50,20);
00042 }
00043
00044
00045
00046 void CalibrationAlgorithm::extract( const std::vector<TH1*>& histos) {
00047 LogDebug(mlCommissioning_) << "##CD## in extract" ;
00048 if ( !anal() ) {
00049 edm::LogWarning(mlCommissioning_)
00050 << "[CalibrationAlgorithm::" << __func__ << "]"
00051 << " NULL pointer to base Analysis object!";
00052 return;
00053 }
00054
00055 CommissioningAnalysis* tmp = const_cast<CommissioningAnalysis*>( anal() );
00056 cal_ = dynamic_cast<CalibrationAnalysis*>( tmp );
00057 if ( !cal_ ) {
00058 edm::LogWarning(mlCommissioning_)
00059 << "[CalibrationAlgorithm::" << __func__ << "]"
00060 << " NULL pointer to derived Analysis object!";
00061 return;
00062 }
00063
00064
00065 if ( histos.size() != 32 && histos.size() !=2 ) {
00066 cal_->addErrorCode(sistrip::numberOfHistos_);
00067 }
00068
00069
00070 if ( !histos.empty() ) { cal_->fedKey( extractFedKey( histos.front() ) ); }
00071
00072
00073 std::vector<TH1*>::const_iterator ihis = histos.begin();
00074 unsigned int cnt = 0;
00075 for ( ; ihis != histos.end(); ihis++,cnt++ ) {
00076
00077
00078 if ( !(*ihis) ) { continue; }
00079
00080
00081 SiStripHistoTitle title( (*ihis)->GetName() );
00082 if ( title.runType() != sistrip::CALIBRATION &&
00083 title.runType() != sistrip::CALIBRATION_DECO &&
00084 title.runType() != sistrip::CALIBRATION_SCAN &&
00085 title.runType() != sistrip::CALIBRATION_SCAN_DECO ) {
00086 cal_->addErrorCode(sistrip::unexpectedTask_);
00087 continue;
00088 }
00089
00090 cal_->isScan_ = ( title.runType() == sistrip::CALIBRATION_SCAN ||
00091 title.runType() == sistrip::CALIBRATION_SCAN_DECO );
00092
00093
00094 histo_[cnt].first = *ihis;
00095 histo_[cnt].second = (*ihis)->GetTitle();
00096 LogDebug(mlCommissioning_) << "##CD## in extract got histogram " << cnt;
00097 }
00098
00099 }
00100
00101
00102
00103 void CalibrationAlgorithm::analyse() {
00104
00105 if ( !cal_ ) {
00106 edm::LogWarning(mlCommissioning_)
00107 << "[CalibrationAlgorithm::" << __func__ << "]"
00108 << " NULL pointer to derived Analysis object!";
00109 return;
00110 }
00111
00112 float Amean[2] = {0.,0.};
00113 float Amin[2] = {2000.,2000.};
00114 float Amax[2] = {0.,0.};
00115 float Aspread[2] = {0.,0.};
00116 float Tmean[2] = {0.,0.};
00117 float Tmin[2] = {2000.,2000.};
00118 float Tmax[2] = {0.,0.};
00119 float Tspread[2] = {0.,0.};
00120 float Rmean[2] = {0.,0.};
00121 float Rmin[2] = {2000.,2000.};
00122 float Rmax[2] = {0.,0.};
00123 float Rspread[2] = {0.,0.};
00124 float Cmean[2] = {0.,0.};
00125 float Cmin[2] = {2000.,2000.};
00126 float Cmax[2] = {0.,0.};
00127 float Cspread[2] = {0.,0.};
00128 float Smean[2] = {0.,0.};
00129 float Smin[2] = {2000.,2000.};
00130 float Smax[2] = {0.,0.};
00131 float Sspread[2] = {0.,0.};
00132 float Kmean[2] = {0.,0.};
00133 float Kmin[2] = {2000000.,2000000.};
00134 float Kmax[2] = {0.,0.};
00135 float Kspread[2] = {0.,0.};
00136
00137 unsigned int upperLimit = cal_->isScan_ ? 2 : 32;
00138 float nStrips = cal_->isScan_ ? 1. : 16.;
00139 for(unsigned int i=0;i<upperLimit;++i) {
00140 if ( !histo_[i].first ) {
00141 edm::LogWarning(mlCommissioning_)
00142 << " NULL pointer to histogram " << i << "!";
00143 return;
00144 }
00145
00146
00147 std::string title = histo_[i].first->GetName();
00148 int apv = 0;
00149 int strip = 0;
00150 if(title.find("STRIP")!=std::string::npos && title.find("Apv")!=std::string::npos) {
00151 strip = atoi(title.c_str()+title.find("STRIP")+6);
00152 apv = (atoi(title.c_str()+title.find("Apv")+3))%2;
00153 } else {
00154 strip = (atoi(title.c_str()+title.find_last_of("_")+1))%16;
00155 apv = (atoi(title.c_str()+title.find_last_of("_")+1))/16;
00156 if(title.find("Apv")!=std::string::npos) {
00157 apv = (atoi(title.c_str()+title.find("Apv")+3))%2;
00158 strip = strip*8 + cal_->calchan_;
00159 } else {
00160 edm::LogWarning(mlCommissioning_)
00161 << " Malformed histogram title! Strip/APV not retreived: "
00162 << title;
00163 }
00164 }
00165
00166 LogDebug(mlCommissioning_)
00167 << "start the calibration analysis for APV " << apv << " strip " << strip;
00168
00169 correctDistribution(histo_[i].first);
00170
00171
00172 cal_->amplitude_[apv][strip] = histo_[i].first->GetMaximum();
00173
00174
00175 cal_->riseTime_[apv][strip] = maximum(histo_[i].first) - turnOn(histo_[i].first);
00176
00177
00178 int lastBin = histo_[i].first->FindBin(histo_[i].first->GetBinCenter(histo_[i].first->GetMaximumBin())+125);
00179 if(lastBin>histo_[i].first->GetNbinsX()-4) lastBin = histo_[i].first->GetNbinsX()-4;
00180 if(histo_[i].first->GetMaximum()!=0)
00181 cal_->tail_[apv][strip] = 100*histo_[i].first->GetBinContent(lastBin)/histo_[i].first->GetMaximum();
00182 else
00183 cal_->tail_[apv][strip] = 100;
00184
00185
00186 TF1* fit = fitPulse(histo_[i].first);
00187
00188
00189 cal_->timeConstant_[apv][strip] = fit->GetParameter(3);
00190
00191
00192 cal_->smearing_[apv][strip] = fit->GetParameter(4);
00193
00194
00195 cal_->chi2_[apv][strip] = fit->GetChisquare();
00196
00197 LogDebug(mlCommissioning_)
00198 << "Results: "
00199 << cal_->chi2_[apv][strip] << " "
00200 << cal_->smearing_[apv][strip] << " "
00201 << cal_->timeConstant_[apv][strip] << " "
00202 << cal_->riseTime_[apv][strip] << " "
00203 << cal_->tail_[apv][strip] << " "
00204 << cal_->amplitude_[apv][strip];
00205
00206
00207 Amean[apv] += cal_->amplitude_[apv][strip]/nStrips;
00208 Amin[apv] = Amin[apv]<cal_->amplitude_[apv][strip] ? Amin[apv] : cal_->amplitude_[apv][strip];
00209 Amax[apv] = Amax[apv]>cal_->amplitude_[apv][strip] ? Amax[apv] : cal_->amplitude_[apv][strip];
00210 Aspread[apv] += cal_->amplitude_[apv][strip]*cal_->amplitude_[apv][strip]/nStrips;
00211 Tmean[apv] += cal_->tail_[apv][strip]/nStrips;
00212 Tmin[apv] = Tmin[apv]<cal_->tail_[apv][strip] ? Tmin[apv] : cal_->tail_[apv][strip];
00213 Tmax[apv] = Tmax[apv]>cal_->tail_[apv][strip] ? Tmax[apv] : cal_->tail_[apv][strip];
00214 Tspread[apv] += cal_->tail_[apv][strip]*cal_->tail_[apv][strip]/nStrips;
00215 Rmean[apv] += cal_->riseTime_[apv][strip]/nStrips;
00216 Rmin[apv] = Rmin[apv]<cal_->riseTime_[apv][strip] ? Rmin[apv] : cal_->riseTime_[apv][strip];
00217 Rmax[apv] = Rmax[apv]>cal_->riseTime_[apv][strip] ? Rmax[apv] : cal_->riseTime_[apv][strip];
00218 Rspread[apv] += cal_->riseTime_[apv][strip]*cal_->riseTime_[apv][strip]/nStrips;
00219 Cmean[apv] += cal_->timeConstant_[apv][strip]/nStrips;
00220 Cmin[apv] = Cmin[apv]<cal_->timeConstant_[apv][strip] ? Cmin[apv] : cal_->timeConstant_[apv][strip];
00221 Cmax[apv] = Cmax[apv]>cal_->timeConstant_[apv][strip] ? Cmax[apv] : cal_->timeConstant_[apv][strip];
00222 Cspread[apv] += cal_->timeConstant_[apv][strip]*cal_->timeConstant_[apv][strip]/nStrips;
00223 Smean[apv] += cal_->smearing_[apv][strip]/nStrips;
00224 Smin[apv] = Smin[apv]<cal_->smearing_[apv][strip] ? Smin[apv] : cal_->smearing_[apv][strip];
00225 Smax[apv] = Smax[apv]>cal_->smearing_[apv][strip] ? Smax[apv] : cal_->smearing_[apv][strip];
00226 Sspread[apv] += cal_->smearing_[apv][strip]*cal_->smearing_[apv][strip]/nStrips;
00227 Kmean[apv] += cal_->chi2_[apv][strip]/nStrips;
00228 Kmin[apv] = Kmin[apv]<cal_->chi2_[apv][strip] ? Kmin[apv] : cal_->chi2_[apv][strip];
00229 Kmax[apv] = Kmax[apv]>cal_->chi2_[apv][strip] ? Kmax[apv] : cal_->chi2_[apv][strip];
00230 Kspread[apv] += cal_->chi2_[apv][strip]*cal_->chi2_[apv][strip]/nStrips;
00231 }
00232
00233
00234 for(int i=0;i<2;++i) {
00235 cal_->mean_amplitude_[i] = Amean[i];
00236 cal_->mean_tail_[i] = Tmean[i];
00237 cal_->mean_riseTime_[i] = Rmean[i];
00238 cal_->mean_timeConstant_[i] = Cmean[i];
00239 cal_->mean_smearing_[i] = Smean[i];
00240 cal_->mean_chi2_[i] = Kmean[i];
00241 cal_->min_amplitude_[i] = Amin[i];
00242 cal_->min_tail_[i] = Tmin[i];
00243 cal_->min_riseTime_[i] = Rmin[i];
00244 cal_->min_timeConstant_[i] = Cmin[i];
00245 cal_->min_smearing_[i] = Smin[i];
00246 cal_->min_chi2_[i] = Kmin[i];
00247 cal_->max_amplitude_[i] = Amax[i];
00248 cal_->max_tail_[i] = Tmax[i];
00249 cal_->max_riseTime_[i] = Rmax[i];
00250 cal_->max_timeConstant_[i] = Cmax[i];
00251 cal_->max_smearing_[i] = Smax[i];
00252 cal_->max_chi2_[i] = Kmax[i];
00253 cal_->spread_amplitude_[i] = sqrt(fabs(Aspread[i]-Amean[i]*Amean[i]));
00254 cal_->spread_tail_[i] = sqrt(fabs(Tspread[i]-Tmean[i]*Tmean[i]));
00255 cal_->spread_riseTime_[i] = sqrt(fabs(Rspread[i]-Rmean[i]*Rmean[i]));
00256 cal_->spread_timeConstant_[i] = sqrt(fabs(Cspread[i]-Cmean[i]*Cmean[i]));
00257 cal_->spread_smearing_[i] = sqrt(fabs(Sspread[i]-Smean[i]*Smean[i]));
00258 cal_->spread_chi2_[i] = sqrt(fabs(Kspread[i]-Kmean[i]*Kmean[i]));
00259 }
00260 }
00261
00262
00263
00264 void CalibrationAlgorithm::correctDistribution( TH1* histo ) const
00265 {
00266
00267 histo->Scale(-1);
00268 if ( cal_ ) { if( cal_->isScan_ ) histo->Scale(1/16.); }
00269 }
00270
00271
00272
00273 TF1* CalibrationAlgorithm::fitPulse( TH1* histo,
00274 float rangeLow,
00275 float rangeHigh )
00276 {
00277 if(!cal_) return 0;
00278 if(!histo) return 0;
00279 float noise = 4.;
00280 float N = round(histo->GetMaximum()/125.);
00281 float error = sqrt(2*N)*noise;
00282
00283 for(int i=1;i<=histo->GetNbinsX();++i) {
00284 histo->SetBinError(i,error);
00285 }
00286 if (cal_->deconv_) {
00287 if(rangeLow>rangeHigh)
00288 histo->Fit(deconv_fitter_,"Q");
00289 else
00290 histo->Fit(deconv_fitter_,"0Q","",rangeLow,rangeHigh);
00291 return deconv_fitter_;
00292 } else {
00293 if(rangeLow>rangeHigh)
00294 histo->Fit(peak_fitter_,"Q");
00295 else
00296 histo->Fit(peak_fitter_,"0Q","",rangeLow,rangeHigh);
00297 return peak_fitter_;
00298 }
00299 }
00300
00301
00302
00303 float CalibrationAlgorithm::maximum( TH1* h ) {
00304 int bin = h->GetMaximumBin();
00305
00306 TF1* fit = fitPulse(h,h->GetBinCenter(bin)-25,h->GetBinCenter(bin)+25);
00307 return fit->GetMaximumX();
00308 }
00309
00310
00311
00312 float CalibrationAlgorithm::turnOn( TH1* h ) {
00313
00314 int bin=1;
00315 float amplitude = h->GetMaximum();
00316 for(; bin<= h->GetNbinsX() && h->GetBinContent(bin)<0.4*amplitude; ++bin) {}
00317 float end = h->GetBinLowEdge(bin);
00318
00319 TF1* sigmoid = new TF1("sigmoid","[0]/(1+exp(-[1]*(x-[2])))+[3]",0,end);
00320 sigmoid->SetParLimits(0,amplitude/10.,amplitude);
00321 sigmoid->SetParLimits(1,0.05,0.5);
00322 sigmoid->SetParLimits(2,end-10,end+10);
00323 sigmoid->SetParLimits(3,-amplitude/10.,amplitude/10.);
00324 sigmoid->SetParameters(amplitude/2.,0.1,end,0.);
00325 h->Fit(sigmoid,"0QR");
00326
00327 float time = 0.;
00328 float base = sigmoid->GetMinimum(0,end);
00329 for(;time<end && (sigmoid->Eval(time)-base)<0.03*(amplitude-base);time += 0.1) {}
00330 delete sigmoid;
00331 return time-0.05;
00332 }
00333