00001 #include "CalibTracker/SiPixelSCurveCalibration/interface/SiPixelSCurveCalibrationAnalysis.h"
00002 #include "TMath.h"
00003
00004 #include <iostream>
00005 #include <fstream>
00006
00007 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00008 #include "DataFormats/SiPixelDigi/interface/SiPixelCalibDigiError.h"
00009 #include "CondFormats/SiPixelObjects/interface/SiPixelFrameConverter.h"
00010 #include "CondFormats/SiPixelObjects/interface/ElectronicIndex.h"
00011 #include "CondFormats/SiPixelObjects/interface/DetectorIndex.h"
00012 #include "CondFormats/SiPixelObjects/interface/LocalPixel.h"
00013 #include <sstream>
00014
00015
00016 std::vector<float> SiPixelSCurveCalibrationAnalysis::efficiencies_(0);
00017 std::vector<float> SiPixelSCurveCalibrationAnalysis::effErrors_(0);
00018
00019
00020 void SiPixelSCurveCalibrationAnalysis::calibrationEnd(){
00021 if(printoutthresholds_)
00022 makeThresholdSummary();
00023 }
00024
00025 void SiPixelSCurveCalibrationAnalysis::makeThresholdSummary(void){
00026 ofstream myfile;
00027 myfile.open (thresholdfilename_.c_str());
00028 for(detIDHistogramMap::iterator thisDetIdHistoGrams= histograms_.begin(); thisDetIdHistoGrams != histograms_.end(); ++thisDetIdHistoGrams){
00029
00030 const MonitorElement *sigmahist = (*thisDetIdHistoGrams).second[kSigmas];
00031 const MonitorElement *thresholdhist = (*thisDetIdHistoGrams).second[kThresholds];
00032 uint32_t detid = (*thisDetIdHistoGrams).first;
00033 std::string name = sigmahist->getTitle();
00034 std::string rocname = name.substr(0,name.size()-7);
00035 rocname+="_ROC";
00036 int total_rows = sigmahist ->getNbinsY();
00037 int total_columns = sigmahist->getNbinsX();
00038
00039 for (int irow=0; irow<total_rows; ++irow){
00040 for (int icol=0; icol<total_columns; ++icol){
00041 float threshold_error = sigmahist->getBinContent(icol+1,irow+1);
00042 if(writeZeroes_ ||(!writeZeroes_ && threshold_error>0)){
00043
00044 int realfedID=-1;
00045 for(int fedid=0; fedid<=40; ++fedid){
00046 SiPixelFrameConverter converter(theCablingMap_.product(),fedid);
00047 if(converter.hasDetUnit(detid)){
00048 realfedID=fedid;
00049 break;
00050 }
00051 }
00052 if (realfedID==-1){
00053 std::cout<<"error: could not obtain real fed ID"<<std::endl;
00054 }
00055 sipixelobjects::DetectorIndex detector ={detid,irow,icol};
00056 sipixelobjects::ElectronicIndex cabling;
00057 SiPixelFrameConverter formatter(theCablingMap_.product(),realfedID);
00058 formatter.toCabling(cabling,detector);
00059
00060
00061 sipixelobjects::LocalPixel::DcolPxid loc;
00062 loc.dcol = cabling.dcol;
00063 loc.pxid = cabling.pxid;
00064 const sipixelobjects::PixelFEDCabling *theFed= theCablingMap_.product()->fed(realfedID);
00065 const sipixelobjects::PixelFEDLink * link = theFed->link(cabling.link);
00066 const sipixelobjects::PixelROC *theRoc = link->roc(cabling.roc);
00067 sipixelobjects::LocalPixel locpixel(loc);
00068 int newrow= locpixel.rocRow();
00069 int newcol = locpixel.rocCol();
00070 myfile<<rocname<<theRoc->idInDetUnit()<<" "<<newcol<<" "<<newrow<<" "<<thresholdhist->getBinContent(icol+1, irow+1)<<" "<<threshold_error;
00071 myfile<<"\n";
00072 }
00073 }
00074 }
00075 }
00076 myfile.close();
00077 }
00078
00079
00080 void chi2toMinimize(int &npar, double* grad, double &fcnval, double* xval, int iflag)
00081 {
00082 TF1 * theFormula = SiPixelSCurveCalibrationAnalysis::fitFunction_;
00083
00084 for (int i = 0; i < npar; i++)
00085 theFormula->SetParameter(i, xval[i]);
00086 fcnval = 0;
00087
00088 const std::vector<short>* theVCalValues = SiPixelSCurveCalibrationAnalysis::getVcalValues();
00089 for (uint32_t i = 0; i < theVCalValues->size(); i++)
00090 {
00091 float chi = (SiPixelSCurveCalibrationAnalysis::efficiencies_[i] - theFormula->Eval((*theVCalValues)[i]) );
00092 chi /= SiPixelSCurveCalibrationAnalysis::effErrors_[i];
00093 fcnval += chi*chi;
00094 }
00095 }
00096
00097 void
00098 SiPixelSCurveCalibrationAnalysis::doSetup(const edm::ParameterSet& iConfig)
00099 {
00100 edm::LogInfo("SiPixelSCurveCalibrationAnalysis") << "Setting up calibration paramters.";
00101 std::vector<uint32_t> anEmptyDefaultVectorOfUInts;
00102 std::vector<uint32_t> detIDsToSaveVector_;
00103 useDetectorHierarchyFolders_ = iConfig.getUntrackedParameter<bool>("useDetectorHierarchyFolders", true);
00104 saveCurvesThatFlaggedBad_ = iConfig.getUntrackedParameter<bool>("saveCurvesThatFlaggedBad", false);
00105 detIDsToSaveVector_ = iConfig.getUntrackedParameter<std::vector<uint32_t> >("detIDsToSave", anEmptyDefaultVectorOfUInts);
00106 maxCurvesToSave_ = iConfig.getUntrackedParameter<uint32_t>("maxCurvesToSave", 1000);
00107 write2dHistograms_ = iConfig.getUntrackedParameter<bool>("write2dHistograms", true);
00108 write2dFitResult_ = iConfig.getUntrackedParameter<bool>("write2dFitResult", true);
00109 printoutthresholds_ = iConfig.getUntrackedParameter<bool>("writeOutThresholdSummary",true);
00110 thresholdfilename_ = iConfig.getUntrackedParameter<std::string>("thresholdOutputFileName","thresholds.txt");
00111 minimumChi2prob_ = iConfig.getUntrackedParameter<double>("minimumChi2prob", 0);
00112 minimumThreshold_ = iConfig.getUntrackedParameter<double>("minimumThreshold", -10);
00113 maximumThreshold_ = iConfig.getUntrackedParameter<double>("maximumThreshold", 300);
00114 minimumSigma_ = iConfig.getUntrackedParameter<double>("minimumSigma", 0);
00115 maximumSigma_ = iConfig.getUntrackedParameter<double>("maximumSigma", 100);
00116 minimumEffAsymptote_ = iConfig.getUntrackedParameter<double>("minimumEffAsymptote", 0);
00117 maximumEffAsymptote_ = iConfig.getUntrackedParameter<double>("maximumEffAsymptote", 1000);
00118 maximumSigmaBin_ = iConfig.getUntrackedParameter<double>("maximumSigmaBin", 10);
00119 maximumThresholdBin_ = iConfig.getUntrackedParameter<double>("maximumThresholdBin", 255);
00120
00121 writeZeroes_= iConfig.getUntrackedParameter<bool>("alsoWriteZeroThresholds", false);
00122
00123
00124 for(unsigned int i = 0; i < detIDsToSaveVector_.size(); i++)
00125 detIDsToSave_.insert( std::make_pair(detIDsToSaveVector_[i], true) );
00126 }
00127
00128 SiPixelSCurveCalibrationAnalysis::~SiPixelSCurveCalibrationAnalysis()
00129 {
00130
00131 }
00132
00133 void SiPixelSCurveCalibrationAnalysis::buildACurveHistogram(const uint32_t& detid, const uint32_t& row, const uint32_t& col, sCurveErrorFlag errorFlag, const std::vector<float>& efficiencies, const std::vector<float>& errors)
00134 {
00135 if (curvesSavedCounter_ > maxCurvesToSave_)
00136 {
00137 edm::LogWarning("SiPixelSCurveCalibrationAnalysis") << "WARNING: Request to save curve for [detid](col/row): [" << detid << "](" << col << "/" << row << ") denied. Maximum number of saved curves (defined in .cfi) exceeded.";
00138 return;
00139 }
00140 std::ostringstream rootName;
00141 rootName << "SCurve_row_" << row << "_col_" << col;
00142 std::ostringstream humanName;
00143 humanName << translateDetIdToString(detid) << "_" << rootName.str() << "_ErrorFlag_" << (int)errorFlag;
00144
00145 unsigned int numberOfVCalPoints = vCalPointsAsFloats_.size()-1;
00146 if (efficiencies.size() != numberOfVCalPoints || errors.size() != numberOfVCalPoints)
00147 {
00148 edm::LogError("SiPixelSCurveCalibrationAnalysis") << "Error saving single curve histogram! Number of Vcal values (" << numberOfVCalPoints << ") does not match number of efficiency points or error points!";
00149 return;
00150 }
00151 setDQMDirectory(detid);
00152 float * vcalValuesToPassToCrappyRoot = &vCalPointsAsFloats_[0];
00153 MonitorElement * aBadHisto = bookDQMHistogram1D(detid, rootName.str(), humanName.str(), numberOfVCalPoints, vcalValuesToPassToCrappyRoot);
00154 curvesSavedCounter_++;
00155 for(unsigned int iBin = 0; iBin < numberOfVCalPoints; ++iBin)
00156 {
00157 int rootBin = iBin + 1;
00158 aBadHisto->setBinContent(rootBin, efficiencies[iBin]);
00159 aBadHisto->setBinError(rootBin, errors[iBin]);
00160 }
00161 }
00162
00163 void SiPixelSCurveCalibrationAnalysis::calibrationSetup(const edm::EventSetup& iSetup)
00164 {
00165 edm::LogInfo("SiPixelSCurveCalibrationAnalysis") << "Calibration Settings: VCalLow: " << vCalValues_[0] << " VCalHigh: " << vCalValues_[vCalValues_.size()-1] << " nVCal: " << vCalValues_.size() << " nTriggers: " << nTriggers_;
00166 curvesSavedCounter_ = 0;
00167 if (saveCurvesThatFlaggedBad_)
00168 {
00169
00170 const std::vector<short>* theVCalValues = this->getVcalValues();
00171 unsigned int numberOfVCalPoints = theVCalValues->size();
00172 edm::LogWarning("SiPixelSCurveCalibrationAnalysis") << "WARNING: Option set to save indiviual S-Curves - max number: "
00173 << maxCurvesToSave_ << " This can lead to large memory consumption! (Got " << numberOfVCalPoints << " VCal Points";
00174 for(unsigned int i = 0; i < numberOfVCalPoints; i++)
00175 {
00176 vCalPointsAsFloats_.push_back( static_cast<float>((*theVCalValues)[i]) );
00177 edm::LogInfo("SiPixelSCurveCalibrationAnalysis") << "Adding calibration Vcal: " << (*theVCalValues)[i];
00178 }
00179
00180 vCalPointsAsFloats_.push_back( vCalPointsAsFloats_[numberOfVCalPoints-1] + 1 );
00181 }
00182
00183 fitFunction_ = new TF1("sCurve", "0.5*[2]*(1+TMath::Erf( (x-[0]) / ([1]*sqrt(2)) ) )", vCalValues_[0], vCalValues_[vCalValues_.size()-1]);
00184 }
00185
00186 bool
00187 SiPixelSCurveCalibrationAnalysis::checkCorrectCalibrationType()
00188 {
00189 if(calibrationMode_=="SCurve")
00190 return true;
00191 else if(calibrationMode_=="unknown"){
00192 edm::LogInfo("SiPixelSCurveCalibrationAnalysis") << "calibration mode is: " << calibrationMode_ << ", continuing anyway..." ;
00193 return true;
00194 }
00195 else{
00196
00197 }
00198 return false;
00199 }
00200
00201 sCurveErrorFlag SiPixelSCurveCalibrationAnalysis::estimateSCurveParameters(const std::vector<float>& eff, float& threshold, float& sigma)
00202 {
00203 sCurveErrorFlag output = errAllZeros;
00204 bool allZeroSoFar = true;
00205 int turnOnBin = -1;
00206 int saturationBin = -1;
00207 for (uint32_t iVcalPt = 0; iVcalPt < eff.size(); iVcalPt++)
00208 {
00209 if (allZeroSoFar && eff[iVcalPt] != 0 ) {
00210 turnOnBin = iVcalPt;
00211 allZeroSoFar = false;
00212 output = errNoTurnOn;
00213 } else if (eff[iVcalPt] > 0.90)
00214 {
00215 saturationBin = iVcalPt;
00216 short turnOnVcal = vCalValues_[turnOnBin];
00217 short saturationVcal = vCalValues_[saturationBin];
00218 short delta = saturationVcal - turnOnVcal;
00219 sigma = delta * 0.682;
00220 if (sigma < 1)
00221 sigma = 1;
00222 threshold = turnOnVcal + (0.5 * delta);
00223 return errOK;
00224 }
00225 }
00226 return output;
00227 }
00228
00229 sCurveErrorFlag SiPixelSCurveCalibrationAnalysis::fittedSCurveSanityCheck(float threshold, float sigma, float amplitude)
00230 {
00231
00232 if (threshold > vCalValues_[vCalValues_.size()-1] || threshold < vCalValues_[0] ||
00233 sigma > vCalValues_[vCalValues_.size()-1] - vCalValues_[0] )
00234 return errFitNonPhysical;
00235
00236 if (threshold < minimumThreshold_ || threshold > maximumThreshold_ ||
00237 sigma < minimumSigma_ || sigma > maximumSigma_ ||
00238 amplitude < minimumEffAsymptote_ || amplitude > maximumEffAsymptote_)
00239 return errFlaggedBadByUser;
00240
00241 return errOK;
00242 }
00243
00244 void calculateEffAndError(int nADCResponse, int nTriggers, float& eff, float& error)
00245 {
00246 eff = (float)nADCResponse / (float)nTriggers;
00247 double effForErrorCalculation = eff;
00248 if (eff <= 0 || eff >= 1)
00249 effForErrorCalculation = 0.5 / (double)nTriggers;
00250 error = TMath::Sqrt(effForErrorCalculation*(1-effForErrorCalculation) / (double)nTriggers);
00251 }
00252
00253
00254 void SiPixelSCurveCalibrationAnalysis::newDetID(uint32_t detid)
00255 {
00256 edm::LogInfo("SiPixelSCurveCalibrationAnalysis") << "Found a new DetID (" << detid << ")! Checking to make sure it has not been added.";
00257
00258 sCurveHistogramHolder tempMap;
00259 std::pair<detIDHistogramMap::iterator, bool> insertResult;
00260 insertResult = histograms_.insert(std::make_pair(detid, tempMap));
00261 if (insertResult.second)
00262 {
00263 edm::LogInfo("SiPixelSCurveCalibrationAnalysisHistogramReport") << "Histogram Map.insert() returned true! Booking new histogrames for detID: " << detid;
00264
00265 if (useDetectorHierarchyFolders_)
00266 setDQMDirectory(detid);
00267
00268 std::string detIdName = translateDetIdToString(detid);
00269 if (write2dHistograms_){
00270 MonitorElement * D2sigma = bookDQMHistoPlaquetteSummary2D(detid,"ScurveSigmas", detIdName + " Sigmas");
00271 MonitorElement * D2thresh = bookDQMHistoPlaquetteSummary2D(detid,"ScurveThresholds", detIdName + " Thresholds");
00272 MonitorElement * D2chi2 = bookDQMHistoPlaquetteSummary2D(detid,"ScurveChi2Prob",detIdName + " Chi2Prob");
00273 insertResult.first->second.insert(std::make_pair(kSigmas, D2sigma));
00274 insertResult.first->second.insert(std::make_pair(kThresholds, D2thresh));
00275 insertResult.first->second.insert(std::make_pair(kChi2s, D2chi2));
00276 }
00277 if (write2dFitResult_){
00278 MonitorElement * D2FitResult = bookDQMHistoPlaquetteSummary2D(detid,"ScurveFitResult", detIdName + " Fit Result");
00279 insertResult.first->second.insert(std::make_pair(kFitResults, D2FitResult));
00280 }
00281 MonitorElement * D1sigma = bookDQMHistogram1D(detid,"ScurveSigmasSummary", detIdName + " Sigmas Summary", 100, 0, maximumSigmaBin_);
00282 MonitorElement * D1thresh = bookDQMHistogram1D(detid,"ScurveThresholdSummary", detIdName + " Thresholds Summary", 255, 0, maximumThresholdBin_);
00283 MonitorElement * D1chi2 = bookDQMHistogram1D(detid,"ScurveChi2ProbSummary", detIdName + " Chi2Prob Summary", 101, 0, 1.01);
00284 MonitorElement * D1FitResult = bookDQMHistogram1D(detid,"ScurveFitResultSummary", detIdName + " Fit Result Summary", 10, -0.5, 9.5);
00285 insertResult.first->second.insert(std::make_pair(kSigmaSummary, D1sigma));
00286 insertResult.first->second.insert(std::make_pair(kThresholdSummary, D1thresh));
00287 insertResult.first->second.insert(std::make_pair(kChi2Summary, D1chi2));
00288 insertResult.first->second.insert(std::make_pair(kFitResultSummary, D1FitResult));
00289 }
00290 }
00291
00292 bool SiPixelSCurveCalibrationAnalysis::doFits(uint32_t detid, std::vector<SiPixelCalibDigi>::const_iterator calibDigi)
00293 {
00294 sCurveErrorFlag errorFlag = errOK;
00295 uint32_t nVCalPts = calibDigi->getnpoints();
00296
00297 efficiencies_.resize(0);
00298 effErrors_.resize(0);
00299 for (uint32_t iVcalPt = 0; iVcalPt < nVCalPts; iVcalPt++)
00300 {
00301 float eff;
00302 float error;
00303 calculateEffAndError(calibDigi->getnentries(iVcalPt), nTriggers_, eff, error);
00304 edm::LogInfo("SiPixelSCurveCalibrationAnalysis") << "Eff: " << eff << " Error: " << error << " nEntries: " << calibDigi->getnentries(iVcalPt) << " nTriggers: " << nTriggers_ << " VCalPt " << vCalValues_[iVcalPt];
00305 efficiencies_.push_back(eff);
00306 effErrors_.push_back(error);
00307 }
00308
00309
00310 float thresholdGuess;
00311 float sigmaGuess;
00312 errorFlag = estimateSCurveParameters(efficiencies_, thresholdGuess, sigmaGuess);
00313
00314
00315 Double_t sigma = -1.0;
00316 Double_t sigmaError = -1.0;
00317 Double_t threshold = -1.0;
00318 Double_t thresholdError = -1.0;
00319 Double_t amplitude = -1.0;
00320 Double_t amplitudeError = -1.0;
00321 Double_t chi2 = -1.0;
00322
00323 Int_t nDOF = vCalValues_.size() - 3;
00324 Double_t chi2probability = 0;
00325
00326 if (errorFlag == errOK)
00327 {
00328
00329 TMinuit *gMinuit = new TMinuit(3);
00330 gMinuit->SetPrintLevel(-1);
00331 gMinuit->SetFCN(chi2toMinimize);
00332
00333
00334 gMinuit->DefineParameter(0, "Threshold", (Double_t)thresholdGuess, 1, -50, 300);
00335
00336 gMinuit->DefineParameter(1, "Sigma", (Double_t)sigmaGuess, 0.1, 0, 255);
00337
00338 gMinuit->DefineParameter(2, "Amplitude", 1, 0.1, -0.001, 200);
00339
00340
00341 gMinuit->Migrad();
00342 gMinuit->GetParameter(0, threshold, thresholdError);
00343 gMinuit->GetParameter(1, sigma, sigmaError);
00344 gMinuit->GetParameter(2, amplitude, amplitudeError);
00345
00346
00347 Double_t params[3] = {threshold, sigma, amplitude};
00348 gMinuit->Eval(3, NULL, chi2, params, 0);
00349
00350 if (nDOF <= 0)
00351 chi2probability = 0;
00352 else
00353 chi2probability = TMath::Prob(chi2, nDOF);
00354
00355
00356 if (chi2probability > minimumChi2prob_)
00357 errorFlag = fittedSCurveSanityCheck(threshold, sigma, amplitude);
00358 else
00359 errorFlag = errBadChi2Prob;
00360
00361 edm::LogInfo("SiPixelSCurveCalibrationAnalysis") << "Fit finished with errorFlag: " << errorFlag << " - threshold: " << threshold << " sigma: " << sigma << " chi2: " << chi2 << " nDOF: " << nDOF << " chi2Prob: " << chi2probability << " chi2MinUser: " << minimumChi2prob_;
00362
00363 delete gMinuit;
00364 }
00365
00366 uint32_t row = calibDigi->row();
00367 uint32_t col = calibDigi->col();
00368
00369
00370 detIDHistogramMap::iterator thisDetIdHistoGrams;
00371 thisDetIdHistoGrams = histograms_.find(detid);
00372 if (thisDetIdHistoGrams != histograms_.end())
00373 {
00374 edm::LogInfo("SiPixelSCurveCalibrationAnalysisHistogramReport") << "Filling histograms for [detid](col/row): [" << detid << "](" << col << "/" << row << ") ErrorFlag: " << errorFlag;
00375
00376 (*thisDetIdHistoGrams).second[kFitResultSummary]->Fill(errorFlag);
00377 if (write2dFitResult_)
00378 (*thisDetIdHistoGrams).second[kFitResults]->setBinContent(col+1, row+1, errorFlag);
00379
00380
00381 (*thisDetIdHistoGrams).second[kSigmaSummary]->Fill(sigma);
00382 (*thisDetIdHistoGrams).second[kThresholdSummary]->Fill(threshold);
00383 if (write2dHistograms_)
00384 {
00385 (*thisDetIdHistoGrams).second[kSigmas]->setBinContent(col+1, row+1, sigma);
00386 (*thisDetIdHistoGrams).second[kThresholds]->setBinContent(col+1, row+1, threshold);
00387 }
00388
00389 (*thisDetIdHistoGrams).second[kChi2Summary]->Fill(chi2probability);
00390 if (write2dHistograms_)
00391 (*thisDetIdHistoGrams).second[kChi2s]->Fill(col, row, chi2probability);
00392 }
00393
00394 if (saveCurvesThatFlaggedBad_)
00395 {
00396 bool thisDetIDinList = false;
00397 if (detIDsToSave_.find(detid) != detIDsToSave_.end())
00398 thisDetIDinList = true;
00399
00400 if (errorFlag != errOK || thisDetIDinList)
00401 {
00402 edm::LogError("SiPixelSCurveCalibrationAnalysis") << "Saving error histogram for [detid](col/row): [" << detid << "](" << col << "/" << row << ") ErrorFlag: " << errorFlag;
00403 buildACurveHistogram(detid, row, col, errorFlag, efficiencies_, effErrors_);
00404 }
00405 }
00406
00407 return true;
00408
00409 }