Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:59:34

0001 // -*- C++ -*-
0002 //
0003 // Package:    SiPixelGainCalibrationAnalysis
0004 // Class:      SiPixelGainCalibrationAnalysis
0005 //
0006 /**\class SiPixelGainCalibrationAnalysis SiPixelGainCalibrationAnalysis.cc CalibTracker/SiPixelGainCalibrationAnalysis/src/SiPixelGainCalibrationAnalysis.cc
0007 
0008 Description: <one line class summary>
0009 
0010 Implementation:
0011 <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Freya Blekman
0015 //         Created:  Wed Nov 14 15:02:06 CET 2007
0016 //
0017 //
0018 
0019 // user include files
0020 #include "SiPixelGainCalibrationAnalysis.h"
0021 #include <sstream>
0022 #include <vector>
0023 #include <cmath>
0024 #include "TGraphErrors.h"
0025 #include "TMath.h"
0026 
0027 using std::cout;
0028 using std::endl;
0029 //
0030 // constructors and destructor
0031 //
0032 SiPixelGainCalibrationAnalysis::SiPixelGainCalibrationAnalysis(const edm::ParameterSet &iConfig)
0033     : SiPixelOfflineCalibAnalysisBase(iConfig),
0034       conf_(iConfig),
0035       bookkeeper_(),
0036       bookkeeper_pixels_(),
0037       nfitparameters_(iConfig.getUntrackedParameter<int>("numberOfFitParameters", 2)),
0038       fitfunction_(iConfig.getUntrackedParameter<std::string>("fitFunctionRootFormula", "pol1")),
0039       listofdetids_(conf_.getUntrackedParameter<std::vector<uint32_t> >("listOfDetIDs")),
0040       ignoreMode_(conf_.getUntrackedParameter<bool>("ignoreMode", false)),
0041       reject_plateaupoints_(iConfig.getUntrackedParameter<bool>("suppressPlateauInFit", true)),
0042       reject_single_entries_(iConfig.getUntrackedParameter<bool>("suppressPointsWithOneEntryOrLess", true)),
0043       plateau_max_slope_(iConfig.getUntrackedParameter<double>("plateauSlopeMax", 1.0)),
0044       reject_first_point_(iConfig.getUntrackedParameter<bool>("rejectVCalZero", true)),
0045       reject_badpoints_frac_(iConfig.getUntrackedParameter<double>("suppressZeroAndPlateausInFitFrac", 0)),
0046       bookBIGCalibPayload_(iConfig.getUntrackedParameter<bool>("saveFullPayloads", false)),
0047       savePixelHists_(iConfig.getUntrackedParameter<bool>("savePixelLevelHists", false)),
0048       chi2Threshold_(iConfig.getUntrackedParameter<double>("minChi2NDFforHistSave", 10)),
0049       chi2ProbThreshold_(iConfig.getUntrackedParameter<double>("minChi2ProbforHistSave", 0.05)),
0050       maxGainInHist_(iConfig.getUntrackedParameter<double>("maxGainInHist", 10)),
0051       maxChi2InHist_(iConfig.getUntrackedParameter<double>("maxChi2InHist", 25)),
0052       saveALLHistograms_(iConfig.getUntrackedParameter<bool>("saveAllHistograms", false)),
0053 
0054       filldb_(iConfig.getUntrackedParameter<bool>("writeDatabase", false)),
0055       writeSummary_(iConfig.getUntrackedParameter<bool>("writeSummary", true)),
0056       recordName_(conf_.getParameter<std::string>("record")),
0057 
0058       appendMode_(conf_.getUntrackedParameter<bool>("appendMode", true)),
0059       /*theGainCalibrationDbInput_(0),
0060   theGainCalibrationDbInputOffline_(0),
0061   theGainCalibrationDbInputHLT_(0),
0062   theGainCalibrationDbInputService_(iConfig),*/
0063       gainlow_(10.),
0064       gainhi_(0.),
0065       pedlow_(255.),
0066       pedhi_(0.),
0067       useVcalHigh_(conf_.getParameter<bool>("useVCALHIGH")),
0068       scalarVcalHigh_VcalLow_(conf_.getParameter<double>("vcalHighToLowConversionFac")) {
0069   if (reject_single_entries_)
0070     min_nentries_ = 1;
0071   else
0072     min_nentries_ = 0;
0073   ::putenv((char *)"CORAL_AUTH_USER=me");
0074   ::putenv((char *)"CORAL_AUTH_PASSWORD=test");
0075   edm::LogInfo("SiPixelGainCalibrationAnalysis") << "now using fit function " << fitfunction_ << ", which has "
0076                                                  << nfitparameters_ << " free parameters. " << std::endl;
0077   func_ = new TF1("func", fitfunction_.c_str(), 0, 256 * scalarVcalHigh_VcalLow_);
0078   graph_ = new TGraphErrors();
0079   currentDetID_ = 0;
0080   summary_.open("SummaryPerDetID.txt");
0081   statusNumbers_ = new int[10];
0082   for (int ii = 0; ii < 10; ii++)
0083     statusNumbers_[ii] = 0;
0084 }
0085 
0086 SiPixelGainCalibrationAnalysis::~SiPixelGainCalibrationAnalysis() {}
0087 // member functions
0088 //
0089 // ------------ method called once each job just before starting event loop  ------------
0090 
0091 std::vector<float> SiPixelGainCalibrationAnalysis::CalculateAveragePerColumn(uint32_t detid, std::string label) {
0092   std::vector<float> result;
0093   int ncols = bookkeeper_[detid][label]->getNbinsX();
0094   int nrows = bookkeeper_[detid][label]->getNbinsY();
0095   for (int icol = 1; icol <= ncols; ++icol) {
0096     float val = 0;
0097     float ntimes = 0;
0098     for (int irow = 1; irow <= nrows; ++irow) {
0099       val += bookkeeper_[detid][label]->getBinContent(icol, irow);
0100       ntimes++;
0101     }
0102     val /= ntimes;
0103     result.push_back(val);
0104   }
0105   return result;
0106 }
0107 
0108 bool SiPixelGainCalibrationAnalysis::checkCorrectCalibrationType() {
0109   if (calibrationMode_ == "GainCalibration")
0110     return true;
0111   else if (ignoreMode_ == true)
0112     return true;
0113   else if (calibrationMode_ == "unknown") {
0114     edm::LogInfo("SiPixelGainCalibrationAnalysis")
0115         << "calibration mode is: " << calibrationMode_ << ", continuing anyway...";
0116     return true;
0117   } else {
0118     //    edm::LogError("SiPixelGainCalibrationAnalysis") << "unknown calibration mode for Gain calibration, should be \"Gain\" and is \"" << calibrationMode_ << "\"";
0119   }
0120   return false;
0121 }
0122 
0123 void SiPixelGainCalibrationAnalysis::calibrationSetup(const edm::EventSetup &) {}
0124 //------- summary printing method. Very verbose.
0125 void SiPixelGainCalibrationAnalysis::printSummary() {
0126   uint32_t detid = 0;
0127   for (std::map<uint32_t, std::map<std::string, MonitorElement *> >::const_iterator idet = bookkeeper_.begin();
0128        idet != bookkeeper_.end();
0129        ++idet) {
0130     if (detid == idet->first)
0131       continue;  // only do things once per detid
0132     detid = idet->first;
0133     std::vector<float> gainvec = CalculateAveragePerColumn(detid, "gain_2d");
0134     std::vector<float> pedvec = CalculateAveragePerColumn(detid, "ped_2d");
0135     std::vector<float> chi2vec = CalculateAveragePerColumn(detid, "chi2_2d");
0136     std::ostringstream summarytext;
0137 
0138     summarytext << "Summary for det ID " << detid << "(" << translateDetIdToString(detid) << ")\n";
0139     summarytext << "\t Following: values per column: column #, gain, pedestal, chi2\n";
0140     for (uint32_t i = 0; i < gainvec.size(); i++)
0141       summarytext << "\t " << i << " \t" << gainvec[i] << " \t" << pedvec[i] << " \t" << chi2vec[i] << "\n";
0142     summarytext << "\t list of pixels with high chi2 (chi2> " << chi2Threshold_ << "): \n";
0143 
0144     for (std::map<std::string, MonitorElement *>::const_iterator ipix = bookkeeper_pixels_[detid].begin();
0145          ipix != bookkeeper_pixels_[detid].end();
0146          ++ipix)
0147       summarytext << "\t " << ipix->first << "\n";
0148     edm::LogInfo("SiPixelGainCalibrationAnalysis") << summarytext.str() << std::endl;
0149   }
0150   if (summary_.is_open()) {
0151     summary_.close();
0152     summary_.open("Summary.txt");
0153     summary_ << "Total Number of Pixel computed :" << statusNumbers_[9] << endl;
0154     summary_ << "Number of pixel tagged with status :" << endl;
0155     for (int ii = 0; ii < 9; ii++)
0156       summary_ << ii << " -> " << statusNumbers_[ii] << " ~ "
0157                << double(statusNumbers_[ii]) / double(statusNumbers_[9]) * 100. << " %" << endl;
0158 
0159     summary_.close();
0160   }
0161 }
0162 
0163 // ------------ method called once each job just after ending the event loop  ------------
0164 
0165 void SiPixelGainCalibrationAnalysis::calibrationEnd() {
0166   if (writeSummary_)
0167     printSummary();
0168 
0169   // this is where we loop over all histograms and save the database objects
0170   if (filldb_)
0171     fillDatabase();
0172 }
0173 //-----------method to fill the database
0174 void SiPixelGainCalibrationAnalysis::fillDatabase() {
0175   // only create when necessary.
0176   // process the minimum and maximum gain & ped values...
0177   edm::LogError("SiPixelGainCalibration::fillDatabase()")
0178       << "PLEASE do not fill the database directly from the gain calibration analyzer. This function is currently "
0179          "disabled and no DB payloads will be produced!"
0180       << std::endl;
0181 }
0182 // ------------ method called to do fits to all objects available  ------------
0183 bool SiPixelGainCalibrationAnalysis::doFits(uint32_t detid, std::vector<SiPixelCalibDigi>::const_iterator ipix) {
0184   float lowmeanval = 255;
0185   float highmeanval = 0;
0186   bool makehistopersistent = saveALLHistograms_;
0187   std::vector<uint32_t>::const_iterator detidfinder = find(listofdetids_.begin(), listofdetids_.end(), detid);
0188   if (detidfinder != listofdetids_.end())
0189     makehistopersistent = true;
0190   // first, fill the input arrays to the TLinearFitter.
0191   double xvals[257];
0192   double yvals[256];
0193   double yerrvals[256];
0194   double xvalsall[257];
0195   double yvalsall[256];
0196   double yerrvalsall[256];
0197   int npoints = 0;
0198   int nallpoints = 0;
0199   bool use_point = true;
0200   int status = 0;
0201   statusNumbers_[9]++;
0202 
0203   bookkeeper_[detid]["status_2d"]->setBinContent(ipix->col() + 1, ipix->row() + 1, 0);
0204   if (writeSummary_ && detid != currentDetID_) {
0205     currentDetID_ = detid;
0206     summary_ << endl << "DetId_" << currentDetID_ << endl;
0207   }
0208 
0209   for (uint32_t ii = 0; ii < ipix->getnpoints() && ii < 200; ii++) {
0210     //    std::cout << ipix->getsum(ii) << " " << ipix->getnentries(ii) << " " << ipix->getsumsquares(ii) << std::endl;
0211     nallpoints++;
0212     use_point = true;
0213     if (useVcalHigh_) {
0214       xvalsall[ii] = vCalValues_[ii] * scalarVcalHigh_VcalLow_;
0215     } else
0216       xvalsall[ii] = vCalValues_[ii];
0217     yerrvalsall[ii] = yvalsall[ii] = 0;
0218 
0219     if (ipix->getnentries(ii) > min_nentries_) {
0220       yvalsall[ii] = ipix->getsum(ii) / (float)ipix->getnentries(ii);
0221       yerrvalsall[ii] = ipix->getsumsquares(ii) / (float)(ipix->getnentries(ii));
0222       yerrvalsall[ii] -= pow(yvalsall[ii], 2);
0223       yerrvalsall[ii] = sqrt(yerrvalsall[ii]) / sqrt(ipix->getnentries(ii));
0224 
0225       if (yvalsall[ii] < lowmeanval)
0226         lowmeanval = yvalsall[ii];
0227       if (yvalsall[ii] > highmeanval)
0228         highmeanval = yvalsall[ii];
0229     }
0230   }
0231 
0232   // calculate plateau value from last 4 entries
0233   double plateauval = 0;
0234   bool noPlateau = false;
0235   if (nallpoints >= 4) {
0236     for (int ii = nallpoints - 1; ii > nallpoints - 5; --ii)
0237       plateauval += yvalsall[ii];
0238     plateauval /= 4;
0239     for (int ii = nallpoints - 1; ii > nallpoints - 5; --ii) {
0240       if (fabs(yvalsall[ii] - plateauval) > 5) {
0241         plateauval = 255;
0242         noPlateau = true;
0243         continue;
0244       }
0245     }
0246 
0247     int NbofPointsInPlateau = 0;
0248     for (int ii = 0; ii < nallpoints; ++ii)
0249       if (fabs(yvalsall[ii] - plateauval) < 10 || yvalsall[ii] == 0)
0250         NbofPointsInPlateau++;
0251     //summary_<<"row_"<<ipix->row()<<" col_"<<ipix->col()<<"   "<<plateauval<<"  "<<NbofPointsInPlateau<<"  "<<nallpoints<<endl;
0252     if (NbofPointsInPlateau >= (nallpoints - 2)) {
0253       status = -2;
0254       bookkeeper_[detid]["status_2d"]->setBinContent(ipix->col() + 1, ipix->row() + 1, status);
0255       if (writeSummary_) {
0256         summary_ << "row_" << ipix->row() << " col_" << ipix->col() << " status_" << status << endl;
0257         statusNumbers_[abs(status)]++;
0258       }
0259       return false;
0260     }
0261   } else
0262     plateauval = 255;
0263 
0264   double maxgoodvalinfit = plateauval * (1. - reject_badpoints_frac_);
0265   npoints = 0;
0266   for (int ii = 0; ii < nallpoints; ++ii) {
0267     // now selecting the appropriate points for the fit.
0268     use_point = true;
0269     if (reject_first_point_ && xvalsall[ii] < 0.1)
0270       use_point = false;
0271     if (ipix->getnentries(ii) <= min_nentries_ && reject_single_entries_)
0272       use_point = false;
0273     if (ipix->getnentries(ii) == 0 && reject_badpoints_)
0274       use_point = false;
0275     if (yvalsall[ii] > maxgoodvalinfit && !noPlateau)
0276       use_point = false;
0277     if (ii > 1 && fabs(yvalsall[ii] - yvalsall[ii - 1]) < 5. && yvalsall[ii] > 0.8 * maxgoodvalinfit &&
0278         reject_plateaupoints_) {
0279       use_point = false;
0280       break;
0281     }
0282 
0283     if (use_point) {
0284       xvals[npoints] = xvalsall[ii];
0285       yvals[npoints] = yvalsall[ii];
0286       yerrvals[npoints] = yerrvalsall[ii];
0287       npoints++;
0288     }
0289   }
0290 
0291   float chi2, slope, intercept, prob, slopeerror, intercepterror;
0292   prob = chi2 = -1;
0293   slope = intercept = slopeerror = intercepterror = 0;
0294 
0295   // now check on number of points. If bad just start taking the first 4:
0296 
0297   if (npoints < 4) {
0298     npoints = 0;
0299     for (int ii = 0; ii < nallpoints && npoints < 4 && yvalsall[ii] < plateauval * 0.97; ++ii) {
0300       if (yvalsall[ii] > 0) {
0301         if (ii > 0 && yvalsall[ii] - yvalsall[ii - 1] < 0.1)
0302           continue;
0303         xvals[npoints] = xvalsall[ii];
0304         yvals[npoints] = yvalsall[ii];
0305         yerrvals[npoints] = yerrvalsall[ii];
0306         npoints++;
0307       }
0308     }
0309   }
0310   if (npoints < 2) {
0311     status = -7;
0312     bookkeeper_[detid]["status_2d"]->setBinContent(ipix->col() + 1, ipix->row() + 1, status);
0313     if (writeSummary_) {
0314       summary_ << "row_" << ipix->row() << " col_" << ipix->col() << " status_" << status << endl;
0315       statusNumbers_[abs(status)]++;
0316     }
0317     std::ostringstream pixelinfo;
0318     pixelinfo << "GainCurve_row_" << ipix->row() << "_col_" << ipix->col();
0319     std::string tempname = translateDetIdToString(detid);
0320     tempname += "_";
0321     tempname += pixelinfo.str();
0322     setDQMDirectory(detid);
0323     bookkeeper_pixels_[detid][pixelinfo.str()] = bookDQMHistogram1D(
0324         detid, pixelinfo.str(), tempname, 105 * nallpoints, xvalsall[0], xvalsall[nallpoints - 1] * 1.05);
0325     for (int ii = 0; ii < nallpoints; ++ii)
0326       bookkeeper_pixels_[detid][pixelinfo.str()]->Fill(xvalsall[ii], yvalsall[ii]);
0327     return false;
0328   }
0329 
0330   //  std::cout << "starting fit!" << std::endl;
0331   graph_->Set(npoints);
0332 
0333   func_->SetParameter(0, 50.);
0334   func_->SetParameter(1, 0.25);
0335   for (int ipointtemp = 0; ipointtemp < npoints; ++ipointtemp) {
0336     graph_->SetPoint(ipointtemp, xvals[ipointtemp], yvals[ipointtemp]);
0337     graph_->SetPointError(ipointtemp, 0, yerrvals[ipointtemp]);
0338   }
0339   Int_t tempresult = graph_->Fit(func_, "FQ0N");
0340   slope = func_->GetParameter(1);
0341   slopeerror = func_->GetParError(1);
0342   intercept = func_->GetParameter(0);
0343   intercepterror = func_->GetParError(0);
0344   chi2 = func_->GetChisquare() / ((float)npoints - func_->GetNpar());
0345   prob = TMath::Prob(func_->GetChisquare(), npoints - func_->GetNpar());
0346   size_t ntimes = 0;
0347   while ((std::isnan(slope) || std::isnan(intercept)) && ntimes < 10) {
0348     ntimes++;
0349     makehistopersistent = true;
0350     //    std::cout << slope << " " << intercept << " " << prob << std::endl;
0351     edm::LogWarning("SiPixelGainCalibrationAnalysis") << "impossible to fit values, try " << ntimes << ": ";
0352     for (int ii = 0; ii < npoints; ++ii) {
0353       edm::LogWarning("SiPixelGainCalibrationAnalysis")
0354           << "vcal " << xvals[ii] << " response: " << yvals[ii] << "+/-" << yerrvals[ii] << std::endl;
0355     }
0356     tempresult = graph_->Fit(func_, "FQ0NW");
0357     slope = func_->GetParameter(1);
0358     slopeerror = func_->GetParError(1);
0359     intercept = func_->GetParameter(0);
0360     intercepterror = func_->GetParError(0);
0361     chi2 = func_->GetChisquare() / ((float)npoints - func_->GetNpar());
0362     prob = TMath::Prob(func_->GetChisquare(), npoints - func_->GetNpar());
0363   }
0364 
0365   if (tempresult == 0)
0366     status = 1;
0367   else
0368     status = 0;
0369   if (slope != 0)
0370     slope = 1. / slope;
0371   if (std::isnan(slope) || std::isnan(intercept)) {
0372     status = -6;
0373     bookkeeper_[detid]["status_2d"]->setBinContent(ipix->col() + 1, ipix->row() + 1, status);
0374     if (writeSummary_) {
0375       summary_ << "row_" << ipix->row() << " col_" << ipix->col() << " status_" << status << endl;
0376       statusNumbers_[abs(status)]++;
0377     }
0378     //return false;
0379   }
0380   if (chi2 > chi2Threshold_ && chi2Threshold_ >= 0)
0381     status = 5;
0382   if (prob < chi2ProbThreshold_)
0383     status = 5;
0384   if (noPlateau)
0385     status = 3;
0386   if (nallpoints < 4)
0387     status = -7;
0388   if (TMath::Abs(slope > maxGainInHist_) || slope < 0)
0389     status = -8;
0390   if (status != 1)
0391     makehistopersistent = true;
0392   statusNumbers_[abs(status)]++;
0393 
0394   if (slope < gainlow_)
0395     gainlow_ = slope;
0396   if (slope > gainhi_)
0397     gainhi_ = slope;
0398   if (intercept > pedhi_)
0399     pedhi_ = intercept;
0400   if (intercept < pedlow_)
0401     pedlow_ = intercept;
0402   bookkeeper_[detid]["gain_1d"]->Fill(slope);
0403   if (slope > maxGainInHist_) {
0404     makehistopersistent = true;
0405     edm::LogWarning("SiPixelGainCalibration")
0406         << "For DETID " << detid << "pixel row,col " << ipix->row() << "," << ipix->col() << " Gain was measured to be "
0407         << slope << " which is outside the range of the summary plot (" << maxGainInHist_ << ") !!!! " << std::endl;
0408   }
0409   bookkeeper_[detid]["dynamicrange_2d"]->setBinContent(ipix->col() + 1, ipix->row() + 1, highmeanval - lowmeanval);
0410   bookkeeper_[detid]["plateau_2d"]->setBinContent(ipix->col() + 1, ipix->row() + 1, highmeanval);
0411   bookkeeper_[detid]["gain_2d"]->setBinContent(ipix->col() + 1, ipix->row() + 1, slope);
0412   bookkeeper_[detid]["errorgain_2d"]->setBinContent(ipix->col() + 1, ipix->row() + 1, slopeerror);
0413   bookkeeper_[detid]["ped_1d"]->Fill(intercept);
0414   bookkeeper_[detid]["ped_2d"]->setBinContent(ipix->col() + 1, ipix->row() + 1, intercept);
0415   bookkeeper_[detid]["errorped_2d"]->setBinContent(ipix->col() + 1, ipix->row() + 1, intercepterror);
0416   bookkeeper_[detid]["chi2_1d"]->Fill(chi2);
0417   bookkeeper_[detid]["chi2_2d"]->setBinContent(ipix->col() + 1, ipix->row() + 1, chi2);
0418   bookkeeper_[detid]["prob_1d"]->Fill(prob);
0419   bookkeeper_[detid]["prob_2d"]->setBinContent(ipix->col() + 1, ipix->row() + 1, prob);
0420   bookkeeper_[detid]["lowpoint_1d"]->Fill(xvals[0]);
0421   bookkeeper_[detid]["lowpoint_2d"]->setBinContent(ipix->col() + 1, ipix->row() + 1, xvals[0]);
0422   bookkeeper_[detid]["highpoint_1d"]->Fill(xvals[npoints - 1]);
0423   bookkeeper_[detid]["highpoint_2d"]->setBinContent(ipix->col() + 1, ipix->row() + 1, xvals[npoints - 1]);
0424   bookkeeper_[detid]["nfitpoints_1d"]->Fill(npoints);
0425   bookkeeper_[detid]["endpoint_1d"]->Fill((255 - intercept) * slope);
0426   bookkeeper_[detid]["status_2d"]->setBinContent(ipix->col() + 1, ipix->row() + 1, status);
0427 
0428   if (!savePixelHists_)
0429     return true;
0430   if (detidfinder == listofdetids_.end() && !listofdetids_.empty())
0431     return true;
0432   if (makehistopersistent) {
0433     std::ostringstream pixelinfo;
0434     pixelinfo << "GainCurve_row_" << ipix->row() << "_col_" << ipix->col();
0435     std::string tempname = translateDetIdToString(detid);
0436     tempname += "_";
0437     tempname += pixelinfo.str();
0438 
0439     // and book the histo
0440     // fill the last value of the vcal array...
0441 
0442     setDQMDirectory(detid);
0443     bookkeeper_pixels_[detid][pixelinfo.str()] = bookDQMHistogram1D(
0444         detid, pixelinfo.str(), tempname, 105 * nallpoints, xvalsall[0], xvalsall[nallpoints - 1] * 1.05);
0445 
0446     edm::LogInfo("SiPixelGainCalibrationAnalysis")
0447         << "now saving histogram for pixel " << tempname << ", gain = " << slope << ", pedestal = " << intercept
0448         << ", chi2/NDF=" << chi2 << "(prob:" << prob << "), fit status " << status;
0449     for (int ii = 0; ii < nallpoints; ++ii) {
0450       //      std::cout << xvalsall[ii]<<","<<yvalsall[ii]<< " " << tempfloats[ii] << std::endl;
0451       bookkeeper_pixels_[detid][pixelinfo.str()]->Fill(xvalsall[ii], yvalsall[ii]);
0452     }
0453 
0454     //    addTF1ToDQMMonitoringElement(bookkeeper_pixels_[detid][pixelinfo.str()],func_);
0455 
0456     if (writeSummary_) {
0457       summary_ << "row_" << ipix->row() << " col_" << ipix->col();
0458       summary_ << " status_" << status;
0459       summary_ << endl;
0460 
0461       //std::cout<<detid<<"  "<<"row " <<ipix->row()<<" col "<<ipix->col()<<"  "<<status<<"  "<<chi2<<"  "<<prob<<"  "<<npoints<<"  "<<xvals[0]<<"  "<<xvals[npoints-1]<<"  "<<plateauval<<std::endl;
0462     }
0463   }
0464   return true;
0465 }
0466 // ------------ method called to do fill new detids  ------------
0467 void SiPixelGainCalibrationAnalysis::newDetID(uint32_t detid) {
0468   setDQMDirectory(detid);
0469   std::string tempname = translateDetIdToString(detid);
0470   bookkeeper_[detid]["gain_1d"] = bookDQMHistogram1D(detid, "Gain1d", "gain for " + tempname, 100, 0., maxGainInHist_);
0471   bookkeeper_[detid]["gain_2d"] = bookDQMHistoPlaquetteSummary2D(detid, "Gain2d", "gain for " + tempname);
0472   bookkeeper_[detid]["errorgain_2d"] =
0473       bookDQMHistoPlaquetteSummary2D(detid, "ErrorGain2d", "error on gain for " + tempname);
0474   bookkeeper_[detid]["ped_1d"] = bookDQMHistogram1D(detid, "Pedestal1d", "pedestal for " + tempname, 256, 0., 256.0);
0475   bookkeeper_[detid]["ped_2d"] = bookDQMHistoPlaquetteSummary2D(detid, "Pedestal2d", "pedestal for " + tempname);
0476   bookkeeper_[detid]["errorped_2d"] =
0477       bookDQMHistoPlaquetteSummary2D(detid, "ErrorPedestal2d", "error on pedestal for " + tempname);
0478   bookkeeper_[detid]["chi2_1d"] =
0479       bookDQMHistogram1D(detid, "GainChi2NDF1d", "#chi^{2}/NDOF for " + tempname, 100, 0., maxChi2InHist_);
0480   bookkeeper_[detid]["chi2_2d"] =
0481       bookDQMHistoPlaquetteSummary2D(detid, "GainChi2NDF2d", "#chi^{2}/NDOF for " + tempname);
0482   bookkeeper_[detid]["prob_1d"] =
0483       bookDQMHistogram1D(detid, "GainChi2Prob1d", "P(#chi^{2},NDOF) for " + tempname, 100, 0., 1.0);
0484   bookkeeper_[detid]["prob_2d"] =
0485       bookDQMHistoPlaquetteSummary2D(detid, "GainChi2Prob2d", "P(#chi^{2},NDOF) for " + tempname);
0486   bookkeeper_[detid]["status_2d"] =
0487       bookDQMHistoPlaquetteSummary2D(detid, "GainFitResult2d", "Fit result for " + tempname);
0488   bookkeeper_[detid]["endpoint_1d"] = bookDQMHistogram1D(
0489       detid, "GainEndPoint1d", "point where fit meets ADC=255 for " + tempname, 256, 0., 256. * scalarVcalHigh_VcalLow_);
0490   bookkeeper_[detid]["lowpoint_1d"] = bookDQMHistogram1D(
0491       detid, "GainLowPoint1d", "lowest fit point for " + tempname, 256, 0., 256. * scalarVcalHigh_VcalLow_);
0492   bookkeeper_[detid]["highpoint_1d"] = bookDQMHistogram1D(
0493       detid, "GainHighPoint1d", "highest fit point for " + tempname, 256, 0., 256. * scalarVcalHigh_VcalLow_);
0494   bookkeeper_[detid]["nfitpoints_1d"] =
0495       bookDQMHistogram1D(detid, "GainNPoints1d", "number of fit point for " + tempname, 20, 0., 20);
0496   bookkeeper_[detid]["dynamicrange_2d"] = bookDQMHistoPlaquetteSummary2D(
0497       detid, "GainDynamicRange2d", "Difference lowest and highest points on gain curve for " + tempname);
0498   bookkeeper_[detid]["lowpoint_2d"] =
0499       bookDQMHistoPlaquetteSummary2D(detid, "GainLowPoint2d", "lowest fit point for " + tempname);
0500   bookkeeper_[detid]["highpoint_2d"] =
0501       bookDQMHistoPlaquetteSummary2D(detid, "GainHighPoint2d", "highest fit point for " + tempname);
0502   bookkeeper_[detid]["plateau_2d"] =
0503       bookDQMHistoPlaquetteSummary2D(detid, "GainSaturate2d", "Highest points on gain curve for " + tempname);
0504 }
0505 //define this as a plug-in
0506 DEFINE_FWK_MODULE(SiPixelGainCalibrationAnalysis);