CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/CalibTracker/SiPixelSCurveCalibration/src/SiPixelSCurveCalibrationAnalysis.cc

Go to the documentation of this file.
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 //initialize static members
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    // loop over det id (det id = number (unsigned int) of pixel module 
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     //loop over all rows on columns on all ROCs 
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); // +1 because root bins start at 1
00042         if(writeZeroes_ ||(!writeZeroes_ && threshold_error>0)){             
00043           //changing from offline to online numbers
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           // cabling should now contain cabling.roc and cabling.dcol  and cabling.pxid
00060           // however, the coordinates now need to be converted from dcl,pxid to the row,col coordinates used in the calibration info 
00061           sipixelobjects::LocalPixel::DcolPxid loc;
00062           loc.dcol = cabling.dcol;
00063           loc.pxid = cabling.pxid;
00064           // FIX to adhere to new cabling map. To be replaced with CalibTracker/SiPixelTools detid - > hardware id classes ASAP.
00065           //        const sipixelobjects::PixelFEDCabling *theFed= theCablingMap.product()->fed(realfedID);
00066           //        const sipixelobjects::PixelFEDLink * link = theFed->link(cabling.link);
00067           //        const sipixelobjects::PixelROC *theRoc = link->roc(cabling.roc);
00068           sipixelobjects::LocalPixel locpixel(loc);
00069           sipixelobjects::CablingPathToDetUnit path = {static_cast<unsigned int>(realfedID),
00070                                                        static_cast<unsigned int>(cabling.link),
00071                                                        static_cast<unsigned int>(cabling.roc)};
00072           const sipixelobjects::PixelROC *theRoc = theCablingMap_->findItem(path);
00073           // END of FIX
00074           int newrow= locpixel.rocRow();
00075           int newcol = locpixel.rocCol();
00076           myfile<<rocname<<theRoc->idInDetUnit()<<" "<<newcol<<" "<<newrow<<" "<<thresholdhist->getBinContent(icol+1, irow+1)<<" "<<threshold_error;  // +1 because root bins start at 1
00077           myfile<<"\n";
00078         }
00079       }
00080     }
00081   }
00082   myfile.close();
00083 }
00084 
00085 //used for TMinuit fitting
00086 void chi2toMinimize(int &npar, double* grad, double &fcnval, double* xval, int iflag)
00087 {
00088    TF1 * theFormula = SiPixelSCurveCalibrationAnalysis::fitFunction_;
00089    //setup function parameters
00090    for (int i = 0; i < npar; i++)
00091       theFormula->SetParameter(i, xval[i]);
00092    fcnval = 0;
00093    //compute Chi2 of all points
00094    const std::vector<short>* theVCalValues = SiPixelSCurveCalibrationAnalysis::getVcalValues();
00095    for (uint32_t i = 0; i < theVCalValues->size(); i++)
00096    {
00097       float chi = (SiPixelSCurveCalibrationAnalysis::efficiencies_[i] - theFormula->Eval((*theVCalValues)[i]) );
00098       chi       /= SiPixelSCurveCalibrationAnalysis::effErrors_[i];
00099       fcnval += chi*chi;
00100    }
00101 }
00102 
00103 void
00104 SiPixelSCurveCalibrationAnalysis::doSetup(const edm::ParameterSet& iConfig)
00105 {
00106    edm::LogInfo("SiPixelSCurveCalibrationAnalysis") << "Setting up calibration paramters.";
00107    std::vector<uint32_t>        anEmptyDefaultVectorOfUInts;
00108    std::vector<uint32_t>        detIDsToSaveVector_;
00109    useDetectorHierarchyFolders_ = iConfig.getUntrackedParameter<bool>("useDetectorHierarchyFolders", true);
00110    saveCurvesThatFlaggedBad_    = iConfig.getUntrackedParameter<bool>("saveCurvesThatFlaggedBad", false);
00111    detIDsToSaveVector_          = iConfig.getUntrackedParameter<std::vector<uint32_t> >("detIDsToSave", anEmptyDefaultVectorOfUInts);
00112    maxCurvesToSave_             = iConfig.getUntrackedParameter<uint32_t>("maxCurvesToSave", 1000);
00113    write2dHistograms_           = iConfig.getUntrackedParameter<bool>("write2dHistograms", true);
00114    write2dFitResult_            = iConfig.getUntrackedParameter<bool>("write2dFitResult", true);
00115    printoutthresholds_          = iConfig.getUntrackedParameter<bool>("writeOutThresholdSummary",true);
00116    thresholdfilename_           = iConfig.getUntrackedParameter<std::string>("thresholdOutputFileName","thresholds.txt");  
00117    minimumChi2prob_             = iConfig.getUntrackedParameter<double>("minimumChi2prob", 0);
00118    minimumThreshold_            = iConfig.getUntrackedParameter<double>("minimumThreshold", -10);
00119    maximumThreshold_            = iConfig.getUntrackedParameter<double>("maximumThreshold", 300);
00120    minimumSigma_                = iConfig.getUntrackedParameter<double>("minimumSigma", 0);
00121    maximumSigma_                = iConfig.getUntrackedParameter<double>("maximumSigma", 100);
00122    minimumEffAsymptote_         = iConfig.getUntrackedParameter<double>("minimumEffAsymptote", 0);
00123    maximumEffAsymptote_         = iConfig.getUntrackedParameter<double>("maximumEffAsymptote", 1000);
00124    maximumSigmaBin_             = iConfig.getUntrackedParameter<double>("maximumSigmaBin", 10);
00125    maximumThresholdBin_         = iConfig.getUntrackedParameter<double>("maximumThresholdBin", 255);
00126 
00127    writeZeroes_= iConfig.getUntrackedParameter<bool>("alsoWriteZeroThresholds", false);
00128 
00129    // convert the vector into a map for quicker lookups.
00130    for(unsigned int i = 0; i < detIDsToSaveVector_.size(); i++)
00131       detIDsToSave_.insert( std::make_pair(detIDsToSaveVector_[i], true) );
00132 }
00133 
00134 SiPixelSCurveCalibrationAnalysis::~SiPixelSCurveCalibrationAnalysis()
00135 {
00136    //do nothing
00137 }
00138 
00139 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)
00140 {
00141    if (curvesSavedCounter_ > maxCurvesToSave_)
00142    {
00143       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.";
00144       return;
00145    }
00146    std::ostringstream rootName;
00147    rootName << "SCurve_row_" << row << "_col_" << col;
00148    std::ostringstream humanName;
00149    humanName << translateDetIdToString(detid) << "_" << rootName.str() << "_ErrorFlag_" << (int)errorFlag;
00150 
00151    unsigned int numberOfVCalPoints = vCalPointsAsFloats_.size()-1; //minus one is necessary since the lower edge of the last bin must be added
00152    if (efficiencies.size() != numberOfVCalPoints || errors.size() != numberOfVCalPoints)
00153    {
00154       edm::LogError("SiPixelSCurveCalibrationAnalysis") << "Error saving single curve histogram!  Number of Vcal values (" << numberOfVCalPoints << ") does not match number of efficiency points or error points!";
00155       return;
00156    }
00157    setDQMDirectory(detid);
00158    float * vcalValuesToPassToCrappyRoot = &vCalPointsAsFloats_[0];
00159    MonitorElement * aBadHisto = bookDQMHistogram1D(detid, rootName.str(), humanName.str(), numberOfVCalPoints, vcalValuesToPassToCrappyRoot);  //ROOT only takes an input as array. :(  HOORAY FOR CINT!
00160    curvesSavedCounter_++;
00161    for(unsigned int iBin = 0; iBin < numberOfVCalPoints; ++iBin)
00162    {
00163       int rootBin = iBin + 1;  //root bins start at 1
00164       aBadHisto->setBinContent(rootBin, efficiencies[iBin]);
00165       aBadHisto->setBinError(rootBin, errors[iBin]);
00166    }
00167 }
00168 
00169 void SiPixelSCurveCalibrationAnalysis::calibrationSetup(const edm::EventSetup& iSetup)
00170 {
00171    edm::LogInfo("SiPixelSCurveCalibrationAnalysis") << "Calibration Settings: VCalLow: " << vCalValues_[0] << "  VCalHigh: " << vCalValues_[vCalValues_.size()-1] << " nVCal: " << vCalValues_.size() << "  nTriggers: " << nTriggers_;
00172    curvesSavedCounter_ = 0;
00173    if (saveCurvesThatFlaggedBad_)
00174    {
00175       //build the vCal values as a vector of floats if we want to save single curves
00176       const std::vector<short>* theVCalValues = this->getVcalValues();
00177       unsigned int numberOfVCalPoints = theVCalValues->size();
00178       edm::LogWarning("SiPixelSCurveCalibrationAnalysis") << "WARNING: Option set to save indiviual S-Curves - max number: " 
00179                                                           << maxCurvesToSave_ << " This can lead to large memory consumption! (Got " << numberOfVCalPoints << " VCal Points";
00180       for(unsigned int i = 0; i < numberOfVCalPoints; i++)
00181       {
00182          vCalPointsAsFloats_.push_back( static_cast<float>((*theVCalValues)[i]) );
00183          edm::LogInfo("SiPixelSCurveCalibrationAnalysis") << "Adding calibration Vcal: " << (*theVCalValues)[i];
00184       }
00185       // must add lower edge of last bin to the vector
00186       vCalPointsAsFloats_.push_back( vCalPointsAsFloats_[numberOfVCalPoints-1] + 1 );
00187    }
00188 
00189    fitFunction_ = new TF1("sCurve", "0.5*[2]*(1+TMath::Erf( (x-[0]) / ([1]*sqrt(2)) ) )", vCalValues_[0], vCalValues_[vCalValues_.size()-1]);
00190 }
00191 
00192 bool
00193 SiPixelSCurveCalibrationAnalysis::checkCorrectCalibrationType()
00194 {
00195   if(calibrationMode_=="SCurve")
00196     return true;
00197   else if(calibrationMode_=="unknown"){
00198     edm::LogInfo("SiPixelSCurveCalibrationAnalysis") <<  "calibration mode is: " << calibrationMode_ << ", continuing anyway..." ;
00199     return true;
00200   }
00201   else{
00202     //    edm::LogDebug("SiPixelSCurveCalibrationAnalysis") << "unknown calibration mode for SCurves, should be \"SCurve\" and is \"" << calibrationMode_ << "\"";
00203   }
00204   return false;
00205 }
00206 
00207 sCurveErrorFlag SiPixelSCurveCalibrationAnalysis::estimateSCurveParameters(const std::vector<float>& eff, float& threshold, float& sigma)
00208 {
00209    sCurveErrorFlag output = errAllZeros;
00210    bool allZeroSoFar    = true;
00211    int turnOnBin        = -1;
00212    int saturationBin    = -1;
00213    for (uint32_t iVcalPt = 0; iVcalPt < eff.size(); iVcalPt++)
00214    {
00215       if (allZeroSoFar && eff[iVcalPt] != 0 ) {
00216          turnOnBin = iVcalPt;
00217          allZeroSoFar = false;
00218          output = errNoTurnOn;
00219       } else if (eff[iVcalPt] > 0.90)
00220       {
00221          saturationBin  = iVcalPt;
00222          short turnOnVcal       = vCalValues_[turnOnBin];
00223          short saturationVcal   = vCalValues_[saturationBin];
00224          short delta            = saturationVcal - turnOnVcal;
00225          sigma                  = delta * 0.682;
00226          if (sigma < 1)         //check to make sure sigma guess is larger than our X resolution.  Hopefully prevents Minuit from getting stuck at boundary
00227             sigma = 1;
00228          threshold              = turnOnVcal + (0.5 * delta);
00229          return errOK;
00230       }
00231    }
00232    return output;
00233 }
00234 
00235 sCurveErrorFlag SiPixelSCurveCalibrationAnalysis::fittedSCurveSanityCheck(float threshold, float sigma, float amplitude)
00236 {
00237    //check if nonsensical
00238    if (threshold > vCalValues_[vCalValues_.size()-1] || threshold < vCalValues_[0] ||
00239          sigma > vCalValues_[vCalValues_.size()-1] - vCalValues_[0] )
00240       return errFitNonPhysical;
00241 
00242    if (threshold < minimumThreshold_ || threshold > maximumThreshold_ ||
00243          sigma < minimumSigma_ || sigma > maximumSigma_ ||
00244          amplitude < minimumEffAsymptote_ || amplitude > maximumEffAsymptote_)
00245       return errFlaggedBadByUser;
00246 
00247    return errOK;
00248 }
00249 
00250 void calculateEffAndError(int nADCResponse, int nTriggers, float& eff, float& error)
00251 {
00252    eff = (float)nADCResponse / (float)nTriggers;
00253    double effForErrorCalculation = eff;
00254    if (eff <= 0 || eff >= 1)
00255       effForErrorCalculation = 0.5 / (double)nTriggers;
00256    error = TMath::Sqrt(effForErrorCalculation*(1-effForErrorCalculation) / (double)nTriggers);
00257 }
00258 
00259 //book histograms when new DetID is encountered in Event Record
00260 void SiPixelSCurveCalibrationAnalysis::newDetID(uint32_t detid)
00261 {
00262    edm::LogInfo("SiPixelSCurveCalibrationAnalysis") << "Found a new DetID (" << detid << ")!  Checking to make sure it has not been added.";
00263    //ensure that this DetID has not been added yet
00264    sCurveHistogramHolder tempMap;
00265    std::pair<detIDHistogramMap::iterator, bool> insertResult; 
00266    insertResult = histograms_.insert(std::make_pair(detid, tempMap));
00267    if (insertResult.second)     //indicates successful insertion
00268    {
00269       edm::LogInfo("SiPixelSCurveCalibrationAnalysisHistogramReport") << "Histogram Map.insert() returned true!  Booking new histogrames for detID: " << detid;
00270       // use detector hierarchy folders if desired
00271       if (useDetectorHierarchyFolders_)
00272          setDQMDirectory(detid);
00273 
00274       std::string detIdName = translateDetIdToString(detid);
00275       if (write2dHistograms_){
00276         MonitorElement * D2sigma       = bookDQMHistoPlaquetteSummary2D(detid,"ScurveSigmas", detIdName + " Sigmas");
00277         MonitorElement * D2thresh      = bookDQMHistoPlaquetteSummary2D(detid,"ScurveThresholds", detIdName + " Thresholds");
00278         MonitorElement * D2chi2        = bookDQMHistoPlaquetteSummary2D(detid,"ScurveChi2Prob",detIdName + " Chi2Prob");
00279          insertResult.first->second.insert(std::make_pair(kSigmas, D2sigma));
00280          insertResult.first->second.insert(std::make_pair(kThresholds, D2thresh));
00281          insertResult.first->second.insert(std::make_pair(kChi2s, D2chi2));
00282       }
00283       if (write2dFitResult_){
00284         MonitorElement * D2FitResult = bookDQMHistoPlaquetteSummary2D(detid,"ScurveFitResult", detIdName + " Fit Result");
00285          insertResult.first->second.insert(std::make_pair(kFitResults, D2FitResult));
00286       }
00287       MonitorElement * D1sigma       = bookDQMHistogram1D(detid,"ScurveSigmasSummary", detIdName + " Sigmas Summary", 100, 0, maximumSigmaBin_);
00288       MonitorElement * D1thresh      = bookDQMHistogram1D(detid,"ScurveThresholdSummary", detIdName + " Thresholds Summary", 255, 0, maximumThresholdBin_);
00289       MonitorElement * D1chi2        = bookDQMHistogram1D(detid,"ScurveChi2ProbSummary", detIdName + " Chi2Prob Summary", 101, 0, 1.01);
00290       MonitorElement * D1FitResult   = bookDQMHistogram1D(detid,"ScurveFitResultSummary", detIdName + " Fit Result Summary", 10, -0.5, 9.5);
00291       insertResult.first->second.insert(std::make_pair(kSigmaSummary, D1sigma));
00292       insertResult.first->second.insert(std::make_pair(kThresholdSummary, D1thresh));
00293       insertResult.first->second.insert(std::make_pair(kChi2Summary, D1chi2));
00294       insertResult.first->second.insert(std::make_pair(kFitResultSummary, D1FitResult));
00295    }
00296 }
00297 
00298 bool SiPixelSCurveCalibrationAnalysis::doFits(uint32_t detid, std::vector<SiPixelCalibDigi>::const_iterator calibDigi)
00299 {
00300    sCurveErrorFlag errorFlag = errOK;
00301    uint32_t nVCalPts = calibDigi->getnpoints();
00302    //reset and fill static datamembers with vector of points and errors
00303    efficiencies_.resize(0);
00304    effErrors_.resize(0);
00305    for (uint32_t iVcalPt = 0; iVcalPt < nVCalPts; iVcalPt++)
00306    {
00307       float eff;
00308       float error;
00309       calculateEffAndError(calibDigi->getnentries(iVcalPt), nTriggers_, eff, error);
00310       edm::LogInfo("SiPixelSCurveCalibrationAnalysis") << "Eff: " << eff << " Error:  " << error << "  nEntries: " << calibDigi->getnentries(iVcalPt) << "  nTriggers: " << nTriggers_ << " VCalPt " << vCalValues_[iVcalPt];
00311       efficiencies_.push_back(eff);
00312       effErrors_.push_back(error);
00313    } 
00314 
00315    //estimate the S-Curve parameters
00316    float thresholdGuess = -1.0;
00317    float sigmaGuess = -1.0;
00318    errorFlag = estimateSCurveParameters(efficiencies_, thresholdGuess, sigmaGuess);
00319 
00320    // these -1.0 default values will only be filled if the curve is all zeroes, or doesn't turn on, WHICH INDICATES A SERIOUS PROBLEM
00321    Double_t sigma                       = -1.0;
00322    Double_t sigmaError                  = -1.0;
00323    Double_t threshold                   = -1.0;
00324    Double_t thresholdError              = -1.0;
00325    Double_t amplitude                   = -1.0;
00326    Double_t amplitudeError              = -1.0;
00327    Double_t chi2                        = -1.0;
00328    //calculate NDF
00329    Int_t nDOF                           = vCalValues_.size() - 3;
00330    Double_t chi2probability             = 0;
00331 
00332    if (errorFlag == errOK)          //only do fit if curve is fittable
00333    {
00334       //set up minuit fit
00335       TMinuit *gMinuit = new TMinuit(3);
00336       gMinuit->SetPrintLevel(-1);  //save ourselves from gigabytes of stdout
00337       gMinuit->SetFCN(chi2toMinimize);
00338 
00339       //define threshold parameters - choose step size 1, max 300, min -50
00340       gMinuit->DefineParameter(0, "Threshold", (Double_t)thresholdGuess, 1, -50, 300);
00341       //sigma
00342       gMinuit->DefineParameter(1, "Sigma", (Double_t)sigmaGuess, 0.1, 0, 255); 
00343       //amplitude
00344       gMinuit->DefineParameter(2, "Amplitude", 1, 0.1, -0.001, 200);
00345 
00346       //Do Chi2 minimazation
00347       gMinuit->Migrad();
00348       gMinuit->GetParameter(0, threshold, thresholdError);
00349       gMinuit->GetParameter(1, sigma, sigmaError);
00350       gMinuit->GetParameter(2, amplitude, amplitudeError);
00351 
00352       //get Chi2
00353       Double_t params[3]   = {threshold, sigma, amplitude};
00354       gMinuit->Eval(3, NULL, chi2, params, 0);
00355       //calculate Chi2 proability
00356       if (nDOF <= 0)
00357          chi2probability = 0;
00358       else
00359          chi2probability = TMath::Prob(chi2, nDOF);
00360       
00361       //check to make sure output makes sense (i.e. threshold > 0)
00362       if (chi2probability > minimumChi2prob_)
00363          errorFlag = fittedSCurveSanityCheck(threshold, sigma, amplitude);
00364       else
00365          errorFlag = errBadChi2Prob;
00366 
00367       edm::LogInfo("SiPixelSCurveCalibrationAnalysis") << "Fit finished with errorFlag: " << errorFlag << " - threshold: " << threshold << "  sigma: " << sigma << "  chi2: " << chi2 << "  nDOF: " << nDOF << " chi2Prob: " << chi2probability << " chi2MinUser: " << minimumChi2prob_;
00368 
00369       delete gMinuit;
00370    }
00371    //get row and column for this pixel
00372    uint32_t row = calibDigi->row();
00373    uint32_t col = calibDigi->col();
00374 
00375    //get iterator to histogram holder for this detid
00376    detIDHistogramMap::iterator thisDetIdHistoGrams;
00377    thisDetIdHistoGrams = histograms_.find(detid);
00378    if (thisDetIdHistoGrams != histograms_.end())
00379    {
00380       edm::LogInfo("SiPixelSCurveCalibrationAnalysisHistogramReport") << "Filling histograms for [detid](col/row):  [" << detid << "](" << col << "/" << row << ") ErrorFlag: " << errorFlag;
00381       //always fill fit result
00382       (*thisDetIdHistoGrams).second[kFitResultSummary]->Fill(errorFlag);
00383       if (write2dFitResult_)
00384          (*thisDetIdHistoGrams).second[kFitResults]->setBinContent(col+1, row+1, errorFlag); // +1 because root bins start at 1
00385 
00386       // fill sigma/threshold result
00387       (*thisDetIdHistoGrams).second[kSigmaSummary]->Fill(sigma);
00388       (*thisDetIdHistoGrams).second[kThresholdSummary]->Fill(threshold);
00389       if (write2dHistograms_)
00390       {
00391          (*thisDetIdHistoGrams).second[kSigmas]->setBinContent(col+1, row+1, sigma); // +1 because root bins start at 1
00392          (*thisDetIdHistoGrams).second[kThresholds]->setBinContent(col+1, row+1, threshold); // +1 because root bins start at 1
00393       }
00394       // fill chi2
00395       (*thisDetIdHistoGrams).second[kChi2Summary]->Fill(chi2probability);
00396       if (write2dHistograms_)
00397          (*thisDetIdHistoGrams).second[kChi2s]->Fill(col, row, chi2probability);
00398    }
00399    // save individual curves, if requested
00400    if (saveCurvesThatFlaggedBad_)
00401    {
00402       bool thisDetIDinList = false;
00403       if (detIDsToSave_.find(detid) != detIDsToSave_.end()) //see if we want to save this histogram
00404          thisDetIDinList = true;
00405 
00406       if (errorFlag != errOK || thisDetIDinList)
00407       {
00408          edm::LogError("SiPixelSCurveCalibrationAnalysis") << "Saving error histogram for [detid](col/row):  [" << detid << "](" << col << "/" << row << ") ErrorFlag: " << errorFlag;
00409          buildACurveHistogram(detid, row, col, errorFlag, efficiencies_, effErrors_);
00410       }
00411    }
00412 
00413    return true;
00414    
00415 }